Fix reading from AOD and changed warnings to AliWarning
[u/mrichter/AliRoot.git] / PWG4 / JetTasks / AliAnalysisHelperJetTasks.cxx
1
2 #include "TROOT.h"
3 #include "TList.h"
4 #include "TH1F.h"
5 #include "TProfile.h"
6 #include "THnSparse.h"
7 #include "TFile.h"
8 #include "AliMCEvent.h"
9 #include "AliLog.h"
10 #include "AliAODJet.h"
11 #include "AliStack.h"
12 #include "AliGenEventHeader.h"
13 #include "AliGenCocktailEventHeader.h"
14 #include "AliGenPythiaEventHeader.h"
15 #include <fstream>
16 #include <iostream>
17 #include "AliAnalysisHelperJetTasks.h"
18
19
20 ClassImp(AliAnalysisHelperJetTasks)
21
22
23
24  
25 AliGenPythiaEventHeader*  AliAnalysisHelperJetTasks::GetPythiaEventHeader(AliMCEvent *mcEvent){
26   
27   AliGenEventHeader* genHeader = mcEvent->GenEventHeader();
28   AliGenPythiaEventHeader* pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(genHeader);
29   if(!pythiaGenHeader){
30     // cocktail ??
31     AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast<AliGenCocktailEventHeader*>(genHeader);
32     
33     if (!genCocktailHeader) {
34       AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Unknown header type (not Pythia or Cocktail)");
35       //      AliWarning(Form("%s %d: Unknown header type (not Pythia or Cocktail)",(char*)__FILE__,__LINE__));
36       return 0;
37     }
38     TList* headerList = genCocktailHeader->GetHeaders();
39     for (Int_t i=0; i<headerList->GetEntries(); i++) {
40       pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(headerList->At(i));
41       if (pythiaGenHeader)
42         break;
43     }
44     if(!pythiaGenHeader){
45       AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Pythia event header not found");
46       return 0;
47     }
48   }
49   return pythiaGenHeader;
50
51 }
52
53
54 void AliAnalysisHelperJetTasks::PrintStack(AliMCEvent *mcEvent,Int_t iFirst,Int_t iLast,Int_t iMaxPrint){
55
56   AliStack *stack = mcEvent->Stack();
57   if(!stack){
58     Printf("%s%d No Stack available",(char*)__FILE__,__LINE__);
59     return;
60   }
61
62   static Int_t iCount = 0;
63   if(iCount>iMaxPrint)return;
64   Int_t nStack = stack->GetNtrack();
65   if(iLast == 0)iLast = nStack;
66   else if(iLast > nStack)iLast = nStack;
67
68
69   Printf("####################################################################");
70   for(Int_t np = iFirst;np<iLast;++np){
71     TParticle *p = stack->Particle(np);
72     Printf("Nr.%d --- Status %d ---- Mother1 %d Mother2 %d Daughter1 %d Daughter2 %d ",
73            np,p->GetStatusCode(),p->GetMother(0),p->GetMother(1),p->GetDaughter(0),p->GetDaughter(1));
74     Printf("Eta %3.3f Phi %3.3f ",p->Eta(),p->Phi()); 
75     p->Print();    
76     Printf("---------------------------------------");
77   }
78   iCount++;
79 }
80
81
82
83
84 void AliAnalysisHelperJetTasks::GetClosestJets(AliAODJet *genJets,const Int_t &kGenJets,
85                                                AliAODJet *recJets,const Int_t &kRecJets,
86                                                Int_t *iGenIndex,Int_t *iRecIndex,
87                                                Int_t iDebug,Float_t maxDist){
88
89   //
90   // Relate the two input jet Arrays
91   //
92
93   //
94   // The association has to be unique
95   // So check in two directions
96   // find the closest rec to a gen
97   // and check if there is no other rec which is closer
98   // Caveat: Close low energy/split jets may disturb this correlation
99
100
101   // Idea: search in two directions generated e.g (a--e) and rec (1--3)
102   // Fill a matrix with Flags (1 for closest rec jet, 2 for closest rec jet
103   // in the end we have something like this
104   //    1   2   3
105   // ------------
106   // a| 3   2   0
107   // b| 0   1   0
108   // c| 0   0   3
109   // d| 0   0   1
110   // e| 0   0   1
111   // Topology
112   //   1     2
113   //     a         b        
114   //
115   //  d      c
116   //        3     e
117   // Only entries with "3" match from both sides
118
119   // In case we have more jets than kmaxjets only the 
120   // first kmaxjets are searched
121   // all other are -1
122   // use kMaxJets for a test not to fragemnt the memory...
123
124   for(int i = 0;i < kGenJets;++i)iGenIndex[i] = -1;
125   for(int j = 0;j < kRecJets;++j)iRecIndex[j] = -1;
126
127
128   
129   const int kMode = 3;
130
131   const Int_t nGenJets = TMath::Min(kMaxJets,kGenJets);
132   const Int_t nRecJets = TMath::Min(kMaxJets,kRecJets);
133
134   if(nRecJets==0||nGenJets==0)return;
135
136   // UShort_t *iFlag = new UShort_t[nGenJets*nRecJets];
137   UShort_t iFlag[kMaxJets*kMaxJets];
138   for(int i = 0;i < nGenJets;++i){
139     for(int j = 0;j < nRecJets;++j){
140       iFlag[i*nGenJets+j] = 0;
141     }
142   }
143
144
145
146   // find the closest distance to the generated
147   for(int ig = 0;ig<nGenJets;++ig){
148     Float_t dist = maxDist;
149     if(iDebug>1)Printf("Gen (%d) p_T %3.3f eta %3.3f ph %3.3f ",ig,genJets[ig].Pt(),genJets[ig].Eta(),genJets[ig].Phi());
150     for(int ir = 0;ir<nRecJets;++ir){
151       Double_t dR = genJets[ig].DeltaR(&recJets[ir]);
152       if(iDebug>1)Printf("Rec (%d) p_T %3.3f eta %3.3f ph %3.3f ",ir,recJets[ir].Pt(),recJets[ir].Eta(),recJets[ir].Phi());
153       if(iDebug>1)Printf("Distance (%d)--(%d) %3.3f ",ig,ir,dR);
154       if(dR<dist){
155         iRecIndex[ig] = ir;
156         dist = dR;
157       } 
158     }
159     if(iRecIndex[ig]>=0)iFlag[ig*nGenJets+iRecIndex[ig]]+=1;
160     // reset...
161     iRecIndex[ig] = -1;
162   }
163   // other way around
164   for(int ir = 0;ir<nRecJets;++ir){
165     Float_t dist = maxDist;
166     for(int ig = 0;ig<nGenJets;++ig){
167       Double_t dR = genJets[ig].DeltaR(&recJets[ir]);
168       if(dR<dist){
169         iGenIndex[ir] = ig;
170         dist = dR;
171       } 
172     }
173     if(iGenIndex[ir]>=0)iFlag[iGenIndex[ir]*nGenJets+ir]+=2;
174     // reset...
175     iGenIndex[ir] = -1;
176   }
177
178   // check for "true" correlations
179
180   if(iDebug>1)Printf(">>>>>> Matrix");
181
182   for(int ig = 0;ig<nGenJets;++ig){
183     for(int ir = 0;ir<nRecJets;++ir){
184       // Print
185       if(iDebug>1)printf("Flag[%d][%d] %d ",ig,ir,iFlag[ig*nGenJets+ir]);
186
187       if(kMode==3){
188         // we have a uniqie correlation
189         if(iFlag[ig*nGenJets+ir]==3){
190           iGenIndex[ir] = ig;
191           iRecIndex[ig] = ir;
192         }
193       }
194       else{
195         // we just take the correlation from on side
196         if((iFlag[ig*nGenJets+ir]&2)==2){
197           iGenIndex[ir] = ig;
198         }
199         if((iFlag[ig*nGenJets+ir]&1)==1){
200           iRecIndex[ig] = ir;
201         }
202       }
203     }
204     if(iDebug>1)printf("\n");
205   }
206 }
207
208
209
210 void  AliAnalysisHelperJetTasks::MergeOutput(char* cFiles, char* cList){
211
212   // This is used to merge the analysis-output from different 
213   // data samples/pt_hard bins
214   // in case the eventweigth was set to xsection/ntrials already, this
215   // is not needed. Both methods only work in case we do not mix different 
216   // pt_hard bins, and do not have overlapping bins
217
218   const Int_t nMaxBins = 12;
219   // LHC08q jetjet100: Mean = 1.42483e-03, RMS = 6.642e-05
220   // LHC08r jetjet50: Mean = 2.44068e-02, RMS = 1.144e-03
221   // LHC08v jetjet15-50: Mean = 2.168291 , RMS = 7.119e-02
222   // const Float_t xsection[nBins] = {2.168291,2.44068e-02};
223
224   Float_t xsection[nMaxBins];
225   Float_t nTrials[nMaxBins];
226   Float_t sf[nMaxBins];
227   TList *lIn[nMaxBins];
228   TFile *fIn[nMaxBins];
229
230   ifstream in1;
231   in1.open(cFiles);
232
233   char cFile[120];
234   Int_t ibTotal = 0;
235   while(in1>>cFile){
236     fIn[ibTotal] = TFile::Open(cFile);
237     lIn[ibTotal] = (TList*)fIn[ibTotal]->Get(cList);
238     Printf("Merging file %s",cFile);
239     if(!lIn[ibTotal]){
240       Printf("%s:%d No list %s found, exiting...",__FILE__,__LINE__,cList);
241       fIn[ibTotal]->ls();
242       return;
243     }
244     TH1* hTrials = (TH1F*)lIn[ibTotal]->FindObject("fh1Trials");
245     if(!hTrials){
246       Printf("%s:%d fh1PtHard_Trials not found in list, exiting...",__FILE__,__LINE__);
247       return;
248     }
249     TProfile* hXsec = (TProfile*)lIn[ibTotal]->FindObject("fh1Xsec");
250     if(!hXsec){
251       Printf("%s:%d fh1Xsec  not found in list, exiting...",__FILE__,__LINE__);
252       return;
253     }
254     xsection[ibTotal] = hXsec->GetBinContent(1);
255     nTrials[ibTotal] = hTrials->Integral();
256     sf[ibTotal] = xsection[ibTotal]/ nTrials[ibTotal];
257     ibTotal++;
258   }
259
260   if(ibTotal==0){
261     Printf("%s:%d No files found for mergin, exiting",__FILE__,__LINE__);
262     return;
263   }
264
265   TFile *fOut = new TFile("allpt.root","RECREATE");
266   TList *lOut = new TList();
267   lOut->SetName(lIn[0]->GetName());
268   // for the start scale all...
269   for(int ie = 0; ie < lIn[0]->GetEntries();++ie){
270     TH1 *h1Add = 0;
271     THnSparse *hnAdd = 0;
272     for(int ib = 0;ib < ibTotal;++ib){
273       // dynamic cast does not work with cint
274       TObject *h = lIn[ib]->At(ie);
275       if(h->InheritsFrom("TH1")){
276         TH1 *h1 = (TH1*)h;
277         if(ib==0){
278           h1Add = (TH1*)h1->Clone(h1->GetName());
279           h1Add->Scale(sf[ib]);
280         }
281         else{
282           h1Add->Add(h1,sf[ib]);
283         }
284       }
285       else if(h->InheritsFrom("THnSparse")){
286         THnSparse *hn = (THnSparse*)h;
287         if(ib==0){
288           hnAdd = (THnSparse*)hn->Clone(hn->GetName());
289           hnAdd->Scale(sf[ib]);
290         }
291         else{
292           hnAdd->Add(hn,sf[ib]);
293         }
294       }
295       
296
297     }// ib
298     if(h1Add)lOut->Add(h1Add);
299     else if(hnAdd)lOut->Add(hnAdd);
300   }
301   fOut->cd();
302   lOut->Write(lOut->GetName(),TObject::kSingleKey);
303   fOut->Close();
304 }