]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/analysis/AliFMDAnalysisTaskCollector.cxx
9e38c6f68f674545645a8f1bf5cf8c2587d5fc0f
[u/mrichter/AliRoot.git] / FMD / analysis / AliFMDAnalysisTaskCollector.cxx
1  
2 #include <TROOT.h>
3 #include <TSystem.h>
4 #include <TInterpreter.h>
5 #include <TChain.h>
6 #include <TFile.h>
7 #include <TList.h>
8 #include <iostream>
9
10 #include "AliFMDAnalysisTaskCollector.h"
11 #include "AliAnalysisManager.h"
12 #include "AliESDFMD.h"
13 #include "AliESDEvent.h"
14 #include "AliAODEvent.h"
15 #include "AliAODHandler.h"
16 #include "AliMCEventHandler.h"
17 #include "AliStack.h"
18 #include "AliESDVertex.h"
19 #include "AliFMDAnaParameters.h"
20
21 ClassImp(AliFMDAnalysisTaskCollector)
22
23
24 AliFMDAnalysisTaskCollector::AliFMDAnalysisTaskCollector()
25 : fDebug(0),
26   fChain(0x0),
27   fESD(0x0),
28   fOutputList(0),
29   fArray(0),
30   fEdistHist(0),
31   fZvtxDist(0)
32 {
33   // Default constructor
34   DefineInput (0, TChain::Class());
35   DefineOutput(0, TList::Class());
36 }
37 //____________________________________________________________________
38 AliFMDAnalysisTaskCollector::AliFMDAnalysisTaskCollector(const char* name):
39     AliAnalysisTask(name, "AnalysisTaskFMD"),
40     fDebug(0),
41     fChain(0x0),
42     fESD(0x0),
43     fOutputList(0),
44     fArray(0),
45     fEdistHist(0),
46     fZvtxDist(0)
47 {
48   // Default constructor
49   DefineInput (0, TChain::Class());
50   DefineOutput(0, TList::Class());
51 }
52 //____________________________________________________________________
53 void AliFMDAnalysisTaskCollector::CreateOutputObjects()
54 {
55   // Create the output container
56   printf("AnalysisTaskFMD::CreateOutPutData() \n");
57   
58   fOutputList = new TList();//(TList*)GetOutputData(0);
59   AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
60   
61   fArray     = new TObjArray();
62   fArray->SetName("FMD");
63   fArray->SetOwner();
64   TH1F* hEdist = 0;
65  
66   for(Int_t nEta = 0; nEta <= pars->GetNetaBins()+1; nEta++) {
67     TObjArray* etaArray = new TObjArray();
68     fArray->AddAtAndExpand(etaArray,nEta);
69     for(Int_t det =1; det<=3;det++)
70       {
71         TObjArray* detArray = new TObjArray();
72         detArray->SetName(Form("FMD%d",det));
73         etaArray->AddAtAndExpand(detArray,det);
74         Int_t nRings = (det==1 ? 1 : 2);
75         for(Int_t ring = 0;ring<nRings;ring++)
76           {
77             Char_t ringChar = (ring == 0 ? 'I' : 'O');
78             hEdist = new TH1F(Form("FMD%d%c_etabin%d",det,ringChar,nEta),Form("FMD%d%c_etabin%d",det,ringChar,nEta),200,0,6);
79             hEdist->SetXTitle("#Delta E / E_{MIP}");
80             fOutputList->Add(hEdist);
81             detArray->AddAtAndExpand(hEdist,ring);
82           } 
83       }
84     
85   }
86   
87   fZvtxDist  = new TH1F("ZvtxDist","Vertex distribution",100,-30,30);
88   fZvtxDist->SetXTitle("z vertex");
89   //fOutputList->Add(fArray);
90   fOutputList->Add(fZvtxDist);
91 }
92 //____________________________________________________________________
93 void AliFMDAnalysisTaskCollector::Init()
94 {
95   // Initialization
96   printf("AnalysisTaskFMD::Init() \n");
97   
98
99 }
100 //____________________________________________________________________
101 void AliFMDAnalysisTaskCollector::ConnectInputData(Option_t */*option*/)
102 {
103   fChain = (TChain*)GetInputData(0);
104   fESD = new AliESDEvent();
105   fESD->ReadFromTree(fChain);
106   
107 }
108 //____________________________________________________________________
109 void AliFMDAnalysisTaskCollector::Exec(Option_t */*option*/)
110 {
111   AliESD* old = fESD->GetAliESDOld();
112   if (old) {
113     fESD->CopyFromOldESD();
114   }
115   AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
116   
117   Double_t vertex[3];
118   
119   GetVertex(vertex);
120   //if(vertex[0] == 0 && vertex[1] == 0 && vertex[2] == 0)
121   //  return;
122   fZvtxDist->Fill(vertex[2]);
123   
124   if(TMath::Abs(vertex[2]) > pars->GetVtxCutZ())
125     return;
126   
127   AliESDFMD* fmd = fESD->GetFMDData();
128   if (!fmd) return;
129   
130   
131   
132   for(UShort_t det=1;det<=3;det++) {
133     
134     
135     Int_t nRings = (det==1 ? 1 : 2);
136     for (UShort_t ir = 0; ir < nRings; ir++) {
137   
138       Char_t   ring = (ir == 0 ? 'I' : 'O');
139       UShort_t nsec = (ir == 0 ? 20  : 40);
140       UShort_t nstr = (ir == 0 ? 512 : 256);
141       TH2F* hBg = pars->GetBackgroundCorrection(det,ring,0);
142       
143       for(UShort_t sec =0; sec < nsec;  sec++)  {
144         for(UShort_t strip = 0; strip < nstr; strip++) {
145           Float_t eta = fmd->Eta(det,ring,sec,strip);
146           Int_t nEta = hBg->GetXaxis()->FindBin(eta);
147          
148           TObjArray* etaArray = (TObjArray*)fArray->At(nEta);
149           TObjArray* detArray = (TObjArray*)etaArray->At(det);
150           TH1F* Edist = (TH1F*)detArray->At(ir);
151           Float_t mult = fmd->Multiplicity(det,ring,sec,strip);
152           if(mult == AliESDFMD::kInvalidMult || mult == 0) continue;
153           Edist->Fill(mult);
154           
155         }
156       }
157     }
158   }
159   
160   PostData(0, fOutputList); 
161   
162 }
163 //____________________________________________________________________
164 void AliFMDAnalysisTaskCollector::Terminate(Option_t */*option*/)
165 {
166   /*
167   for(UShort_t det=1;det<=3;det++) {
168     TObjArray* detArray = (TObjArray*)fArray->At(det);
169     Int_t nRings = (det==1 ? 1 : 2);
170     for (UShort_t ir = 0; ir < nRings; ir++) {
171       TH1F* hEdist = (TH1F*)detArray->At(ir);
172       hEdist->SetAxisRange(0.4,hEdist->GetXaxis()->GetXmax());
173       Float_t max = hEdist->GetBinCenter(hEdist->GetMaximumBin());
174       hEdist->Fit("landau","","",max-0.1,2*max);
175     }
176   }
177   */
178 }
179 //_____________________________________________________________________
180 void AliFMDAnalysisTaskCollector::GetVertex(Double_t* vertexXYZ) 
181 {
182   const AliESDVertex* vertex = 0;
183   vertex = fESD->GetPrimaryVertex();
184   if (!vertex || (vertexXYZ[0] == 0 && vertexXYZ[1] == 0 && vertexXYZ[2] == 0))        
185     vertex = fESD->GetPrimaryVertexSPD();
186   if (!vertex || (vertexXYZ[0] == 0 && vertexXYZ[1] == 0 && vertexXYZ[2] == 0))        
187     vertex = fESD->GetPrimaryVertexTPC();
188   
189   if (!vertex || (vertexXYZ[0] == 0 && vertexXYZ[1] == 0 && vertexXYZ[2] == 0))    
190     vertex = fESD->GetVertex();
191   
192   if (vertex && (vertexXYZ[0] != 0 || vertexXYZ[1] != 0 || vertexXYZ[2] != 0)) {
193     vertex->GetXYZ(vertexXYZ);
194     //std::cout<<vertex->GetName()<<"   "<< vertex->GetTitle() <<"   "<< vertex->GetZv()<<std::endl;
195     return;
196   }
197   else if (fESD->GetESDTZERO()) { 
198     vertexXYZ[0] = 0;
199     vertexXYZ[1] = 0;
200     vertexXYZ[2] = fESD->GetT0zVertex();
201     
202     return;
203   }
204   
205   return;
206 }
207 //____________________________________________________________________
208 //
209 // EOF
210 //