]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/SPECTRA/PiKaPr/TestAOD/AliAnalysisTaskSpectraBoth.cxx
Merge branch 'TPCdev' of https://git.cern.ch/reps/AliRoot into TPCdev
[u/mrichter/AliRoot.git] / PWGLF / SPECTRA / PiKaPr / TestAOD / AliAnalysisTaskSpectraBoth.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2009, 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 //-----------------------------------------------------------------
17 //         AliAnalysisTaskSpectraBoth class
18 //-----------------------------------------------------------------
19
20 #include "TChain.h"
21 #include "TTree.h"
22 #include "TLegend.h"
23 #include "TH1F.h"
24 #include "TH2F.h"
25 #include "TH3F.h"
26 #include "TCanvas.h"
27 #include "AliAnalysisTask.h"
28 #include "AliAnalysisManager.h"
29 #include "AliAODTrack.h"
30 #include "AliAODMCParticle.h"
31 #include "AliVParticle.h"
32 #include "AliAODEvent.h"
33 #include "AliAODInputHandler.h"
34 #include "AliAnalysisTaskSpectraBoth.h"
35 #include "AliAnalysisTaskESDfilter.h"
36 #include "AliAnalysisDataContainer.h"
37 #include "AliSpectraBothHistoManager.h"
38 #include "AliSpectraBothTrackCuts.h"
39 #include "AliSpectraBothEventCuts.h"
40 #include "AliCentrality.h"
41 #include "TProof.h"
42 #include "AliPID.h"
43 #include "AliVEvent.h"
44 #include "AliESDEvent.h"
45 #include "AliPIDResponse.h"
46 #include "AliStack.h"
47 #include "AliSpectraBothPID.h"
48 #include "AliGenEventHeader.h"  
49 #include <TMCProcess.h>
50
51 #include <iostream>
52
53
54
55
56 using namespace AliSpectraNameSpaceBoth;
57 using namespace std;
58
59 ClassImp(AliAnalysisTaskSpectraBoth)
60
61 //________________________________________________________________________
62 AliAnalysisTaskSpectraBoth::AliAnalysisTaskSpectraBoth(const char *name) : AliAnalysisTaskSE(name), fAOD(0), fHistMan(0), fTrackCuts(0), fEventCuts(0),  fPID(0), fIsMC(0), fNRebin(0),fUseMinSigma(0),fCuts(0),fdotheMCLoopAfterEventCuts(0)
63 {
64   // Default constructor
65   
66   DefineInput(0, TChain::Class());
67   DefineOutput(1, AliSpectraBothHistoManager::Class());
68   DefineOutput(2, AliSpectraBothEventCuts::Class());
69   DefineOutput(3, AliSpectraBothTrackCuts::Class());
70   DefineOutput(4, AliSpectraBothPID::Class());
71   fNRebin=0;
72   
73 }
74 //________________________________________________________________________
75 //________________________________________________________________________
76 void AliAnalysisTaskSpectraBoth::UserCreateOutputObjects()
77 {
78   // create output objects
79   fHistMan = new AliSpectraBothHistoManager("SpectraHistos",fNRebin);
80
81   if (!fTrackCuts) AliFatal("Track Cuts should be set in the steering macro");
82   if (!fEventCuts) AliFatal("Event Cuts should be set in the steering macro");
83   if (!fPID)       AliFatal("PID object should be set in the steering macro");
84   fTrackCuts->SetAliESDtrackCuts(fCuts);
85   fEventCuts->InitHisto();
86   fTrackCuts->InitHisto();
87
88   PostData(1, fHistMan  );
89   PostData(2, fEventCuts);
90   PostData(3, fTrackCuts);
91   PostData(4, fPID      );
92 }
93 //________________________________________________________________________
94 void AliAnalysisTaskSpectraBoth::UserExec(Option_t *)
95 {
96   // main event loop
97         Int_t ifAODEvent=AliSpectraBothTrackCuts::kotherobject;
98         fAOD = dynamic_cast<AliVEvent*>(InputEvent());
99   //    AliESDEvent* esdevent=0x0;
100  //     AliAODEvent* aodevent=0x0;
101    
102         TString nameoftrack(fAOD->ClassName());  
103         if(!nameoftrack.CompareTo("AliESDEvent"))
104         {
105                 ifAODEvent=AliSpectraBothTrackCuts::kESDobject;
106                 //esdevent=dynamic_cast<AliESDEvent*>(fAOD);
107         }
108         else if(!nameoftrack.CompareTo("AliAODEvent"))
109         {
110                 ifAODEvent=AliSpectraBothTrackCuts::kAODobject;
111                 //aodevent=dynamic_cast<AliAODEvent*>(fAOD);
112         }
113         else
114                 AliFatal("Not processing AODs or ESDS") ;
115         if(fdotheMCLoopAfterEventCuts)
116                 if(!fEventCuts->IsSelected(fAOD,fTrackCuts))return;//event selection
117         TClonesArray *arrayMC = 0;
118         Int_t npar=0;
119         AliStack* stack=0x0;
120          Double_t mcZ=-100;
121
122         if (fIsMC)
123         {
124                 TArrayF mcVertex(3);
125                 mcVertex[0]=9999.; mcVertex[1]=9999.; mcVertex[2]=9999.;
126                 AliMCEvent* mcEvent=(AliMCEvent*)MCEvent();
127                 if (!mcEvent) 
128                 {
129                         AliFatal("Error: MC particles branch not found!\n");
130                 }
131                 AliHeader* header = mcEvent->Header();
132                 if (!header) 
133                 {
134                         AliDebug(AliLog::kError, "Header not available");
135                         return;
136                 }
137         
138                 AliGenEventHeader* genHeader = header->GenEventHeader();
139                 if(genHeader)
140                 {
141                         genHeader->PrimaryVertex(mcVertex);
142                         mcZ=mcVertex[2];
143                 }
144                   if(ifAODEvent==AliSpectraBothTrackCuts::kAODobject)
145                   {
146                           arrayMC = (TClonesArray*) fAOD->GetList()->FindObject(AliAODMCParticle::StdBranchName());
147                           if (!arrayMC) 
148                           {
149                                         AliFatal("Error: MC particles branch not found!\n");
150                           }
151                           Int_t nMC = arrayMC->GetEntries();
152                           for (Int_t iMC = 0; iMC < nMC; iMC++)
153                           {
154                                   AliAODMCParticle *partMC = (AliAODMCParticle*) arrayMC->At(iMC);
155                                   if(!partMC->Charge()) continue;//Skip neutrals
156                                   //if(partMC->Eta() > fTrackCuts->GetEtaMin() && partMC->Eta() < fTrackCuts->GetEtaMax()){//charged hadron are filled inside the eta acceptance
157                                   //Printf("%f     %f-%f",partMC->Eta(),fTrackCuts->GetEtaMin(),fTrackCuts->GetEtaMax());
158                                   if(partMC->Eta() > fTrackCuts->GetEtaMin() && partMC->Eta() < fTrackCuts->GetEtaMax())
159                                                 fHistMan->GetPtHistogram(kHistPtGen)->Fill(partMC->Pt(),partMC->IsPhysicalPrimary());                                                                    
160                                   //rapidity cut
161                                   if(partMC->Y() > fTrackCuts->GetYMax()|| partMC->Y() < fTrackCuts->GetYMin() ) 
162                                         continue;       
163                                   if(partMC->IsPhysicalPrimary())
164                                          npar++;    
165                                   // check for true PID + and fill P_t histos 
166                                   Int_t charge = partMC->Charge() > 0 ? kChPos : kChNeg ;
167                                   Int_t id = fPID->GetParticleSpecie(partMC);
168                                   if(id != kSpUndefined) 
169                                   {
170                                         fHistMan->GetHistogram2D(kHistPtGenTruePrimary,id,charge)->Fill(partMC->Pt(),partMC->IsPhysicalPrimary());
171                                   }
172                           }
173                   }
174                   if(ifAODEvent==AliSpectraBothTrackCuts::kESDobject)
175                   {
176                         stack = mcEvent->Stack();
177                         Int_t nMC = stack->GetNtrack();
178                         for (Int_t iMC = 0; iMC < nMC; iMC++)
179                         {
180                                  
181                                 TParticle *partMC = stack->Particle(iMC);
182                                 
183                                 if(!partMC)     
184                                         continue;       
185         
186                                 if(!partMC->GetPDG(0))
187                                         continue;
188                                 if(TMath::Abs(partMC->GetPDG(0)->Charge()/3.0)<0.01) 
189                                         continue;//Skip neutrals
190                                 if(partMC->Eta() > fTrackCuts->GetEtaMin() && partMC->Eta() < fTrackCuts->GetEtaMax())
191                                         fHistMan->GetPtHistogram(kHistPtGen)->Fill(partMC->Pt(),stack->IsPhysicalPrimary(iMC));
192                                 if(partMC->Y()   > fTrackCuts->GetYMax() ||partMC->Y()   < fTrackCuts->GetYMin()  ) 
193                                         continue;
194                                 if(stack->IsPhysicalPrimary(iMC))
195                                          npar++;    
196                                   // check for true PID + and fill P_t histos 
197                                 Int_t charge = partMC->GetPDG(0)->Charge()/3.0 > 0 ? kChPos : kChNeg ;
198                                 Int_t id = fPID->GetParticleSpecie(partMC);
199                                 if(id != kSpUndefined) 
200                                 {
201                                         fHistMan->GetHistogram2D(kHistPtGenTruePrimary,id,charge)->Fill(partMC->Pt(),stack->IsPhysicalPrimary(iMC));
202                                 }
203                           }
204                   }
205         }
206         if(!fdotheMCLoopAfterEventCuts)
207                 if(!fEventCuts->IsSelected(fAOD,fTrackCuts,fIsMC,mcZ))return;//event selection
208         //main loop on tracks
209         Int_t ntracks=0;
210         //cout<<fAOD->GetNumberOfTracks()<<endl;
211         for (Int_t iTracks = 0; iTracks < fAOD->GetNumberOfTracks(); iTracks++) 
212         {
213                 AliVTrack* track = dynamic_cast<AliVTrack*>(fAOD->GetTrack(iTracks));
214                 AliAODTrack* aodtrack=0;
215                 AliESDtrack* esdtrack=0;
216                 Float_t dca=-999.;
217                 Float_t dcaz=-999.;
218                 Short_t ncls=-1; 
219                 Float_t chi2perndf=-1.0;
220                 if(ifAODEvent==AliSpectraBothTrackCuts::kESDobject)
221                 {
222                         esdtrack=dynamic_cast<AliESDtrack*>(track);
223                         if(!esdtrack)
224                                 continue;
225                         esdtrack->GetImpactParameters(dca,dcaz);
226                         ncls=esdtrack->GetTPCNcls();
227                         if ( ncls > 5) 
228                                 chi2perndf=(esdtrack->GetTPCchi2()/Float_t(ncls - 5));
229   
230                         else 
231                                 chi2perndf=-1;
232                         
233                 }
234                 else if (ifAODEvent==AliSpectraBothTrackCuts::kAODobject)
235                 {
236                         aodtrack=dynamic_cast<AliAODTrack*>(track);
237                         if(!aodtrack)
238                                 continue;               
239                         dca=aodtrack->DCA();
240                         dcaz=aodtrack->ZAtDCA();
241                         ncls=aodtrack->GetTPCNcls();
242                         chi2perndf=aodtrack->Chi2perNDF();
243                 }
244                 else
245                         continue;
246                 if (!fTrackCuts->IsSelected(track,kTRUE)) 
247                         continue;
248                 ntracks++;
249                 fPID->FillQAHistos(fHistMan, track, fTrackCuts);
250                 
251                 //calculate DCA for AOD track
252                 if(dca==-999.)
253                 {// track->DCA() does not work in old AOD production
254                         Double_t d[2], covd[3];
255                         AliVTrack* track_clone=(AliVTrack*)track->Clone("track_clone"); // need to clone because PropagateToDCA updates the track parameters
256                         Bool_t isDCA = track_clone->PropagateToDCA(fAOD->GetPrimaryVertex(),fAOD->GetMagneticField(),9999.,d,covd);
257                         delete track_clone;
258                         if(!isDCA)
259                                 d[0]=-999.;
260                         dca=d[0];
261                         dcaz=d[1];      
262                 }
263         
264                 fHistMan->GetPtHistogram(kHistPtRec)->Fill(track->Pt(),dca);  // PT histo
265                 // get identity and charge
266                 Bool_t rec[3]={false,false,false};
267                 Int_t idRec  = fPID->GetParticleSpecie(fHistMan,track, fTrackCuts,rec);
268                 for(int irec=kSpPion;irec<kNSpecies;irec++)
269                 {
270    
271                         if(fUseMinSigma)
272                         {
273                                 if(irec>kSpPion)
274                                         break;
275                         }
276                         else
277                         {       
278                                 if(!rec[irec]) 
279                                         idRec = kSpUndefined;
280                                 else    
281                                         idRec=irec;
282                         }               
283    
284                         Int_t charge = track->Charge() > 0 ? kChPos : kChNeg;
285                 
286                         // Fill histograms, only if inside y and nsigma acceptance
287                         if(idRec != kSpUndefined && fTrackCuts->CheckYCut ((BothParticleSpecies_t)idRec))
288                         {
289                                 fHistMan->GetHistogram2D(kHistPtRecSigma,idRec,charge)->Fill(track->Pt(),dca);
290                                 fTrackCuts->GetHistoDCAzQA()->Fill(idRec,track->Pt(),dcaz);
291                                 fTrackCuts->GetHistoNclustersQA()->Fill(idRec,track->Pt(),ncls);
292                                 fTrackCuts->GetHistochi2perNDFQA()->Fill(idRec,track->Pt(),chi2perndf);
293                         }
294                         //can't put a continue because we still have to fill allcharged primaries, done later
295                 
296                         /* MC Part */
297                         if (arrayMC||stack) 
298                         {
299                                 Bool_t isPrimary           = kFALSE;
300                                 Bool_t isSecondaryMaterial = kFALSE; 
301                                 Bool_t isSecondaryWeak     = kFALSE; 
302                                 Int_t idGen     =kSpUndefined;
303                                 Int_t pdgcode=0;
304                                 if (ifAODEvent==AliSpectraBothTrackCuts::kAODobject)
305                                 {
306                                         AliAODMCParticle *partMC = (AliAODMCParticle*) arrayMC->At(TMath::Abs(track->GetLabel()));
307                                         if (!partMC) 
308                                         { 
309                                                 AliError("Cannot get MC particle");
310                                                 continue; 
311                                         }
312                                         // Check if it is primary, secondary from material or secondary from weak decay
313                                         isPrimary           = partMC->IsPhysicalPrimary();
314                                         isSecondaryWeak     = partMC->IsSecondaryFromWeakDecay();
315                                         isSecondaryMaterial      = partMC->IsSecondaryFromMaterial();
316                                         //cout<<"AOD tagging "<<isPrimary<<" "<<isSecondaryWeak<<isSecondaryMaterial<<" "<<partMC->GetMCProcessCode()<<endl;
317
318                                         if(!isPrimary&&!isSecondaryWeak&&!isSecondaryMaterial) 
319                                         {
320                                                 AliError("old tagging");
321                                                 Int_t mfl=-999,codemoth=-999;
322                                                 Int_t indexMoth=partMC->GetMother(); // FIXME ignore fakes? TO BE CHECKED, on ESD is GetFirstMother()
323                                                 if(indexMoth>=0)
324                                                 {//is not fake
325                                                         AliAODMCParticle* moth = (AliAODMCParticle*) arrayMC->At(indexMoth);
326                                                         codemoth = TMath::Abs(moth->GetPdgCode());
327                                                         mfl = Int_t (codemoth/ TMath::Power(10, Int_t(TMath::Log10(codemoth))));
328                                                 }
329                                                 //Int_t uniqueID = partMC->GetUniqueID();
330                                                 //cout<<"uniqueID: "<<partMC->GetUniqueID()<<"       "<<kPDecay<<endl;
331                                                 //cout<<"status: "<<partMC->GetStatus()<<"       "<<kPDecay<<endl;
332                                                 // if(uniqueID == kPDecay)Printf("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
333                                                 if(mfl==3) 
334                                                         isSecondaryWeak     = kTRUE; // add if(partMC->GetStatus() & kPDecay)? FIXME
335                                                 else       
336                                                         isSecondaryMaterial = kTRUE;
337                                         }
338                                         //cout<<"AOD 2 tagging "<<isPrimary<<" "<<isSecondaryWeak<<isSecondaryMaterial<<" "<<partMC->GetMCProcessCode()<<endl;
339
340
341                                         idGen     = fPID->GetParticleSpecie(partMC);
342                                         pdgcode=partMC->GetPdgCode(); 
343                                 }
344                                 else if (ifAODEvent==AliSpectraBothTrackCuts::kESDobject)
345                                 {
346                                         TParticle *partMC =stack->Particle(TMath::Abs(track->GetLabel()));
347                                         if (!partMC) 
348                                         { 
349                                                 AliError("Cannot get MC particle");
350                                                 continue; 
351                                         }
352                                         isPrimary           = stack->IsPhysicalPrimary(TMath::Abs(track->GetLabel()));
353                                         isSecondaryWeak     = stack->IsSecondaryFromWeakDecay(TMath::Abs(track->GetLabel()));
354                                         isSecondaryMaterial      = stack->IsSecondaryFromMaterial(TMath::Abs(track->GetLabel()));
355                                         //cout<<"ESD tagging "<<isPrimary<<" "<<isSecondaryWeak<<isSecondaryMaterial<<endl;
356                                         
357                                         idGen     = fPID->GetParticleSpecie(partMC);
358                                         pdgcode=partMC->GetPdgCode(); 
359                                 }
360                                 else
361                                         return;
362                         
363                         //        cout<<isPrimary<<" "<<isSecondaryWeak<<" "<<isSecondaryMaterial<<endl;
364                         //        cout<<" functions "<<partMC->IsPhysicalPrimary()<<" "<<partMC->IsSecondaryFromWeakDecay()<<" "<<partMC->IsSecondaryFromMaterial()<<endl;
365                   
366                                 if (isPrimary&&irec==kSpPion)
367                                         fHistMan->GetPtHistogram(kHistPtRecPrimaryAll)->Fill(track->Pt(),dca);  // PT histo of reconstrutsed primaries in defined eta
368                   
369                         //nsigma cut (reconstructed nsigma)
370                                 if(idRec == kSpUndefined) 
371                                         continue;
372                   
373                         // rapidity cut (reconstructed pt and identity)
374                                  if(!fTrackCuts->CheckYCut ((BothParticleSpecies_t)idRec)) continue;
375                   
376                   // Get true ID
377                   
378                   
379                                  if (idRec == idGen) fHistMan->GetHistogram2D(kHistPtRecTrue,  idGen, charge)->Fill(track->Pt(),dca); 
380                   
381                                 if (isPrimary) 
382                                 {
383                                         fHistMan->GetHistogram2D(kHistPtRecSigmaPrimary, idRec, charge)->Fill(track->Pt(),dca); 
384                                         if(idGen != kSpUndefined) 
385                                         {
386                                                 fHistMan->GetHistogram2D(kHistPtRecPrimary,      idGen, charge)->Fill(track->Pt(),dca);
387                                                 if (idRec == idGen) 
388                                                         fHistMan->GetHistogram2D(kHistPtRecTruePrimary,  idGen, charge)->Fill(track->Pt(),dca); 
389                                         }
390                                 }
391                                  //25th Apr - Muons are added to Pions -- FIXME
392                                 if ( pdgcode == 13 && idRec == kSpPion) 
393                                 { 
394                                         fHistMan->GetPtHistogram(kHistPtRecTrueMuonPlus)->Fill(track->Pt(),dca); 
395                                         if(isPrimary)
396                                                 fHistMan->GetPtHistogram(kHistPtRecTruePrimaryMuonPlus)->Fill(track->Pt(),dca); 
397                                 }
398                                 if ( pdgcode == -13 && idRec == kSpPion) 
399                                 { 
400                                         fHistMan->GetPtHistogram(kHistPtRecTrueMuonMinus)->Fill(track->Pt(),dca); 
401                                         if (isPrimary) 
402                                         {
403                                                 fHistMan->GetPtHistogram(kHistPtRecTruePrimaryMuonMinus)->Fill(track->Pt(),dca); 
404                                         }
405                                 }
406                   
407                                  ///..... END FIXME
408                   
409                                 // Fill secondaries
410                                 if(isSecondaryWeak    )  
411                                         fHistMan->GetHistogram2D(kHistPtRecSigmaSecondaryWeakDecay, idRec, charge)->Fill(track->Pt(),dca);
412                                 if(isSecondaryMaterial)  
413                                         fHistMan->GetHistogram2D(kHistPtRecSigmaSecondaryMaterial , idRec, charge)->Fill(track->Pt(),dca);
414                   
415                         }//end if(arrayMC)
416                 }
417         
418         
419         } // end loop on tracks
420
421  // cout<< ntracks<<endl;
422   fHistMan->GetGenMulvsRawMulHistogram("hHistGenMulvsRawMul")->Fill(npar,ntracks);
423   PostData(1, fHistMan  );
424   PostData(2, fEventCuts);
425   PostData(3, fTrackCuts);
426   PostData(4, fPID      );
427 }
428
429 //_________________________________________________________________
430 void   AliAnalysisTaskSpectraBoth::Terminate(Option_t *)
431 {
432   // Terminate
433 }