]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/SPECTRA/PiKaPr/TestAOD/AliAnalysisTaskSpectraAOD.cxx
track->DCA() intead of PropagateToDCA
[u/mrichter/AliRoot.git] / PWGLF / SPECTRA / PiKaPr / TestAOD / AliAnalysisTaskSpectraAOD.cxx
1 /**************************************************************************\r
2  * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *\r
3  *                                                                        *\r
4  * Author: The ALICE Off-line Project.                                    *\r
5  * Contributors are mentioned in the code where appropriate.              *\r
6  *                                                                        *\r
7  * Permission to use, copy, modify and distribute this software and its   *\r
8  * documentation strictly for non-commercial purposes is hereby granted   *\r
9  * without fee, provided that the above copyright notice appears in all   *\r
10  * copies and that both the copyright notice and this permission notice   *\r
11  * appear in the supporting documentation. The authors make no claims     *\r
12  * about the suitability of this software for any purpose. It is          *\r
13  * provided "as is" without express or implied warranty.                  *\r
14  **************************************************************************/\r
15 \r
16 //-----------------------------------------------------------------\r
17 //         AliAnalysisTaskSpectraAOD class\r
18 //-----------------------------------------------------------------\r
19 \r
20 #include "TChain.h"\r
21 #include "TTree.h"\r
22 #include "TLegend.h"\r
23 #include "TH1F.h"\r
24 #include "TH2F.h"\r
25 #include "TCanvas.h"\r
26 #include "AliAnalysisTask.h"\r
27 #include "AliAnalysisManager.h"\r
28 #include "AliAODTrack.h"\r
29 #include "AliAODMCParticle.h"\r
30 #include "AliVParticle.h"\r
31 #include "AliAODEvent.h"\r
32 #include "AliAODInputHandler.h"\r
33 #include "AliAnalysisTaskSpectraAOD.h"\r
34 #include "AliAnalysisTaskESDfilter.h"\r
35 #include "AliAnalysisDataContainer.h"\r
36 #include "AliSpectraAODHistoManager.h"\r
37 #include "AliSpectraAODTrackCuts.h"\r
38 #include "AliSpectraAODEventCuts.h"\r
39 #include "AliCentrality.h"\r
40 #include "TProof.h"\r
41 #include "AliPID.h"\r
42 #include "AliVEvent.h"\r
43 #include "AliPIDResponse.h"\r
44 #include "AliStack.h"\r
45 #include "AliSpectraAODPID.h"\r
46 #include <TMCProcess.h>\r
47 \r
48 #include <iostream>\r
49 \r
50 \r
51 \r
52 \r
53 using namespace AliSpectraNameSpace;\r
54 using namespace std;\r
55 \r
56 ClassImp(AliAnalysisTaskSpectraAOD) // EX1 This stuff tells root to implement the streamer, inspector methods etc (we discussed about it today)\r
57 \r
58 //________________________________________________________________________\r
59 AliAnalysisTaskSpectraAOD::AliAnalysisTaskSpectraAOD(const char *name) : AliAnalysisTaskSE(name), fAOD(0), fHistMan(0), fTrackCuts(0), fEventCuts(0),  fPID(0), fIsMC(0)\r
60 {\r
61   // Default constructor\r
62 \r
63   DefineInput(0, TChain::Class());\r
64   DefineOutput(1, AliSpectraAODHistoManager::Class());\r
65   DefineOutput(2, AliSpectraAODEventCuts::Class());\r
66   DefineOutput(3, AliSpectraAODTrackCuts::Class());\r
67   DefineOutput(4, AliSpectraAODPID::Class());\r
68 \r
69 }\r
70 //________________________________________________________________________\r
71 //________________________________________________________________________\r
72 void AliAnalysisTaskSpectraAOD::UserCreateOutputObjects()\r
73 {\r
74   // create output objects\r
75   fHistMan = new AliSpectraAODHistoManager("SpectraHistos");\r
76 \r
77   if (!fTrackCuts) AliFatal("Track Cuts should be set in the steering macro");\r
78   if (!fEventCuts) AliFatal("Event Cuts should be set in the steering macro");\r
79   if (!fPID)       AliFatal("PID object should be set in the steering macro");\r
80 \r
81   PostData(1, fHistMan  );\r
82   PostData(2, fEventCuts);\r
83   PostData(3, fTrackCuts);\r
84   PostData(4, fPID      );\r
85 \r
86 }\r
87 //________________________________________________________________________\r
88 void AliAnalysisTaskSpectraAOD::UserExec(Option_t *)\r
89 {\r
90   // main event loop\r
91   fAOD = dynamic_cast<AliAODEvent*>(fInputEvent);\r
92   if (strcmp(fAOD->ClassName(), "AliAODEvent"))\r
93     {\r
94       AliFatal("Not processing AODs");\r
95     }\r
96     \r
97   //check on centrality distribution\r
98   fHistMan->GetPtHistogram("CentCheck")->Fill(fAOD->GetCentrality()->GetCentralityPercentile("V0M"),fAOD->GetHeader()->GetCentralityP()->GetCentralityPercentileUnchecked("V0M"));\r
99   \r
100   if(!fEventCuts->IsSelected(fAOD,fTrackCuts))return;//event selection\r
101   \r
102   //AliCentrality fAliCentral*;\r
103   //   if ((fAOD->GetCentrality())->GetCentralityPercentile("V0M") > 5.) return;\r
104   \r
105   // First do MC to fill up the MC particle array, such that we can use it later\r
106   TClonesArray *arrayMC = 0;\r
107   if (fIsMC)\r
108     {\r
109       arrayMC = (TClonesArray*) fAOD->GetList()->FindObject(AliAODMCParticle::StdBranchName());\r
110       if (!arrayMC) {\r
111         AliFatal("Error: MC particles branch not found!\n");\r
112       }\r
113       Int_t nMC = arrayMC->GetEntries();\r
114       for (Int_t iMC = 0; iMC < nMC; iMC++)\r
115         {\r
116           AliAODMCParticle *partMC = (AliAODMCParticle*) arrayMC->At(iMC);\r
117           if(!partMC->Charge()) continue;//Skip neutrals\r
118           if(TMath::Abs(partMC->Eta()) < fTrackCuts->GetEta()){//charged hadron are filled inside the eta acceptance\r
119             fHistMan->GetPtHistogram(kHistPtGen)->Fill(partMC->Pt(),partMC->IsPhysicalPrimary());\r
120           }\r
121           //rapidity cut\r
122           if(TMath::Abs(partMC->Y())   > fTrackCuts->GetY()  ) continue;            \r
123           // check for true PID + and fill P_t histos \r
124           Int_t charge = partMC->Charge() > 0 ? kChPos : kChNeg ;\r
125           Int_t id = fPID->GetParticleSpecie(partMC);\r
126           if(id != kSpUndefined) {\r
127             fHistMan->GetHistogram2D(kHistPtGenTruePrimary,id,charge)->Fill(partMC->Pt(),partMC->IsPhysicalPrimary());\r
128           }\r
129         }\r
130     }\r
131   \r
132   //main loop on tracks\r
133   for (Int_t iTracks = 0; iTracks < fAOD->GetNumberOfTracks(); iTracks++) {\r
134     AliAODTrack* track = fAOD->GetTrack(iTracks);\r
135     if (!fTrackCuts->IsSelected(track)) continue;\r
136     \r
137     fPID->FillQAHistos(fHistMan, track, fTrackCuts);\r
138     \r
139     // //cut on q vectors track-by-track FIXME\r
140     // if ((qPos>fTrackCuts->GetQvecMin() && qPos<fTrackCuts->GetQvecMax() && track->Eta()<0) || (qNeg>fTrackCuts->GetQvecMin() && qNeg<fTrackCuts->GetQvecMax() && track->Eta()>=0)){\r
141     \r
142     //calculate DCA for AOD track\r
143     //Double_t d[2], covd[3];\r
144     // AliAODTrack* track_clone=(AliAODTrack*)track->Clone("track_clone"); // need to clone because PropagateToDCA updates the track parameters\r
145     // Bool_t isDCA = track_clone->PropagateToDCA(fAOD->GetPrimaryVertex(),fAOD->GetMagneticField(),9999.,d,covd);\r
146     // delete track_clone;\r
147     // if(!isDCA)d[0]=-999;\r
148     Double_t dca=track->DCA();\r
149     \r
150     fHistMan->GetPtHistogram(kHistPtRec)->Fill(track->Pt(),dca);  // PT histo\r
151     \r
152     // get identity and charge\r
153     Int_t idRec  = fPID->GetParticleSpecie(track, fTrackCuts);\r
154     \r
155     Int_t charge = track->Charge() > 0 ? kChPos : kChNeg;\r
156     \r
157     // Fill histograms, only if inside y and nsigma acceptance\r
158     if(idRec != kSpUndefined && fTrackCuts->CheckYCut ((AODParticleSpecies_t)idRec))fHistMan->GetHistogram2D(kHistPtRecSigma,idRec,charge)->Fill(track->Pt(),dca);\r
159     //can't put a continue because we still have to fill allcharged primaries, done later\r
160     \r
161     /* MC Part */\r
162     if (arrayMC) {\r
163       AliAODMCParticle *partMC = (AliAODMCParticle*) arrayMC->At(TMath::Abs(track->GetLabel()));\r
164       if (!partMC) { \r
165         AliError("Cannot get MC particle");\r
166         continue; \r
167       }\r
168       // Check if it is primary, secondary from material or secondary from weak decay\r
169       Bool_t isPrimary           = partMC->IsPhysicalPrimary();\r
170       Bool_t isSecondaryMaterial = kFALSE; \r
171       Bool_t isSecondaryWeak     = kFALSE; \r
172       if(!isPrimary) {\r
173         Int_t mfl=-999,codemoth=-999;\r
174         Int_t indexMoth=partMC->GetMother(); // FIXME ignore fakes? TO BE CHECKED, on ESD is GetFirstMother()\r
175         if(indexMoth>=0){//is not fake\r
176           AliAODMCParticle* moth = (AliAODMCParticle*) arrayMC->At(indexMoth);\r
177           codemoth = TMath::Abs(moth->GetPdgCode());\r
178           mfl = Int_t (codemoth/ TMath::Power(10, Int_t(TMath::Log10(codemoth))));\r
179         }\r
180         // Int_t uniqueID = partMC->GetUniqueID();\r
181         // cout<<"uniqueID: "<<partMC->GetUniqueID()<<"       "<<kPDecay<<endl;\r
182         // cout<<"status: "<<partMC->GetStatus()<<"       "<<kPDecay<<endl;\r
183         // if(uniqueID == kPDecay)Printf("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");\r
184         if(mfl==3) isSecondaryWeak     = kTRUE; // add if(partMC->GetStatus() & kPDecay)? FIXME\r
185         else       isSecondaryMaterial = kTRUE;\r
186       }\r
187       \r
188       if (isPrimary)fHistMan->GetPtHistogram(kHistPtRecPrimary)->Fill(track->Pt(),dca);  // PT histo of primaries\r
189       \r
190       //nsigma cut (reconstructed nsigma)\r
191       if(idRec == kSpUndefined) continue;\r
192       \r
193       // rapidity cut (reconstructed pt and identity)\r
194       if(!fTrackCuts->CheckYCut ((AODParticleSpecies_t)idRec)) continue;\r
195       \r
196       // Get true ID\r
197       Int_t idGen     = fPID->GetParticleSpecie(partMC);\r
198       //if(TMath::Abs(partMC->Y())   > fTrackCuts->GetY()  ) continue;      // FIXME: do we need a rapidity cut on the generated?\r
199       // Fill histograms for primaries\r
200       \r
201       if (idRec == idGen) fHistMan->GetHistogram2D(kHistPtRecTrue,  idGen, charge)->Fill(track->Pt(),dca); \r
202       \r
203       if (isPrimary) {\r
204         fHistMan                    ->GetHistogram2D(kHistPtRecSigmaPrimary, idRec, charge)->Fill(track->Pt(),dca); \r
205         if(idGen != kSpUndefined) {\r
206           fHistMan                    ->GetHistogram2D(kHistPtRecPrimary,      idGen, charge)->Fill(track->Pt(),dca);\r
207           if (idRec == idGen) fHistMan->GetHistogram2D(kHistPtRecTruePrimary,  idGen, charge)->Fill(track->Pt(),dca); \r
208         }\r
209       }\r
210       //25th Apr - Muons are added to Pions -- FIXME\r
211       if ( partMC->PdgCode() == 13 && idRec == kSpPion) { \r
212         fHistMan->GetPtHistogram(kHistPtRecTrueMuonPlus)->Fill(track->Pt(),dca); \r
213         if(isPrimary)\r
214           fHistMan->GetPtHistogram(kHistPtRecTruePrimaryMuonPlus)->Fill(track->Pt(),dca); \r
215       }\r
216       if ( partMC->PdgCode() == -13 && idRec == kSpPion) { \r
217         fHistMan->GetPtHistogram(kHistPtRecTrueMuonMinus)->Fill(track->Pt(),dca); \r
218         if (isPrimary) {\r
219           fHistMan->GetPtHistogram(kHistPtRecTruePrimaryMuonMinus)->Fill(track->Pt(),dca); \r
220         }\r
221       }\r
222       \r
223       ///..... END FIXME\r
224       \r
225       // Fill secondaries\r
226       if(isSecondaryWeak    )  fHistMan->GetHistogram2D(kHistPtRecSigmaSecondaryWeakDecay, idRec, charge)->Fill(track->Pt(),dca);\r
227       if(isSecondaryMaterial)  fHistMan->GetHistogram2D(kHistPtRecSigmaSecondaryMaterial , idRec, charge)->Fill(track->Pt(),dca);\r
228       \r
229     }//end if(arrayMC)\r
230   } // end loop on tracks\r
231   \r
232   PostData(1, fHistMan  );\r
233   PostData(2, fEventCuts);\r
234   PostData(3, fTrackCuts);\r
235   PostData(4, fPID      );\r
236 }\r
237 \r
238 //_________________________________________________________________\r
239 void   AliAnalysisTaskSpectraAOD::Terminate(Option_t *)\r
240 {\r
241   // Terminate\r
242 }\r