a04038c11a50a64537219aefac1deac888e245ea
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis / AliFMDAnalysisTaskGenerateCorrection.cxx
1 #include "AliFMDAnalysisTaskGenerateCorrection.h"
2 #include "AliESDEvent.h"
3 #include "iostream"
4 #include "AliESDFMD.h"
5 #include "TH2F.h"
6 #include "AliTrackReference.h"
7 #include "AliStack.h"
8 #include "AliFMDAnaParameters.h"
9 #include "AliFMDStripIndex.h"
10 #include "AliStack.h"
11 #include "AliMCParticle.h"
12 #include "AliMCEvent.h"
13 //#include "AliFMDGeometry.h"
14 #include "TArray.h"
15 #include "AliGenEventHeader.h"
16 #include "AliMultiplicity.h"
17 #include "AliHeader.h"
18 #include "AliFMDAnaCalibBackgroundCorrection.h"
19 #include "AliFMDAnaCalibEventSelectionEfficiency.h"
20 #include "AliGenPythiaEventHeader.h"
21 //#include "AliCDBManager.h"
22 //#include "AliCDBId.h"
23 //#include "AliCDBMetaData.h"
24 #include "TSystem.h"
25 #include "TROOT.h"
26 #include "TAxis.h"
27 ClassImp(AliFMDAnalysisTaskGenerateCorrection)
28
29 //_____________________________________________________________________
30 AliFMDAnalysisTaskGenerateCorrection::AliFMDAnalysisTaskGenerateCorrection():
31 AliAnalysisTaskSE(),
32   fListOfHits(), 
33   fListOfPrimaries(),
34   fListOfCorrection(),
35   fVertexBins(),
36   fLastTrackByStrip(0),
37   fHitsByStrip(0),
38   fZvtxCut(10),
39   fNvtxBins(10),
40   fNbinsEta(200),
41   fBackground(0),
42   fEventSelectionEff(0)
43 {
44   // Default constructor
45 }
46 //_____________________________________________________________________
47 AliFMDAnalysisTaskGenerateCorrection::AliFMDAnalysisTaskGenerateCorrection(const char* name):
48   AliAnalysisTaskSE(name),
49   fListOfHits(), 
50   fListOfPrimaries(),
51   fListOfCorrection(),
52   fVertexBins(),
53   fLastTrackByStrip(0),
54   fHitsByStrip(0),
55   fZvtxCut(10),
56   fNvtxBins(10),
57   fNbinsEta(200),
58   fBackground(0),
59   fEventSelectionEff(0)
60 {
61  
62   DefineOutput(1, TList::Class());
63   DefineOutput(2, TList::Class());
64   DefineOutput(3, TH1F::Class());
65   DefineOutput(4, TList::Class());
66 }
67 //_____________________________________________________________________
68 void AliFMDAnalysisTaskGenerateCorrection::UserCreateOutputObjects()
69 {
70 // Create the output containers
71 //
72   
73   std::cout<<"Creating output objects"<<std::endl;
74   for(Int_t v=0; v<fNvtxBins;v++) {
75     
76     TH2F* hSPDhits       = new TH2F(Form("hSPDhits_vtx%d",v),
77                                     Form("hSPDhits_vtx%d",v),
78                                     fNbinsEta, -4,6, 20, 0,2*TMath::Pi());
79     hSPDhits->Sumw2();
80     fListOfHits.Add(hSPDhits);
81     
82     for(Int_t iring = 0; iring<2;iring++) {
83       Char_t ringChar = (iring == 0 ? 'I' : 'O');
84       Int_t nSec = (iring == 1 ? 40 : 20);
85       
86       TH2F* hPrimary       = new TH2F(Form("hPrimary_FMD_%c_vtx%d",ringChar,v),
87                                       Form("hPrimary_FMD_%c_vtx%d",ringChar,v),
88                                       fNbinsEta, -4,6, nSec, 0,2*TMath::Pi());
89       hPrimary->Sumw2();
90       fListOfPrimaries.Add(hPrimary);
91       
92       
93     }
94   }
95   
96  
97   
98   
99   for(Int_t det =1; det<=3;det++) {
100     Int_t nRings = (det==1 ? 1 : 2);
101     for(Int_t ring = 0;ring<nRings;ring++) {
102       Int_t nSec = (ring == 1 ? 40 : 20);
103       Char_t ringChar = (ring == 0 ? 'I' : 'O');
104       TH1F* doubleHits = new TH1F(Form("DoubleHits_FMD%d%c",det,ringChar),
105                                   Form("DoubleHits_FMD%d%c",det,ringChar),fNbinsEta, -4,6);
106       TH1F* allHits = new TH1F(Form("allHits_FMD%d%c",det,ringChar),
107                                Form("allHits_FMD%d%c",det,ringChar), fNbinsEta, -4,6);
108       
109       doubleHits->Sumw2();
110       allHits->Sumw2();
111       fListOfHits.Add(allHits);
112       fListOfHits.Add(doubleHits);
113         
114       for(Int_t v=0; v<fNvtxBins;v++) {
115         TH2F* hHits = new TH2F(Form("hHits_FMD%d%c_vtx%d", det,ringChar,v),
116                                Form("hHits_FMD%d%c_vtx%d", det,ringChar,v),
117                                fNbinsEta, -4,6, nSec, 0, 2*TMath::Pi());
118         hHits->Sumw2();
119         fListOfHits.Add(hHits);
120         
121       } 
122     }
123   }
124   
125   TH1F* hEventsSelected  = new TH1F("EventsSelected","EventsSelected",fNvtxBins,0,fNvtxBins);
126   TH1F* hEventsAll    = new TH1F("EventsAll","EventsAll",fNvtxBins,0,fNvtxBins);
127   TH1F* hEventsAllNSD    = new TH1F("EventsAllNSD","EventsAllNSD",fNvtxBins,0,fNvtxBins);
128   TH1F* hEventsSelectedVtx  = new TH1F("EventsSelectedVtx","EventsSelectedVtx",fNvtxBins,0,fNvtxBins);  
129   TH1F* hEventsSelectedTrigger  = new TH1F("EventsSelectedTrigger","EventsSelectedTrigger",fNvtxBins,0,fNvtxBins);
130   
131   TH1F* hEventsSelectedNSDVtx  = new TH1F("EventsSelectedNSDVtx","EventsSelectedNSDVtx",fNvtxBins,0,fNvtxBins);  
132   TH1F* hEventsSelectedNSD     = new TH1F("EventsSelectedNSD","EventsSelectedNSD",fNvtxBins,0,fNvtxBins);
133   
134   TH1F* hXvtx = new TH1F("hXvtx","x vertex distribution",100,-2,2);
135   TH1F* hYvtx = new TH1F("hYvtx","y vertex distribution",100,-2,2);
136   TH1F* hZvtx = new TH1F("hZvtx","z vertex distribution",4*fNvtxBins,-4*fZvtxCut,4*fZvtxCut);
137   
138   fListOfPrimaries.Add(hXvtx);
139   fListOfPrimaries.Add(hYvtx);
140   fListOfPrimaries.Add(hZvtx);
141   
142   hEventsSelected->Sumw2();
143   hEventsAll->Sumw2();
144   hEventsAllNSD->Sumw2();
145   
146   fListOfHits.Add(hEventsSelected);
147   fListOfHits.Add(hEventsSelectedVtx);
148   fListOfHits.Add(hEventsSelectedTrigger);
149   fListOfHits.Add(hEventsSelectedNSDVtx);
150   fListOfHits.Add(hEventsSelectedNSD);
151   
152   
153   fListOfPrimaries.Add(hEventsAll);
154   fListOfPrimaries.Add(hEventsAllNSD);
155   
156   for(Int_t v=0; v<fNvtxBins;v++) {
157     
158     for(Int_t iring = 0; iring<2;iring++) {
159       Char_t ringChar = (iring == 0 ? 'I' : 'O');
160       Int_t nSec = (iring == 1 ? 40 : 20);
161       TH2F* hAnalysed       = new TH2F(Form("Analysed_FMD%c_vtx%d",ringChar,v),
162                                        Form("Analysed_FMD%c_vtx%d",ringChar,v),
163                                        fNbinsEta, -4,6, nSec, 0,2*TMath::Pi());
164       
165       hAnalysed->Sumw2();
166       fListOfPrimaries.Add(hAnalysed);
167       
168       TH2F* hAnalysedNSD       = new TH2F(Form("AnalysedNSD_FMD%c_vtx%d",ringChar,v),
169                                           Form("AnalysedNSD_FMD%c_vtx%d",ringChar,v),
170                                           fNbinsEta, -4,6, nSec, 0,2*TMath::Pi());
171       
172       hAnalysedNSD->Sumw2();
173       fListOfPrimaries.Add(hAnalysedNSD);
174       
175       TH2F* hInel           = new TH2F(Form("Inel_FMD%c_vtx%d",ringChar,v),
176                                        Form("Inel_FMD%c_vtx%d",ringChar,v),
177                                        fNbinsEta, -4,6, nSec, 0,2*TMath::Pi());
178       
179       hInel->Sumw2();
180       fListOfPrimaries.Add(hInel);
181       
182       TH2F* hNSD           = new TH2F(Form("NSD_FMD%c_vtx%d",ringChar,v),
183                                       Form("NSD_FMD%c_vtx%d",ringChar,v),
184                                       fNbinsEta, -4,6, nSec, 0,2*TMath::Pi());
185       
186       hNSD->Sumw2();
187       fListOfPrimaries.Add(hNSD);
188       
189     }
190   }
191   
192   fVertexBins.SetName("VertexBins");
193   fVertexBins.GetXaxis()->Set(fNvtxBins,-1*fZvtxCut,fZvtxCut);
194   
195 }
196 //_____________________________________________________________________
197 void AliFMDAnalysisTaskGenerateCorrection::Init()
198 {
199   fLastTrackByStrip.Reset(-1);
200   
201   
202 }
203 //_____________________________________________________________________
204 void AliFMDAnalysisTaskGenerateCorrection::UserExec(Option_t */*option*/)
205 {
206   
207   fLastTrackByStrip.Reset(-1);
208   fHitsByStrip.Reset(0);
209   AliMCEvent* mcevent = MCEvent();
210   
211   AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
212   
213   
214   AliESDEvent* esdevent = (AliESDEvent*)InputEvent();
215   
216   pars->SetTriggerStatus(esdevent);
217   
218   Double_t esdvertex[3];
219   Bool_t vtxStatus =  pars->GetVertex(esdevent,esdvertex);
220   
221   AliMCParticle* particle = 0;
222   AliStack* stack = mcevent->Stack();
223   
224   UShort_t det,sec,strip;
225   Char_t   ring;
226   
227   Int_t nTracks                = mcevent->GetNumberOfTracks();
228   AliHeader* header            = mcevent->Header();
229   AliGenEventHeader* genHeader = header->GenEventHeader();
230   
231   AliGenPythiaEventHeader* pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(genHeader);
232   if (!pythiaGenHeader) {
233     std::cout<<" no pythia header!"<<std::endl;
234     return; 
235   }
236   
237   Int_t pythiaType = pythiaGenHeader->ProcessType();
238   Bool_t nsd = kTRUE;
239   if(pythiaType==92 || pythiaType==93)
240     nsd = kFALSE;
241   
242   //if(pythiaType==94){
243   //  std::cout<<"double diffractive"<<std::endl;
244   //  return;
245   //}
246   
247   TArrayF vertex;
248   genHeader->PrimaryVertex(vertex);
249   
250   TH1F* hXvtx = (TH1F*)fListOfPrimaries.FindObject("hXvtx");
251   hXvtx->Fill(vertex.At(0));
252   TH1F* hYvtx = (TH1F*)fListOfPrimaries.FindObject("hYvtx");
253   hYvtx->Fill(vertex.At(1));
254   TH1F* hZvtx = (TH1F*)fListOfPrimaries.FindObject("hZvtx");
255   hZvtx->Fill(vertex.At(2));
256
257   
258   
259   
260   Double_t delta           = 2*fZvtxCut/fNvtxBins;
261   Double_t vertexBinDouble = (vertex.At(2) + fZvtxCut) / delta;
262   Int_t    vertexBin       = (Int_t)vertexBinDouble;
263   
264   // Vertex determination correction
265   TH1F* hEventsSelected           = (TH1F*)fListOfHits.FindObject("EventsSelected");
266   TH1F* hEventsSelectedVtx        = (TH1F*)fListOfHits.FindObject("EventsSelectedVtx");
267   TH1F* hEventsSelectedTrigger    = (TH1F*)fListOfHits.FindObject("EventsSelectedTrigger");
268   TH1F* hEventsSelectedNSDVtx     = (TH1F*)fListOfHits.FindObject("EventsSelectedNSDVtx");
269   TH1F* hEventsSelectedNSD        = (TH1F*)fListOfHits.FindObject("EventsSelectedNSD");
270   TH1F* hEventsAll                = (TH1F*)fListOfPrimaries.FindObject("EventsAll");
271   TH1F* hEventsAllNSD             = (TH1F*)fListOfPrimaries.FindObject("EventsAllNSD");
272   
273   // TH1F* hTriggered      = (TH1F*)fListOfHits.FindObject("Triggered");
274   //  TH1F* hTriggeredAll   = (TH1F*)fListOfPrimaries.FindObject("TriggeredAll");
275   
276   Bool_t vtxFound = kTRUE;
277   if(!vtxStatus)
278     vtxFound = kFALSE;
279   
280   if(TMath::Abs(vertex.At(2)) > fZvtxCut) {
281     vtxFound = kFALSE;
282     
283   }
284   Bool_t isTriggered    = pars->IsEventTriggered(AliFMDAnaParameters::kMB1);
285   Bool_t isTriggeredNSD = pars->IsEventTriggered(AliFMDAnaParameters::kNSD);
286   if(vtxFound && isTriggered) hEventsSelected->Fill(vertexBin);
287   
288   if(vtxFound) hEventsSelectedVtx->Fill(vertexBin);
289   if(isTriggered) hEventsSelectedTrigger->Fill(vertexBin);
290
291   if(vtxFound && isTriggeredNSD) hEventsSelectedNSDVtx->Fill(vertexBin);
292   if(isTriggeredNSD) hEventsSelectedNSD->Fill(vertexBin);
293
294   
295   hEventsAll->Fill(vertexBin);
296   if(nsd) hEventsAllNSD->Fill(vertexBin);
297   
298   //  if(!vtxFound || !isTriggered) return;
299   
300   for(Int_t i = 0 ;i<nTracks;i++) {
301     particle = (AliMCParticle*) mcevent->GetTrack(i);
302     
303     if(!particle)
304       continue;
305         
306     if(stack->IsPhysicalPrimary(i) && particle->Charge() != 0) {
307       
308       
309       TH2F* hPrimaryInner = (TH2F*)fListOfPrimaries.FindObject( Form("hPrimary_FMD_%c_vtx%d",'I',vertexBin));
310       TH2F* hPrimaryOuter = (TH2F*)fListOfPrimaries.FindObject( Form("hPrimary_FMD_%c_vtx%d",'O',vertexBin));
311       hPrimaryInner->Fill(particle->Eta(),particle->Phi());
312       hPrimaryOuter->Fill(particle->Eta(),particle->Phi());
313       TH2F* hAnalysedInner = (TH2F*)fListOfPrimaries.FindObject( Form("Analysed_FMD%c_vtx%d",'I',vertexBin));
314       TH2F* hAnalysedOuter = (TH2F*)fListOfPrimaries.FindObject( Form("Analysed_FMD%c_vtx%d",'O',vertexBin));
315       TH2F* hAnalysedNSDInner = (TH2F*)fListOfPrimaries.FindObject( Form("AnalysedNSD_FMD%c_vtx%d",'I',vertexBin));
316       TH2F* hAnalysedNSDOuter = (TH2F*)fListOfPrimaries.FindObject( Form("AnalysedNSD_FMD%c_vtx%d",'O',vertexBin));
317       TH2F* hInelInner = (TH2F*)fListOfPrimaries.FindObject( Form("Inel_FMD%c_vtx%d",'I',vertexBin));
318       TH2F* hInelOuter = (TH2F*)fListOfPrimaries.FindObject( Form("Inel_FMD%c_vtx%d",'O',vertexBin));
319       TH2F* hNSDInner = (TH2F*)fListOfPrimaries.FindObject( Form("NSD_FMD%c_vtx%d",'I',vertexBin));
320       TH2F* hNSDOuter = (TH2F*)fListOfPrimaries.FindObject( Form("NSD_FMD%c_vtx%d",'O',vertexBin));
321       
322       //if(vtxFound && isTriggered) {
323       if(isTriggeredNSD) {
324         hAnalysedNSDInner->Fill(particle->Eta(),particle->Phi());
325         hAnalysedNSDOuter->Fill(particle->Eta(),particle->Phi());
326       }
327       if(isTriggered) {
328         hAnalysedInner->Fill(particle->Eta(),particle->Phi());
329         hAnalysedOuter->Fill(particle->Eta(),particle->Phi());
330       }
331       hInelInner->Fill(particle->Eta(),particle->Phi());
332       hInelOuter->Fill(particle->Eta(),particle->Phi());
333       
334       if(nsd) {
335         hNSDInner->Fill(particle->Eta(),particle->Phi());
336         hNSDOuter->Fill(particle->Eta(),particle->Phi());
337       }
338     }
339     
340     for(Int_t j=0; j<particle->GetNumberOfTrackReferences();j++) {
341       
342       AliTrackReference* ref = particle->GetTrackReference(j);
343       
344       if(ref->DetectorId() != AliTrackReference::kFMD)
345         continue;
346
347       AliFMDStripIndex::Unpack(ref->UserId(),det,ring,sec,strip);
348       Float_t thisStripTrack = fLastTrackByStrip(det,ring,sec,strip);
349       if(particle->Charge() != 0 && i != thisStripTrack ) {
350         
351         Float_t phi = pars->GetPhiFromSector(det,ring,sec);
352         Float_t eta = pars->GetEtaFromStrip(det,ring,sec,strip,vertex.At(2));
353         
354         TH2F* hHits = (TH2F*)fListOfHits.FindObject(Form("hHits_FMD%d%c_vtx%d", det,ring,vertexBin));
355         hHits->Fill(eta,phi);
356         Float_t nstrips = (ring =='O' ? 256 : 512);
357         fHitsByStrip(det,ring,sec,strip) +=1;
358         TH1F* allHits = (TH1F*)fListOfHits.FindObject(Form("allHits_FMD%d%c",det,ring));
359         TH1F* doubleHits = (TH1F*)fListOfHits.FindObject(Form("DoubleHits_FMD%d%c",det,ring));
360         
361         if(fHitsByStrip(det,ring,sec,strip) == 1)
362           allHits->Fill(eta);
363         
364         doubleHits->Fill(eta);
365                 
366         fLastTrackByStrip(det,ring,sec,strip) = (Float_t)i;
367         if(strip >0)
368           fLastTrackByStrip(det,ring,sec,strip-1) = (Float_t)i;
369         if(strip < (nstrips - 1))
370           fLastTrackByStrip(det,ring,sec,strip+1) = (Float_t)i;
371       }
372     }
373
374   }
375   
376   //SPD part HHD
377   TH2F* hSPDMult = (TH2F*)fListOfHits.FindObject(Form("hSPDhits_vtx%d", vertexBin));
378   
379   const AliMultiplicity* spdmult = esdevent->GetMultiplicity();
380   for(Int_t j = 0; j< spdmult->GetNumberOfTracklets();j++) 
381     hSPDMult->Fill(spdmult->GetEta(j),spdmult->GetPhi(j));
382   
383   for(Int_t j = 0; j< spdmult->GetNumberOfSingleClusters();j++) 
384     hSPDMult->Fill(-TMath::Log(TMath::Tan(spdmult->GetThetaSingle(j)/2.)),spdmult->GetPhiSingle(j));
385     
386         
387   PostData(1, &fListOfHits);
388   PostData(2, &fListOfPrimaries);
389   PostData(3, &fVertexBins);
390 }
391 //_____________________________________________________________________
392 void AliFMDAnalysisTaskGenerateCorrection::Terminate(Option_t */*option*/)
393 {
394   /*  TH1F* allHits = (TH1F*)fListOfHits.FindObject("allHits");
395   TH1F* doubleHits = (TH1F*)fListOfHits.FindObject("DoubleHits");
396   
397   doubleHits->Divide(allHits);
398   GenerateCorrection();
399   PostData(1, &fListOfHits);
400   PostData(4, &fListOfCorrection);*/
401   
402 }
403 //_____________________________________________________________________
404 void AliFMDAnalysisTaskGenerateCorrection::GenerateCorrection() {
405   
406   fBackground         = new AliFMDAnaCalibBackgroundCorrection();
407   fEventSelectionEff  = new AliFMDAnaCalibEventSelectionEfficiency();
408   
409   //Event selection
410   TH1F* hEventsSelected           = (TH1F*)fListOfHits.FindObject("EventsSelected");
411   TH1F* hEventsSelectedVtx        = (TH1F*)fListOfHits.FindObject("EventsSelectedVtx");
412   TH1F* hEventsSelectedTrigger    = (TH1F*)fListOfHits.FindObject("EventsSelectedTrigger");
413   TH1F* hEventsSelectedNSDVtx     = (TH1F*)fListOfHits.FindObject("EventsSelectedNSDVtx");
414   TH1F* hEventsSelectedNSD        = (TH1F*)fListOfHits.FindObject("EventsSelectedNSD");
415     
416   TH1F* hEventsAll                = (TH1F*)fListOfPrimaries.FindObject("EventsAll");
417   TH1F* hEventsAllNSD             = (TH1F*)fListOfPrimaries.FindObject("EventsAllNSD");
418   TH1F* hEventsSelectedVtxDivByTr = (TH1F*)hEventsSelectedVtx->Clone("hEventsSelectedVtxDivByTr");
419   TH1F* hEventsSelectedNSDVtxDivByNSD = (TH1F*)hEventsSelectedNSDVtx->Clone("hEventsSelectedNSDVtxDivByNSD");
420   
421   fListOfHits.Add(hEventsSelectedVtxDivByTr);
422   fListOfHits.Add(hEventsSelectedNSDVtxDivByNSD);
423   hEventsSelectedNSDVtxDivByNSD->Divide(hEventsSelectedNSD);
424   hEventsSelectedVtxDivByTr->Divide(hEventsSelectedTrigger);
425   
426   for(Int_t v=0;v<fNvtxBins ; v++) {
427     TH2F* hAnalysedInner = (TH2F*)fListOfPrimaries.FindObject(Form("Analysed_FMD%c_vtx%d",'I',v));
428     TH2F* hAnalysedOuter = (TH2F*)fListOfPrimaries.FindObject(Form("Analysed_FMD%c_vtx%d",'O',v));
429     TH2F* hAnalysedNSDInner = (TH2F*)fListOfPrimaries.FindObject( Form("AnalysedNSD_FMD%c_vtx%d",'I',v));
430     TH2F* hAnalysedNSDOuter = (TH2F*)fListOfPrimaries.FindObject( Form("AnalysedNSD_FMD%c_vtx%d",'O',v));
431     
432     TH2F* hInelInner = (TH2F*)fListOfPrimaries.FindObject(Form("Inel_FMD%c_vtx%d",'I',v));
433     TH2F* hInelOuter = (TH2F*)fListOfPrimaries.FindObject(Form("Inel_FMD%c_vtx%d",'O',v));
434     TH2F* hNSDInner  = (TH2F*)fListOfPrimaries.FindObject( Form("NSD_FMD%c_vtx%d",'I',v));
435     TH2F* hNSDOuter  = (TH2F*)fListOfPrimaries.FindObject( Form("NSD_FMD%c_vtx%d",'O',v));
436     // hAnalysedInner->Scale((hEventsSelected->GetBinContent(v+1) == 0 ? 0 : 1/hEventsSelected->GetBinContent(v+1))); 
437     //hAnalysedOuter->Scale((hEventsSelected->GetBinContent(v+1) == 0 ? 0 : 1/hEventsSelected->GetBinContent(v+1))); 
438
439     hAnalysedInner->Scale((hEventsSelectedTrigger->GetBinContent(v+1) == 0 ? 0 : 1/hEventsSelectedTrigger->GetBinContent(v+1))); 
440     
441     hAnalysedOuter->Scale((hEventsSelectedTrigger->GetBinContent(v+1) == 0 ? 0 : 1/hEventsSelectedTrigger->GetBinContent(v+1))); 
442     hAnalysedNSDInner->Scale((hEventsSelectedNSD->GetBinContent(v+1) == 0 ? 0 : 1/hEventsSelectedNSD->GetBinContent(v+1))); 
443     
444     hAnalysedNSDOuter->Scale((hEventsSelectedNSD->GetBinContent(v+1) == 0 ? 0 : 1/hEventsSelectedNSD->GetBinContent(v+1))); 
445     
446     
447
448     hInelInner->Scale((hEventsAll->GetBinContent(v+1) == 0 ? 0 : 1/hEventsAll->GetBinContent(v+1))); 
449     hInelOuter->Scale((hEventsAll->GetBinContent(v+1) == 0 ? 0 : 1/hEventsAll->GetBinContent(v+1)));
450     
451     hNSDInner->Scale((hEventsAllNSD->GetBinContent(v+1) == 0 ? 0 : 1/hEventsAllNSD->GetBinContent(v+1))); 
452     hNSDOuter->Scale((hEventsAllNSD->GetBinContent(v+1) == 0 ? 0 : 1/hEventsAllNSD->GetBinContent(v+1)));
453     
454
455     
456     //hAnalysedInner->Scale(1./0.84);  //numbers for real data, run104892
457     //hAnalysedOuter->Scale(1./0.84);
458     hAnalysedInner->Divide(hInelInner);
459     hAnalysedOuter->Divide(hInelOuter);
460     
461     hAnalysedNSDInner->Divide(hNSDInner);
462     hAnalysedNSDOuter->Divide(hNSDOuter);
463     
464     fEventSelectionEff->SetCorrection("INEL",v,'I',hAnalysedInner);
465     fEventSelectionEff->SetCorrection("INEL",v,'O',hAnalysedOuter);
466     //NSD
467     fEventSelectionEff->SetCorrection("NSD",v,'I',hAnalysedNSDInner);
468     fEventSelectionEff->SetCorrection("NSD",v,'O',hAnalysedNSDOuter);
469     
470     
471 }
472   
473   Float_t vtxEff = 1;
474   if(hEventsSelectedTrigger->GetEntries())
475     vtxEff = hEventsSelectedVtx->GetEntries() / hEventsSelectedTrigger->GetEntries();
476   fEventSelectionEff->SetVtxToTriggerRatio(vtxEff);
477   
478   //  hEventsAll->Divide(hEventsAll,hEventsSelected,1,1,"B");
479   hEventsSelectedVtx->Divide(hEventsAll);
480   hEventsSelectedTrigger->Divide(hEventsAll);
481   
482   hEventsSelectedNSDVtx->Divide(hEventsAllNSD);
483   hEventsSelectedNSD->Divide(hEventsAllNSD);
484   
485   for(Int_t i = 1; i<=hEventsSelected->GetNbinsX(); i++) {
486     if(hEventsSelected->GetBinContent(i) == 0 )
487       continue;
488     Float_t b    = hEventsSelected->GetBinContent(i);
489     // Float_t b    = hEventsSelected->GetBinContent(i)*hEventsSelectedTrigger->GetBinContent(i);
490     Float_t db   = hEventsSelected->GetBinError(i);
491     Float_t sum  = hEventsAll->GetBinContent(i);
492     Float_t dsum = hEventsAll->GetBinError(i);
493     Float_t a    = sum-b;
494     Float_t da   = TMath::Sqrt(TMath::Power(db,2) + TMath::Power(dsum,2));
495     
496     Float_t cor  = sum / b;
497     Float_t ecor = TMath::Sqrt(TMath::Power(da,2) + TMath::Power(a/(b*db),2)) / b;
498     
499     hEventsAll->SetBinContent(i,cor);
500     hEventsAll->SetBinError(i,ecor);
501     
502   }
503   
504   //TH1F* hEventTest = (TH1F*)hEventsAll->Clone("hEventTest");
505   
506   
507   
508   fEventSelectionEff->SetCorrection(hEventsAll);
509   
510   for(Int_t det= 1; det <=3; det++) {
511     Int_t nRings = (det==1 ? 1 : 2);
512     
513     for(Int_t iring = 0; iring<nRings; iring++) {
514       Char_t ring = (iring == 0 ? 'I' : 'O');
515       TH1F* allHits = (TH1F*)fListOfHits.FindObject(Form("allHits_FMD%d%c",det,ring));
516       TH1F* doubleHits = (TH1F*)fListOfHits.FindObject(Form("DoubleHits_FMD%d%c",det,ring));
517       allHits->Divide(doubleHits);
518       
519       fBackground->SetDoubleHitCorrection(det,ring,allHits);
520       
521       for(Int_t vertexBin=0;vertexBin<fNvtxBins  ;vertexBin++) {
522         TH2F* hHits          = (TH2F*)fListOfHits.FindObject(Form("hHits_FMD%d%c_vtx%d", det,ring,vertexBin));
523         TH2F* hPrimary  = (TH2F*)fListOfPrimaries.FindObject( Form("hPrimary_FMD_%c_vtx%d",ring,vertexBin));
524         TH2F* hCorrection = (TH2F*)hHits->Clone(Form("FMD%d%c_vtxbin_%d_correction",det,ring,vertexBin));
525         hCorrection->Divide(hPrimary);
526
527         
528         hCorrection->SetTitle(hCorrection->GetName());
529         fListOfCorrection.Add(hCorrection);
530         fBackground->SetBgCorrection(det,ring,vertexBin,hCorrection);
531         
532         
533       }
534       
535     }
536   }
537   for(Int_t vertexBin=0;vertexBin<fNvtxBins  ;vertexBin++) {
538     TH2F* hPrimary  = (TH2F*)fListOfPrimaries.FindObject( Form("hPrimary_FMD_%c_vtx%d",'I',vertexBin));
539     TH2F* hSPDMult = (TH2F*)fListOfHits.FindObject(Form("hSPDhits_vtx%d", vertexBin));
540     if(!hSPDMult) continue;
541     
542     TH2F* hCorrection = (TH2F*)hSPDMult->Clone(Form("SPD_vtxbin_%d_correction",vertexBin));
543     hCorrection->SetTitle(hCorrection->GetName());
544     fListOfCorrection.Add(hCorrection);
545     hCorrection->Divide(hPrimary);
546     fBackground->SetBgCorrection(0,'Q',vertexBin,hCorrection);
547     
548     TH1F* hAlive = new TH1F(Form("hAliveSPD_vtxbin%d",vertexBin),Form("hAliveSPD_vtxbin%d",vertexBin),hSPDMult->GetNbinsX(),hSPDMult->GetXaxis()->GetXmin(), hSPDMult->GetXaxis()->GetXmax());
549     TH1F* hPresent = new TH1F(Form("hPresentSPD_vtxbin%d",vertexBin),Form("hPresentSPD_vtxbin%d",vertexBin),hSPDMult->GetNbinsX(),hSPDMult->GetXaxis()->GetXmin(), hSPDMult->GetXaxis()->GetXmax());
550     for(Int_t xx = 1; xx <=hSPDMult->GetNbinsX(); xx++) {
551       
552       if(TMath::Abs(hCorrection->GetXaxis()->GetBinCenter(xx)) > 2)
553         continue;
554       for(Int_t yy = 1; yy <=hSPDMult->GetNbinsY(); yy++) {
555         if(hCorrection->GetBinContent(xx,yy) > 0.1)
556           hAlive->Fill(hCorrection->GetXaxis()->GetBinCenter(xx));
557         hPresent->Fill(hCorrection->GetXaxis()->GetBinCenter(xx));
558         
559       }
560     }
561     TH1F* hDeadCorrection = (TH1F*)hAlive->Clone(Form("hSPDDeadCorrection_vtxbin%d",vertexBin));
562     hDeadCorrection->Divide(hPresent);
563     fBackground->SetSPDDeadCorrection(vertexBin,hDeadCorrection);
564     fListOfCorrection.Add(hDeadCorrection);
565   }
566   
567   
568   
569   
570   
571   
572   TAxis refAxis(fNvtxBins,-1*fZvtxCut,fZvtxCut);
573   fBackground->SetRefAxis(&refAxis);
574
575 }
576 //_____________________________________________________________________
577 void AliFMDAnalysisTaskGenerateCorrection::ReadFromFile(const Char_t* filename, Bool_t storeInOCDB, Int_t /*runNo*/) {
578   
579   TFile infile(filename);
580   TH1F* hVertex = (TH1F*)infile.Get("VertexBins");
581   fZvtxCut = hVertex->GetXaxis()->GetXmax();
582   fNvtxBins = hVertex->GetXaxis()->GetNbins();
583   fVertexBins.SetName("VertexBins");
584   fVertexBins.GetXaxis()->Set(fNvtxBins,-1*fZvtxCut,fZvtxCut);
585   
586   TList* listOfHits = (TList*)infile.Get("Hits");
587   TList* listOfPrim = (TList*)infile.Get("Primaries");
588   
589   TH1F* hEventsSelected           = (TH1F*)listOfHits->FindObject("EventsSelected");
590   TH1F* hEventsSelectedVtx        = (TH1F*)listOfHits->FindObject("EventsSelectedVtx");
591   TH1F* hEventsSelectedTrigger    = (TH1F*)listOfHits->FindObject("EventsSelectedTrigger");
592   TH1F* hEventsSelectedNSDVtx     = (TH1F*)listOfHits->FindObject("EventsSelectedNSDVtx");
593   TH1F* hEventsSelectedNSD        = (TH1F*)listOfHits->FindObject("EventsSelectedNSD");  
594   TH1F* hEventsAll                = (TH1F*)listOfPrim->FindObject("EventsAll");
595   TH1F* hEventsAllNSD             = (TH1F*)listOfPrim->FindObject("EventsAllNSD");
596   
597   fListOfHits.Add(hEventsSelected);
598   fListOfHits.Add(hEventsSelectedVtx);
599   fListOfHits.Add(hEventsSelectedNSD);
600   fListOfHits.Add(hEventsSelectedNSDVtx);
601   fListOfHits.Add(hEventsSelectedTrigger);
602   fListOfPrimaries.Add(hEventsAll);
603   fListOfPrimaries.Add(hEventsAllNSD);
604   
605   TH1F* hXvtx = (TH1F*)listOfPrim->FindObject("hXvtx");
606   TH1F* hYvtx = (TH1F*)listOfPrim->FindObject("hYvtx");
607   TH1F* hZvtx = (TH1F*)listOfPrim->FindObject("hZvtx");
608   fListOfPrimaries.Add(hXvtx);
609   fListOfPrimaries.Add(hYvtx);
610   fListOfPrimaries.Add(hZvtx);  
611
612   for(Int_t det =1; det<=3;det++)
613     {
614       Int_t nRings = (det==1 ? 1 : 2);
615       for(Int_t ring = 0;ring<nRings;ring++)
616         {
617           Char_t ringChar = (ring == 0 ? 'I' : 'O');
618           TH1F* allHits = (TH1F*)listOfHits->FindObject(Form("allHits_FMD%d%c",det,ringChar));
619           TH1F* doubleHits = (TH1F*)listOfHits->FindObject(Form("DoubleHits_FMD%d%c",det,ringChar));
620           fListOfHits.Add(allHits);
621           fListOfHits.Add(doubleHits);
622           for(Int_t v=0; v<fNvtxBins;v++)
623             {
624               
625               TH2F* hHits          = (TH2F*)listOfHits->FindObject(Form("hHits_FMD%d%c_vtx%d", det,ringChar,v));
626               fListOfHits.Add(hHits);
627             }
628         }
629     }
630   for(Int_t v=0; v<fNvtxBins;v++) {
631     TH2F* hSPDHits          = (TH2F*)listOfHits->FindObject(Form("hSPDhits_vtx%d", v));   
632     fListOfHits.Add(hSPDHits);
633     
634     for(Int_t iring = 0; iring<2;iring++) {
635       Char_t ringChar = (iring == 0 ? 'I' : 'O');
636       
637       
638       TH2F* hPrimary       = (TH2F*)listOfPrim->FindObject( Form("hPrimary_FMD_%c_vtx%d",ringChar,v));
639       TH2F* hAnalysed      = (TH2F*)listOfPrim->FindObject(Form("Analysed_FMD%c_vtx%d",ringChar,v));
640       TH2F* hAnalysedNSD   = (TH2F*)listOfPrim->FindObject(Form("AnalysedNSD_FMD%c_vtx%d",ringChar,v));
641       TH2F* hInel          = (TH2F*)listOfPrim->FindObject(Form("Inel_FMD%c_vtx%d",ringChar,v));
642       TH2F* hNSD           = (TH2F*)listOfPrim->FindObject(Form("NSD_FMD%c_vtx%d",ringChar,v));
643       
644       fListOfPrimaries.Add(hPrimary);      
645       fListOfPrimaries.Add(hAnalysed);
646       fListOfPrimaries.Add(hInel);      
647       fListOfPrimaries.Add(hAnalysedNSD);
648       fListOfPrimaries.Add(hNSD);      
649     }
650   }
651   GenerateCorrection();
652   
653   TFile fout("backgroundFromFile.root","recreate");
654   fListOfHits.Write();
655   fListOfPrimaries.Write();
656   fListOfCorrection.Write();
657   fVertexBins.Write();
658   
659   fout.Close();
660   AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
661   if(storeInOCDB) {
662     TFile fbg(pars->GetPath(pars->GetBackgroundID()),"RECREATE");
663     fBackground->Write(AliFMDAnaParameters::GetBackgroundID());
664     fbg.Close();
665     TFile feselect(pars->GetPath(pars->GetEventSelectionEffID()),"RECREATE");
666     fEventSelectionEff->Write(AliFMDAnaParameters::GetEventSelectionEffID());
667     feselect.Close();
668     
669   }
670   
671   
672   
673   
674 }
675 //_____________________________________________________________________
676 //
677 // EOF
678 //