]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDBackgroundCorrection.cxx
No need to liknk with lhapdf, pythia6 and microcern libraries
[u/mrichter/AliRoot.git] / FMD / AliFMDBackgroundCorrection.cxx
1 /**************************************************************************
2  * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
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.
19 //
20 // Author: Hans Hjersing Dalsgaard, NBI, hans.dalsgaard@cern.ch
21 //
22 //
23
24 #include "AliSimulation.h"
25 #include "TStopwatch.h"
26 #include "iostream"
27 #include "TGrid.h"
28 #include "AliRunLoader.h"
29 #include "AliGeomManager.h"
30 #include "AliFMDGeometry.h"
31 #include "AliStack.h"
32 #include "TParticle.h"
33 #include "TDatabasePDG.h"
34 #include "TParticlePDG.h"
35 #include "AliRun.h"
36 #include "AliFMDBackgroundCorrection.h"
37 #include "TSystem.h"
38 #include "AliCDBManager.h"
39 #include "TTree.h"
40 #include "TClonesArray.h"
41 #include "TBranch.h"
42 #include "AliFMDHit.h"
43 #include "AliLoader.h" 
44 #include "AliFMD.h"
45 #include "TH2F.h"
46 #include "AliGenEventHeader.h"
47 #include "AliHeader.h"
48 #include "TFile.h"
49 #include "TAxis.h"
50 #include "AliCDBId.h"
51 #include "AliCDBMetaData.h"
52 #include "TROOT.h"
53 #include "AliFMDParameters.h"
54 #include "AliLog.h"
55
56
57
58 ClassImp(AliFMDBackgroundCorrection)
59 //_____________________________________________________________________
60 AliFMDBackgroundCorrection::AliFMDBackgroundCorrection() : 
61   TNamed(),
62   fCorrectionArray()
63 {} 
64
65 //_____________________________________________________________________
66 AliFMDBackgroundCorrection::AliFMDInputBG::AliFMDInputBG() : 
67   AliFMDInput(),
68   fPrim(0),
69   fHits(0),
70   fZvtxCut(0),
71   fNvtxBins(0),
72   fPrevTrack(-1),
73   fPrevDetector(-1),
74   fPrevRing('Q'),
75   fNbinsEta(100)
76 {
77   AddLoad(kTracks); 
78   AddLoad(kHits); 
79   AddLoad(kKinematics); 
80   AddLoad(kHeader);
81 }
82
83 //_____________________________________________________________________
84
85 void AliFMDBackgroundCorrection::GenerateBackgroundCorrection(Int_t nvtxbins,
86                                                               Float_t zvtxcut, 
87                                                               Int_t nBinsEta, 
88                                                               Bool_t storeInAlien, 
89                                                               Int_t runNo,
90                                                               Int_t endRunNo, 
91                                                               const Char_t* filename, 
92                                                               Bool_t simulate,
93                                                               Int_t nEvents) {
94   
95   TGrid::Connect("alien:",0,0,"t");
96   if(simulate)
97     Simulate(nEvents);
98   else {
99     AliCDBManager::Instance()->SetDefaultStorage("alien://Folder=/alice/data/2008/LHC08d/OCDB/");
100     AliCDBManager::Instance()->SetRun(runNo);
101     
102 #if defined(__CINT__)
103     gSystem->Load("liblhapdf");
104     gSystem->Load("libEGPythia6");
105     gSystem->Load("libpythia6");
106     gSystem->Load("libAliPythia6");
107     gSystem->Load("libgeant321");
108 #endif
109   }  
110   
111   //Setting up the geometry
112   //-----------------------------------------------
113   if (AliGeomManager::GetGeometry() == NULL)
114     AliGeomManager::LoadGeometry();
115   
116   AliFMDGeometry* geo = AliFMDGeometry::Instance();
117   geo->Init();
118   geo->InitTransformations();
119   
120   
121   
122   AliInfo("Processing hits and primaries ");
123   
124   AliFMDInputBG input;
125   input.SetVtxCutZ(zvtxcut);
126   input.SetNvtxBins(nvtxbins);
127   input.SetNbinsEta(nBinsEta);
128   input.Run();
129   
130   AliInfo(Form("Found %d primaries and %d hits.", input.GetNprim(),input.GetNhits()));
131   
132   TObjArray* hitArray = input.GetHits();
133   TObjArray* primaryArray = input.GetPrimaries();
134   fCorrectionArray.SetName("FMD_bg_correction");
135   fCorrectionArray.SetOwner();
136   for(Int_t det= 1; det <=3; det++) {
137     Int_t nRings = (det==1 ? 1 : 2);
138     
139     TObjArray* detArrayCorrection = new TObjArray();
140     detArrayCorrection->SetName(Form("FMD%d",det));
141     fCorrectionArray.AddAtAndExpand(detArrayCorrection,det);
142     
143     
144     for(Int_t iring = 0; iring<nRings; iring++) {
145       TObjArray* primRingArray = (TObjArray*)primaryArray->At(iring);
146       Char_t ringChar = (iring == 0 ? 'I' : 'O');
147       TObjArray* vtxArrayCorrection = new TObjArray();
148       vtxArrayCorrection->SetName(Form("FMD%d%c",det,ringChar));
149       detArrayCorrection->AddAtAndExpand(vtxArrayCorrection,iring);
150       
151       for(Int_t vertexBin=0;vertexBin<nvtxbins;vertexBin++) {
152         TObjArray* detArray  = (TObjArray*)hitArray->At(det);
153         TObjArray* vtxArray  = (TObjArray*)detArray->At(iring);
154         TH2F* hHits          = (TH2F*)vtxArray->At(vertexBin);
155         TH2F* hPrimary  = (TH2F*)primRingArray->At(vertexBin);
156         TH2F* hCorrection = (TH2F*)hHits->Clone(Form("FMD%d%c_vtxbin_%d_correction",det,ringChar,vertexBin));
157         hCorrection->Divide(hPrimary);
158         vtxArrayCorrection->AddAtAndExpand(hCorrection,vertexBin);
159       }
160       
161     }
162   }
163   
164   TAxis refAxis(nvtxbins,-1*zvtxcut,zvtxcut);
165   
166   
167
168   
169   TFile fout(filename,"RECREATE");
170   refAxis.Write("vertexbins");
171   fout.mkdir("Hits");
172   fout.cd("Hits");
173   hitArray->Write();
174   fout.mkdir("Primaries");
175   fout.cd("Primaries");
176   primaryArray->Write();
177   fout.mkdir("Correction");
178   fout.cd("Correction");
179   fCorrectionArray.Write();
180   
181   if(storeInAlien) {
182     AliCDBManager* cdb = AliCDBManager::Instance();
183     cdb->SetDefaultStorage("local://$ALICE_ROOT");
184     AliCDBId      id("FMD/AnalysisCalib/Background",runNo,endRunNo);
185     
186     AliCDBMetaData* meta = new AliCDBMetaData;                          
187     meta->SetResponsible(gSystem->GetUserInfo()->fRealName.Data());     
188     meta->SetAliRootVersion(gROOT->GetVersion());                       
189     meta->SetBeamPeriod(1);                                             
190     meta->SetComment("Background Correction for FMD");
191     meta->SetProperty("key1", &fout );
192     cdb->Put(&fout, id, meta);
193     
194   }
195   
196
197   fout.Close();
198   
199 }
200 //_____________________________________________________________________
201 void AliFMDBackgroundCorrection::Simulate(Int_t nEvents) {
202   
203   AliSimulation sim ; 
204   sim.SetRunNumber(0);
205   TGrid::Connect("alien:",0,0,"t");
206   sim.SetDefaultStorage("alien://Folder=/alice/data/2008/LHC08d/OCDB/");
207   sim.SetConfigFile("Config.C");
208   sim.SetRunQA("FMD:");
209   TStopwatch timer;
210   timer.Start();
211   sim.RunSimulation(nEvents);    
212   timer.Stop();
213   timer.Print();
214   
215 }
216
217 //_____________________________________________________________________
218 Bool_t AliFMDBackgroundCorrection::AliFMDInputBG::ProcessTrack(Int_t i , TParticle* p, AliFMDHit* h) {
219   
220  
221   //  if(!h)
222   //  return kTRUE;
223   //if(!p)
224   //  return kTRUE;
225   
226   // if(h) {
227   //  if(h->Detector() == 3)
228   //    std::cout<<"Track "<<i<<" hit "<<"FMD3"<<h->Ring()<<std::endl;
229     
230   //}
231   
232   AliGenEventHeader* genHeader = fLoader->GetHeader()->GenEventHeader();
233   TArrayF vertex;
234   genHeader->PrimaryVertex(vertex);
235   
236   if(TMath::Abs(vertex.At(2)) > fZvtxCut) 
237     return kTRUE;
238   
239   Double_t delta = 2*fZvtxCut/fNvtxBins;
240   Double_t vertexBinDouble = (vertex.At(2) + fZvtxCut) / delta;
241   Int_t vertexBin = (Int_t)vertexBinDouble;
242   
243   if(h)
244     if(h->Q() !=  0 && (i != fPrevTrack || h->Detector() != fPrevDetector || h->Ring() != fPrevRing)) {
245       
246       Int_t det = h->Detector();
247       Char_t ring = h->Ring();
248       //  if(h->Detector() == 3)
249       //        std::cout<<"Detected a hit in " <<"FMD"<<det<<ring<<"    ! "<<std::endl;
250      
251       Int_t iring = (ring == 'I' ? 0 : 1);
252       
253       TObjArray* detArray  = (TObjArray*)fHitArray.At(det);
254       TObjArray* vtxArray  = (TObjArray*)detArray->At(iring);
255       TH2F* hHits          = (TH2F*)vtxArray->At(vertexBin);
256       
257       Float_t phi   = TMath::ATan2(h->Py(),h->Px());
258       if(phi<0)
259         phi = phi+2*TMath::Pi();
260       Float_t theta = TMath::ATan2(TMath::Sqrt(TMath::Power(h->X(),2)+TMath::Power(h->Y(),2)),h->Z()+vertex.At(2));
261       Float_t eta   = -1*TMath::Log(TMath::Tan(0.5*theta));
262       hHits->Fill(eta,phi);
263       fHits++;
264       //fPrevTrack = i;
265       fPrevDetector = det;
266       fPrevRing     = ring;
267     }
268   if(p && i != fPrevTrack) {
269     
270     TDatabasePDG* pdgDB = TDatabasePDG::Instance();
271     TParticlePDG* pdgPart = pdgDB->GetParticle(p->GetPdgCode());
272     Float_t charge = (pdgPart ? pdgPart->Charge() : 0);
273     Float_t phi = TMath::ATan2(p->Py(),p->Px());
274
275     if(phi<0)
276       phi = phi+2*TMath::Pi();
277     
278     
279     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);
280     if(primary) {
281       fPrim++;
282       TObjArray* innerArray = (TObjArray*)fPrimaryArray.At(0);
283       TObjArray* outerArray = (TObjArray*)fPrimaryArray.At(1);
284       
285       TH2F* hPrimaryInner  = (TH2F*)innerArray->At(vertexBin);
286       TH2F* hPrimaryOuter  = (TH2F*)outerArray->At(vertexBin);
287       hPrimaryInner->Fill(p->Eta(),phi);
288       hPrimaryOuter->Fill(p->Eta(),phi);
289     }
290   }
291   
292   
293   fPrevTrack = i;
294   return kTRUE;
295 }
296
297 //_____________________________________________________________________
298 Bool_t AliFMDBackgroundCorrection::AliFMDInputBG::Init() {
299   
300   fPrimaryArray.SetOwner();
301   fPrimaryArray.SetName("FMD_primary");
302   
303   for(Int_t iring = 0; iring<2;iring++) {
304     Char_t ringChar = (iring == 0 ? 'I' : 'O');
305     TObjArray* ringArray = new TObjArray();
306     ringArray->SetName(Form("FMD_%c",ringChar));
307     fPrimaryArray.AddAtAndExpand(ringArray,iring);
308     Int_t nSec = (iring == 1 ? 40 : 20);
309     for(Int_t v=0; v<fNvtxBins;v++) {
310
311       TH2F* hPrimary       = new TH2F(Form("hPrimary_FMD_%c_vtx%d",ringChar,v),
312                                       Form("hPrimary_FMD_%c_vtx%d",ringChar,v),
313                                       fNbinsEta, -6,6, nSec, 0,2*TMath::Pi());
314       ringArray->AddAtAndExpand(hPrimary,v);
315     }
316   }
317   
318   
319   fHitArray.SetOwner();
320   fHitArray.SetName("FMD_hits");
321    
322   for(Int_t det =1; det<=3;det++)
323     {
324       TObjArray* detArrayHits = new TObjArray();
325       detArrayHits->SetName(Form("FMD%d",det));
326       fHitArray.AddAtAndExpand(detArrayHits,det);
327       Int_t nRings = (det==1 ? 1 : 2);
328       for(Int_t ring = 0;ring<nRings;ring++)
329         {
330           Int_t nSec = (ring == 1 ? 40 : 20);
331           Char_t ringChar = (ring == 0 ? 'I' : 'O');
332           TObjArray* vtxArrayHits = new TObjArray();
333           vtxArrayHits->SetName(Form("FMD%d%c",det,ringChar));
334           detArrayHits->AddAtAndExpand(vtxArrayHits,ring);
335           for(Int_t v=0; v<fNvtxBins;v++)
336             {
337               
338               TH2F* hHits          = new TH2F(Form("hHits_FMD%d%c_vtx%d",det,ringChar,v),
339                                               Form("hHits_FMD%d%c_vtx%d",det,ringChar,v),
340                                               fNbinsEta, -6,6, nSec, 0, 2*TMath::Pi());
341               vtxArrayHits->AddAtAndExpand(hHits,v);
342                       
343             } 
344         }
345       
346     }
347
348
349   
350   AliFMDInput::Init();
351   
352   
353   return kTRUE;
354 }
355 // EOF
356
357
358