1 /**************************************************************************
2 * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 // Thil class computes background corrections for the FMD. The correction is computed
17 // in eta,phi cells and the objects stored can be put into alien to use with analysis.
18 // It is based on the AliFMDInput class that is used to loop over hits and primaries.
20 // Author: Hans Hjersing Dalsgaard, NBI, hans.dalsgaard@cern.ch
24 #include "AliSimulation.h"
25 #include "TStopwatch.h"
28 #include "AliRunLoader.h"
29 #include "AliGeomManager.h"
30 #include "AliFMDGeometry.h"
32 #include "TParticle.h"
33 #include "TDatabasePDG.h"
34 #include "TParticlePDG.h"
36 #include "AliFMDBackgroundCorrection.h"
38 #include "AliCDBManager.h"
40 #include "TClonesArray.h"
42 #include "AliFMDHit.h"
43 #include "AliLoader.h"
46 #include "AliGenEventHeader.h"
47 #include "AliHeader.h"
51 #include "AliCDBMetaData.h"
53 #include "AliFMDParameters.h"
55 #include "AliFMDAnaParameters.h"
56 #include "AliFMDAnaCalibBackgroundCorrection.h"
58 ClassImp(AliFMDBackgroundCorrection)
59 //_____________________________________________________________________
60 AliFMDBackgroundCorrection::AliFMDBackgroundCorrection() :
65 //_____________________________________________________________________
66 AliFMDBackgroundCorrection::AliFMDInputBG::AliFMDInputBG() :
85 //_____________________________________________________________________
87 void AliFMDBackgroundCorrection::GenerateBackgroundCorrection(Int_t nvtxbins,
93 const Char_t* filename,
97 TGrid::Connect("alien:",0,0,"t");
101 AliCDBManager::Instance()->SetDefaultStorage("alien://Folder=/alice/data/2008/LHC08d/OCDB/");
102 AliCDBManager::Instance()->SetRun(runNo);
104 #if defined(__CINT__)
105 gSystem->Load("liblhapdf");
106 gSystem->Load("libEGPythia6");
107 gSystem->Load("libpythia6");
108 gSystem->Load("libAliPythia6");
109 gSystem->Load("libgeant321");
113 //Setting up the geometry
114 //-----------------------------------------------
115 if (AliGeomManager::GetGeometry() == NULL)
116 AliGeomManager::LoadGeometry();
118 AliFMDGeometry* geo = AliFMDGeometry::Instance();
120 geo->InitTransformations();
124 AliInfo("Processing hits and primaries ");
127 input.SetVtxCutZ(zvtxcut);
128 input.SetNvtxBins(nvtxbins);
129 input.SetNbinsEta(nBinsEta);
132 AliInfo(Form("Found %d primaries and %d hits.", input.GetNprim(),input.GetNhits()));
134 TObjArray* hitArray = input.GetHits();
135 TObjArray* primaryArray = input.GetPrimaries();
136 fCorrectionArray.SetName("FMD_bg_correction");
137 fCorrectionArray.SetOwner();
139 AliFMDAnaCalibBackgroundCorrection* background = new AliFMDAnaCalibBackgroundCorrection();
141 for(Int_t det= 1; det <=3; det++) {
142 Int_t nRings = (det==1 ? 1 : 2);
144 TObjArray* detArrayCorrection = new TObjArray();
145 detArrayCorrection->SetName(Form("FMD%d",det));
146 fCorrectionArray.AddAtAndExpand(detArrayCorrection,det);
149 for(Int_t iring = 0; iring<nRings; iring++) {
150 TObjArray* primRingArray = (TObjArray*)primaryArray->At(iring);
151 Char_t ringChar = (iring == 0 ? 'I' : 'O');
152 TObjArray* vtxArrayCorrection = new TObjArray();
153 vtxArrayCorrection->SetName(Form("FMD%d%c",det,ringChar));
154 detArrayCorrection->AddAtAndExpand(vtxArrayCorrection,iring);
156 for(Int_t vertexBin=0;vertexBin<nvtxbins;vertexBin++) {
157 TObjArray* detArray = (TObjArray*)hitArray->At(det);
158 TObjArray* vtxArray = (TObjArray*)detArray->At(iring);
159 TH2F* hHits = (TH2F*)vtxArray->At(vertexBin);
160 TH2F* hPrimary = (TH2F*)primRingArray->At(vertexBin);
161 TH2F* hCorrection = (TH2F*)hHits->Clone(Form("FMD%d%c_vtxbin_%d_correction",det,ringChar,vertexBin));
162 hCorrection->SetTitle(hCorrection->GetName());
163 hCorrection->Divide(hPrimary);
164 vtxArrayCorrection->AddAtAndExpand(hCorrection,vertexBin);
165 background->SetBgCorrection(det,ringChar,vertexBin,hCorrection);
171 TAxis refAxis(nvtxbins,-1*zvtxcut,zvtxcut);
173 background->SetRefAxis(&refAxis);
175 TFile* fout = new TFile(filename,"RECREATE");
176 refAxis.Write("vertexbins");
180 fout->mkdir("Primaries");
181 fout->cd("Primaries");
182 primaryArray->Write();
183 fout->mkdir("Correction");
184 fout->cd("Correction");
185 fCorrectionArray.Write();
187 TObjArray* container = new TObjArray();
188 container->SetOwner();
189 container->AddAtAndExpand(&refAxis,0);
190 container->AddAtAndExpand(&fCorrectionArray,1);
191 container->AddAtAndExpand(hitArray,2);
192 container->AddAtAndExpand(primaryArray,3);
196 AliCDBManager* cdb = AliCDBManager::Instance();
197 cdb->SetDefaultStorage("local://$ALICE_ROOT");
198 AliCDBId id(AliFMDAnaParameters::GetBackgroundPath(),runNo,endRunNo);
200 AliCDBMetaData* meta = new AliCDBMetaData;
201 meta->SetResponsible(gSystem->GetUserInfo()->fRealName.Data());
202 meta->SetAliRootVersion(gROOT->GetVersion());
203 meta->SetBeamPeriod(1);
204 meta->SetComment("Background Correction for FMD");
205 meta->SetProperty("key1", background );
206 cdb->Put(background, id, meta);
214 //_____________________________________________________________________
215 void AliFMDBackgroundCorrection::Simulate(Int_t nEvents) {
219 TGrid::Connect("alien:",0,0,"t");
220 sim.SetDefaultStorage("alien://Folder=/alice/data/2008/LHC08d/OCDB/");
221 sim.SetConfigFile("Config.C");
222 sim.SetRunQA("FMD:");
225 sim.RunSimulation(nEvents);
231 //_____________________________________________________________________
232 Bool_t AliFMDBackgroundCorrection::AliFMDInputBG::ProcessTrack(Int_t i , TParticle* p, AliFMDHit* h) {
234 AliGenEventHeader* genHeader = fLoader->GetHeader()->GenEventHeader();
236 genHeader->PrimaryVertex(vertex);
238 if(TMath::Abs(vertex.At(2)) > fZvtxCut)
241 Double_t delta = 2*fZvtxCut/fNvtxBins;
242 Double_t vertexBinDouble = (vertex.At(2) + fZvtxCut) / delta;
243 Int_t vertexBin = (Int_t)vertexBinDouble;
246 if(h->Q() != 0 && (i != fPrevTrack || h->Detector() != fPrevDetector || h->Ring() != fPrevRing)) {
248 Int_t det = h->Detector();
249 Char_t ring = h->Ring();
250 // if(h->Detector() == 3)
251 // std::cout<<"Detected a hit in " <<"FMD"<<det<<ring<<" ! "<<std::endl;
253 Int_t iring = (ring == 'I' ? 0 : 1);
255 TObjArray* detArray = (TObjArray*)fHitArray.At(det);
256 TObjArray* vtxArray = (TObjArray*)detArray->At(iring);
257 TH2F* hHits = (TH2F*)vtxArray->At(vertexBin);
259 Float_t phi = TMath::ATan2(h->Py(),h->Px());
261 phi = phi+2*TMath::Pi();
262 Float_t theta = TMath::ATan2(TMath::Sqrt(TMath::Power(h->X(),2)+TMath::Power(h->Y(),2)),h->Z()+vertex.At(2));
263 Float_t eta = -1*TMath::Log(TMath::Tan(0.5*theta));
264 hHits->Fill(eta,phi);
270 if(p && i != fPrevTrack) {
272 TDatabasePDG* pdgDB = TDatabasePDG::Instance();
273 TParticlePDG* pdgPart = pdgDB->GetParticle(p->GetPdgCode());
274 Float_t charge = (pdgPart ? pdgPart->Charge() : 0);
275 Float_t phi = TMath::ATan2(p->Py(),p->Px());
278 phi = phi+2*TMath::Pi();
281 Bool_t primary = (charge!=0)&&(TMath::Abs(p->Vx() - vertex.At(0))<0.1)&&(TMath::Abs(p->Vy() - vertex.At(1))<0.1)&&(TMath::Abs(p->Vz() - vertex.At(2))<0.1);
284 TObjArray* innerArray = (TObjArray*)fPrimaryArray.At(0);
285 TObjArray* outerArray = (TObjArray*)fPrimaryArray.At(1);
287 TH2F* hPrimaryInner = (TH2F*)innerArray->At(vertexBin);
288 TH2F* hPrimaryOuter = (TH2F*)outerArray->At(vertexBin);
289 hPrimaryInner->Fill(p->Eta(),phi);
290 hPrimaryOuter->Fill(p->Eta(),phi);
299 //_____________________________________________________________________
300 Bool_t AliFMDBackgroundCorrection::AliFMDInputBG::Init() {
302 fPrimaryArray.SetOwner();
303 fPrimaryArray.SetName("FMD_primary");
305 for(Int_t iring = 0; iring<2;iring++) {
306 Char_t ringChar = (iring == 0 ? 'I' : 'O');
307 TObjArray* ringArray = new TObjArray();
308 ringArray->SetName(Form("FMD_%c",ringChar));
309 fPrimaryArray.AddAtAndExpand(ringArray,iring);
310 Int_t nSec = (iring == 1 ? 40 : 20);
311 for(Int_t v=0; v<fNvtxBins;v++) {
313 TH2F* hPrimary = new TH2F(Form("hPrimary_FMD_%c_vtx%d",ringChar,v),
314 Form("hPrimary_FMD_%c_vtx%d",ringChar,v),
315 fNbinsEta, -6,6, nSec, 0,2*TMath::Pi());
317 ringArray->AddAtAndExpand(hPrimary,v);
322 fHitArray.SetOwner();
323 fHitArray.SetName("FMD_hits");
325 for(Int_t det =1; det<=3;det++)
327 TObjArray* detArrayHits = new TObjArray();
328 detArrayHits->SetName(Form("FMD%d",det));
329 fHitArray.AddAtAndExpand(detArrayHits,det);
330 Int_t nRings = (det==1 ? 1 : 2);
331 for(Int_t ring = 0;ring<nRings;ring++)
333 Int_t nSec = (ring == 1 ? 40 : 20);
334 Char_t ringChar = (ring == 0 ? 'I' : 'O');
335 TObjArray* vtxArrayHits = new TObjArray();
336 vtxArrayHits->SetName(Form("FMD%d%c",det,ringChar));
337 detArrayHits->AddAtAndExpand(vtxArrayHits,ring);
338 for(Int_t v=0; v<fNvtxBins;v++)
341 TH2F* hHits = new TH2F(Form("hHits_FMD%d%c_vtx%d",det,ringChar,v),
342 Form("hHits_FMD%d%c_vtx%d",det,ringChar,v),
343 fNbinsEta, -6,6, nSec, 0, 2*TMath::Pi());
345 vtxArrayHits->AddAtAndExpand(hHits,v);