]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/JetTasks/AliAnalysisHelperJetTasks.cxx
6dea29c3856fc6f9b4ce3b77d570886c12e7df25
[u/mrichter/AliRoot.git] / PWG4 / JetTasks / AliAnalysisHelperJetTasks.cxx
1
2 #include "TROOT.h"
3 #include "TKey.h"
4 #include "TList.h"
5 #include "TSystem.h"
6 #include "TH1F.h"
7 #include "TProfile.h"
8 #include "THnSparse.h"
9 #include "TFile.h"
10 #include "TString.h"
11 #include "AliMCEvent.h"
12 #include "AliLog.h"
13 #include "AliAODJet.h"
14 #include "AliStack.h"
15 #include "AliGenEventHeader.h"
16 #include "AliGenCocktailEventHeader.h"
17 #include "AliGenPythiaEventHeader.h"
18 #include <fstream>
19 #include <iostream>
20 #include "AliAnalysisHelperJetTasks.h"
21 #include "TMatrixDSym.h"
22 #include "TMatrixDSymEigen.h"
23 #include "TVector.h"
24
25 ClassImp(AliAnalysisHelperJetTasks)
26
27
28
29  
30 AliGenPythiaEventHeader*  AliAnalysisHelperJetTasks::GetPythiaEventHeader(AliMCEvent *mcEvent){
31   
32   AliGenEventHeader* genHeader = mcEvent->GenEventHeader();
33   AliGenPythiaEventHeader* pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(genHeader);
34   if(!pythiaGenHeader){
35     // cocktail ??
36     AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast<AliGenCocktailEventHeader*>(genHeader);
37     
38     if (!genCocktailHeader) {
39       AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Unknown header type (not Pythia or Cocktail)");
40       //      AliWarning(Form("%s %d: Unknown header type (not Pythia or Cocktail)",(char*)__FILE__,__LINE__));
41       return 0;
42     }
43     TList* headerList = genCocktailHeader->GetHeaders();
44     for (Int_t i=0; i<headerList->GetEntries(); i++) {
45       pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(headerList->At(i));
46       if (pythiaGenHeader)
47         break;
48     }
49     if(!pythiaGenHeader){
50       AliWarningGeneral(Form(" %s:%d",(char*)__FILE__,__LINE__),"Pythia event header not found");
51       return 0;
52     }
53   }
54   return pythiaGenHeader;
55
56 }
57
58
59 void AliAnalysisHelperJetTasks::PrintStack(AliMCEvent *mcEvent,Int_t iFirst,Int_t iLast,Int_t iMaxPrint){
60
61   AliStack *stack = mcEvent->Stack();
62   if(!stack){
63     Printf("%s%d No Stack available",(char*)__FILE__,__LINE__);
64     return;
65   }
66
67   static Int_t iCount = 0;
68   if(iCount>iMaxPrint)return;
69   Int_t nStack = stack->GetNtrack();
70   if(iLast == 0)iLast = nStack;
71   else if(iLast > nStack)iLast = nStack;
72
73
74   Printf("####################################################################");
75   for(Int_t np = iFirst;np<iLast;++np){
76     TParticle *p = stack->Particle(np);
77     Printf("Nr.%d --- Status %d ---- Mother1 %d Mother2 %d Daughter1 %d Daughter2 %d ",
78            np,p->GetStatusCode(),p->GetMother(0),p->GetMother(1),p->GetDaughter(0),p->GetDaughter(1));
79     Printf("Eta %3.3f Phi %3.3f ",p->Eta(),p->Phi()); 
80     p->Print();    
81     Printf("---------------------------------------");
82   }
83   iCount++;
84 }
85
86
87
88
89 void AliAnalysisHelperJetTasks::GetClosestJets(AliAODJet *genJets,const Int_t &kGenJets,
90                                                AliAODJet *recJets,const Int_t &kRecJets,
91                                                Int_t *iGenIndex,Int_t *iRecIndex,
92                                                Int_t iDebug,Float_t maxDist){
93
94   //
95   // Relate the two input jet Arrays
96   //
97
98   //
99   // The association has to be unique
100   // So check in two directions
101   // find the closest rec to a gen
102   // and check if there is no other rec which is closer
103   // Caveat: Close low energy/split jets may disturb this correlation
104
105
106   // Idea: search in two directions generated e.g (a--e) and rec (1--3)
107   // Fill a matrix with Flags (1 for closest rec jet, 2 for closest rec jet
108   // in the end we have something like this
109   //    1   2   3
110   // ------------
111   // a| 3   2   0
112   // b| 0   1   0
113   // c| 0   0   3
114   // d| 0   0   1
115   // e| 0   0   1
116   // Topology
117   //   1     2
118   //     a         b        
119   //
120   //  d      c
121   //        3     e
122   // Only entries with "3" match from both sides
123
124   // In case we have more jets than kmaxjets only the 
125   // first kmaxjets are searched
126   // all other are -1
127   // use kMaxJets for a test not to fragemnt the memory...
128
129   for(int i = 0;i < kGenJets;++i)iGenIndex[i] = -1;
130   for(int j = 0;j < kRecJets;++j)iRecIndex[j] = -1;
131
132
133   
134   const int kMode = 3;
135
136   const Int_t nGenJets = TMath::Min(kMaxJets,kGenJets);
137   const Int_t nRecJets = TMath::Min(kMaxJets,kRecJets);
138
139   if(nRecJets==0||nGenJets==0)return;
140
141   // UShort_t *iFlag = new UShort_t[nGenJets*nRecJets];
142   UShort_t iFlag[kMaxJets*kMaxJets];
143   for(int i = 0;i < nGenJets;++i){
144     for(int j = 0;j < nRecJets;++j){
145       iFlag[i*nGenJets+j] = 0;
146     }
147   }
148
149
150
151   // find the closest distance to the generated
152   for(int ig = 0;ig<nGenJets;++ig){
153     Float_t dist = maxDist;
154     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());
155     for(int ir = 0;ir<nRecJets;++ir){
156       Double_t dR = genJets[ig].DeltaR(&recJets[ir]);
157       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());
158       if(iDebug>1)Printf("Distance (%d)--(%d) %3.3f ",ig,ir,dR);
159       if(dR<dist){
160         iRecIndex[ig] = ir;
161         dist = dR;
162       } 
163     }
164     if(iRecIndex[ig]>=0)iFlag[ig*nGenJets+iRecIndex[ig]]+=1;
165     // reset...
166     iRecIndex[ig] = -1;
167   }
168   // other way around
169   for(int ir = 0;ir<nRecJets;++ir){
170     Float_t dist = maxDist;
171     for(int ig = 0;ig<nGenJets;++ig){
172       Double_t dR = genJets[ig].DeltaR(&recJets[ir]);
173       if(dR<dist){
174         iGenIndex[ir] = ig;
175         dist = dR;
176       } 
177     }
178     if(iGenIndex[ir]>=0)iFlag[iGenIndex[ir]*nGenJets+ir]+=2;
179     // reset...
180     iGenIndex[ir] = -1;
181   }
182
183   // check for "true" correlations
184
185   if(iDebug>1)Printf(">>>>>> Matrix");
186
187   for(int ig = 0;ig<nGenJets;++ig){
188     for(int ir = 0;ir<nRecJets;++ir){
189       // Print
190       if(iDebug>1)printf("Flag[%d][%d] %d ",ig,ir,iFlag[ig*nGenJets+ir]);
191
192       if(kMode==3){
193         // we have a uniqie correlation
194         if(iFlag[ig*nGenJets+ir]==3){
195           iGenIndex[ir] = ig;
196           iRecIndex[ig] = ir;
197         }
198       }
199       else{
200         // we just take the correlation from on side
201         if((iFlag[ig*nGenJets+ir]&2)==2){
202           iGenIndex[ir] = ig;
203         }
204         if((iFlag[ig*nGenJets+ir]&1)==1){
205           iRecIndex[ig] = ir;
206         }
207       }
208     }
209     if(iDebug>1)printf("\n");
210   }
211 }
212
213
214
215 void  AliAnalysisHelperJetTasks::MergeOutput(char* cFiles, char* cList){
216
217   // This is used to merge the analysis-output from different 
218   // data samples/pt_hard bins
219   // in case the eventweigth was set to xsection/ntrials already, this
220   // is not needed. Both methods only work in case we do not mix different 
221   // pt_hard bins, and do not have overlapping bins
222
223   const Int_t nMaxBins = 12;
224   // LHC08q jetjet100: Mean = 1.42483e-03, RMS = 6.642e-05
225   // LHC08r jetjet50: Mean = 2.44068e-02, RMS = 1.144e-03
226   // LHC08v jetjet15-50: Mean = 2.168291 , RMS = 7.119e-02
227   // const Float_t xsection[nBins] = {2.168291,2.44068e-02};
228
229   Float_t xsection[nMaxBins];
230   Float_t nTrials[nMaxBins];
231   Float_t sf[nMaxBins];
232   TList *lIn[nMaxBins];
233   TFile *fIn[nMaxBins];
234
235   ifstream in1;
236   in1.open(cFiles);
237
238   char cFile[120];
239   Int_t ibTotal = 0;
240   while(in1>>cFile){
241     fIn[ibTotal] = TFile::Open(cFile);
242     lIn[ibTotal] = (TList*)fIn[ibTotal]->Get(cList);
243     Printf("Merging file %s",cFile);
244     if(!lIn[ibTotal]){
245       Printf("%s:%d No list %s found, exiting...",__FILE__,__LINE__,cList);
246       fIn[ibTotal]->ls();
247       return;
248     }
249     TH1* hTrials = (TH1F*)lIn[ibTotal]->FindObject("fh1Trials");
250     if(!hTrials){
251       Printf("%s:%d fh1PtHard_Trials not found in list, exiting...",__FILE__,__LINE__);
252       return;
253     }
254     TProfile* hXsec = (TProfile*)lIn[ibTotal]->FindObject("fh1Xsec");
255     if(!hXsec){
256       Printf("%s:%d fh1Xsec  not found in list, exiting...",__FILE__,__LINE__);
257       return;
258     }
259     xsection[ibTotal] = hXsec->GetBinContent(1);
260     nTrials[ibTotal] = hTrials->Integral();
261     sf[ibTotal] = xsection[ibTotal]/ nTrials[ibTotal];
262     ibTotal++;
263   }
264
265   if(ibTotal==0){
266     Printf("%s:%d No files found for mergin, exiting",__FILE__,__LINE__);
267     return;
268   }
269
270   TFile *fOut = new TFile("allpt.root","RECREATE");
271   TList *lOut = new TList();
272   lOut->SetName(lIn[0]->GetName());
273   // for the start scale all...
274   for(int ie = 0; ie < lIn[0]->GetEntries();++ie){
275     TH1 *h1Add = 0;
276     THnSparse *hnAdd = 0;
277     for(int ib = 0;ib < ibTotal;++ib){
278       // dynamic cast does not work with cint
279       TObject *h = lIn[ib]->At(ie);
280       if(h->InheritsFrom("TH1")){
281         TH1 *h1 = (TH1*)h;
282         if(ib==0){
283           h1Add = (TH1*)h1->Clone(h1->GetName());
284           h1Add->Scale(sf[ib]);
285         }
286         else{
287           h1Add->Add(h1,sf[ib]);
288         }
289       }
290       else if(h->InheritsFrom("THnSparse")){
291         THnSparse *hn = (THnSparse*)h;
292         if(ib==0){
293           hnAdd = (THnSparse*)hn->Clone(hn->GetName());
294           hnAdd->Scale(sf[ib]);
295         }
296         else{
297           hnAdd->Add(hn,sf[ib]);
298         }
299       }
300       
301
302     }// ib
303     if(h1Add)lOut->Add(h1Add);
304     else if(hnAdd)lOut->Add(hnAdd);
305   }
306   fOut->cd();
307   lOut->Write(lOut->GetName(),TObject::kSingleKey);
308   fOut->Close();
309 }
310
311 Bool_t AliAnalysisHelperJetTasks::PythiaInfoFromFile(const char* currFile,Float_t &fXsec,Float_t &fTrials){
312   //
313   // get the cross section and the trails either from pyxsec.root or from pysec_hists.root
314   // This is to called in Notify and should provide the path to the AOD/ESD file
315
316   TString file(currFile);  
317   fXsec = 0;
318   fTrials = 1;
319
320   if(file.Contains("root_archive.zip#")){
321     Ssiz_t pos1 = file.Index("root_archive",12,TString::kExact);
322     Ssiz_t pos = file.Index("#",1,pos1,TString::kExact);
323     file.Replace(pos+1,20,"");
324   }
325   else {
326     // not an archive take the basename....
327     file.ReplaceAll(gSystem->BaseName(file.Data()),"");
328   }
329   Printf("%s",file.Data());
330   
331  
332    
333
334   TFile *fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec.root")); // problem that we cannot really test the existance of a file in a archive so we have to lvie with open error message from root
335   if(!fxsec){
336     // next trial fetch the histgram file
337     fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec_hists.root"));
338     if(!fxsec){
339         // not a severe condition but inciate that we have no information
340       return kFALSE;
341     }
342     else{
343       // find the tlist we want to be independtent of the name so use the Tkey
344       TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0); 
345       if(!key){
346         fxsec->Close();
347         return kFALSE;
348       }
349       TList *list = dynamic_cast<TList*>(key->ReadObj());
350       if(!list){
351         fxsec->Close();
352         return kFALSE;
353       }
354       fXsec = ((TProfile*)list->FindObject("h1Xsec"))->GetBinContent(1);
355       fTrials  = ((TH1F*)list->FindObject("h1Trials"))->GetBinContent(1);
356       fxsec->Close();
357     }
358   } // no tree pyxsec.root
359   else {
360     TTree *xtree = (TTree*)fxsec->Get("Xsection");
361     if(!xtree){
362       fxsec->Close();
363       return kFALSE;
364     }
365     UInt_t   ntrials  = 0;
366     Double_t  xsection  = 0;
367     xtree->SetBranchAddress("xsection",&xsection);
368     xtree->SetBranchAddress("ntrials",&ntrials);
369     xtree->GetEntry(0);
370     fTrials = ntrials;
371     fXsec = xsection;
372     fxsec->Close();
373   }
374   return kTRUE;
375 }
376
377
378 Bool_t AliAnalysisHelperJetTasks::GetThrustAxis(TVector3 &n01, TVector3 * pTrack, const Int_t &nTracks)
379 {
380   //
381   // fetch the thrust axis
382   // sona.pochybova@cern.ch
383
384   const Int_t kTracks = 1000;
385   if(nTracks>kTracks)return kFALSE;
386
387   TVector3 psum;
388   Double_t psum1 = 0;
389   Double_t psum2 = 0;
390   Double_t thrust[kTracks];
391   Double_t th = -3;
392   Double_t tpom = -1;
393   Int_t j = 0;
394
395   
396   //  for(Int_t j = 0; j < nTracks; j++)
397   while(TMath::Abs(th-tpom) > 10e-7 && j < nTracks)
398     {
399       th = tpom;
400       psum.SetXYZ(0., 0., 0.);
401       psum1 = 0;
402       psum2 = 0;
403       for(Int_t i = 0; i < nTracks; i++)
404         {
405           psum1 += (TMath::Abs(n01.Dot(pTrack[i])));
406           psum2 += pTrack[i].Mag();
407           
408           if (n01.Dot(pTrack[i]) > 0) psum += pTrack[i];
409           if (n01.Dot(pTrack[i]) < 0) psum -= pTrack[i];
410         }
411       
412       thrust[j] = psum1/psum2;
413       
414       //      if(th == thrust[j]) 
415       //        break;
416       
417       tpom = thrust[j];
418       
419       n01 = psum.Unit();
420       j++;
421     }
422   return kTRUE;
423 }
424
425 //___________________________________________________________________________________________________________
426
427 Bool_t AliAnalysisHelperJetTasks::GetEventShapes(TVector3 &n01, TVector3 * pTrack, Int_t nTracks, Double_t * eventShapes)
428 {       
429   // ***
430   // Event shape calculation
431   // sona.pochybova@cern.ch
432
433   const Int_t kTracks = 1000;
434   if(nTracks>kTracks)return kFALSE;
435
436   TVector3 psum;
437   TVector3 pTrackPerp[kTracks];
438   Double_t psum1 = 0;
439   Double_t psum2 = 0;
440   Double_t thrust[kTracks];
441   Double_t th = -3;
442   Double_t tpom = -1;
443
444   //Sphericity calculation
445   TMatrixDSym m(3);
446   Double_t s00 = 0;
447   Double_t s01 = 0;
448   Double_t s02 = 0;
449   
450   Double_t s10 = 0;
451   Double_t s11 = 0;
452   Double_t s12 = 0;
453   
454   Double_t s20 = 0;
455   Double_t s21 = 0;
456   Double_t s22 = 0;
457   
458   Double_t ptot = 0;
459   
460   Double_t c = -10;
461
462 //
463 //loop for thrust calculation  
464 //
465   Int_t k = 0;
466   while(TMath::Abs(th-tpom) > 10e-7 && k < nTracks)
467     {  
468       th = tpom;
469       psum.SetXYZ(0., 0., 0.);
470       psum1 = 0;
471       psum2 = 0;
472       for(Int_t i = 0; i < nTracks; i++)
473         {
474           pTrackPerp[i].SetXYZ(pTrack[i].X(), pTrack[i].Y(), 0);
475
476           psum1 += (TMath::Abs(n01.Dot(pTrackPerp[i])));
477           psum2 += pTrackPerp[i].Mag();
478
479           if (n01.Dot(pTrackPerp[i]) > 0) psum += pTrackPerp[i];
480           if (n01.Dot(pTrackPerp[i]) < 0) psum -= pTrackPerp[i];
481         }
482       
483       thrust[k] = psum1/psum2;
484       
485       tpom = thrust[k];
486       n01 = psum.Unit();
487       k++;
488 //      Printf("========== %d +++ th :: %f, tpom :: %f=============\n", k, th, tpom);
489     }
490
491   eventShapes[0] = th;
492
493 //
494 //other event shapes variables
495 //
496   for(Int_t j = 0; j < nTracks; j++)
497     {  
498       s00 = s00 + (pTrack[j].Px()*pTrack[j].Px())/pTrack[j].Mag();
499       s01 = s01 + (pTrack[j].Px()*pTrack[j].Py())/pTrack[j].Mag();
500       s02 = s02 + (pTrack[j].Px()*pTrack[j].Pz())/pTrack[j].Mag();
501       
502       s10 = s10 + (pTrack[j].Py()*pTrack[j].Px())/pTrack[j].Mag();
503       s11 = s11 + (pTrack[j].Py()*pTrack[j].Py())/pTrack[j].Mag();
504       s12 = s12 + (pTrack[j].Py()*pTrack[j].Pz())/pTrack[j].Mag();
505       
506       s20 = s20 + (pTrack[j].Pz()*pTrack[j].Px())/pTrack[j].Mag();
507       s21 = s21 + (pTrack[j].Pz()*pTrack[j].Py())/pTrack[j].Mag();
508       s22 = s22 + (pTrack[j].Pz()*pTrack[j].Pz())/pTrack[j].Mag();
509       
510       ptot += pTrack[j].Mag();
511     }
512
513   if(ptot > 0.)
514     {
515       m(0,0) = s00/ptot;
516       m(0,1) = s01/ptot;
517       m(0,2) = s02/ptot;
518
519       m(1,0) = s10/ptot;
520       m(1,1) = s11/ptot;
521       m(1,2) = s12/ptot;
522
523       m(2,0) = s20/ptot;
524       m(2,1) = s21/ptot;
525       m(2,2) = s22/ptot;
526
527       TMatrixDSymEigen eigen(m);
528       TVectorD eigenVal = eigen.GetEigenValues();
529
530       Double_t sphericity = (3/2)*(eigenVal(2)+eigenVal(1));
531       eventShapes[1] = sphericity;
532
533       Double_t aplanarity = (3/2)*(eigenVal(2));
534       eventShapes[2] = aplanarity;
535
536       c = 3*(eigenVal(0)*eigenVal(1)+eigenVal(0)*eigenVal(2)+eigenVal(1)*eigenVal(2));
537       eventShapes[3] = c;
538     }
539   return kTRUE;
540 }
541   
542
543
544  //__________________________________________________________________________________________________________________________