]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/AliAnalysisTaskContMC.cxx
197508061b93749589f1ef73ff1239b9c96865bb
[u/mrichter/AliRoot.git] / PWGCF / Correlations / DPhi / AliAnalysisTaskContMC.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 //         AliAnalysisTaskContMC 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 "TH3F.h"\r
26 #include "TCanvas.h"\r
27 #include "AliAnalysisTask.h"\r
28 #include "AliAnalysisManager.h"\r
29 #include "AliVTrack.h"\r
30 #include "AliAODMCParticle.h"\r
31 #include "AliVParticle.h"\r
32 #include "AliAODEvent.h"\r
33 #include "AliAODInputHandler.h"\r
34 #include "AliAnalysisTaskContMC.h"\r
35 #include "AliAnalysisTaskESDfilter.h"\r
36 #include "AliAnalysisDataContainer.h"\r
37 #include "AliHelperPID.h"\r
38 #include "AliCentrality.h"\r
39 #include "TProof.h"\r
40 #include "AliPID.h"\r
41 #include "AliVEvent.h"\r
42 #include "AliPIDResponse.h"\r
43 #include "AliStack.h"\r
44 #include <TMCProcess.h>\r
45 \r
46 #include <iostream>\r
47 \r
48 using namespace std;\r
49 \r
50 ClassImp(AliAnalysisTaskContMC) \r
51 \r
52 //________________________________________________________________________\r
53 AliAnalysisTaskContMC::AliAnalysisTaskContMC(const char *name) : AliAnalysisTaskSE(name), fAOD(0), fNSigmaPID(0),fIsMC(0),fOutput(0)\r
54 {\r
55   // Default constructor\r
56   \r
57 \r
58   DefineInput(0, TChain::Class());\r
59   //DefineOutput(1, AliHelperPID::Class());\r
60   DefineOutput(1, TList::Class());\r
61   \r
62 }\r
63 //________________________________________________________________________\r
64 //________________________________________________________________________\r
65 void AliAnalysisTaskContMC::UserCreateOutputObjects()\r
66 {\r
67   Printf("\n\n\n\n\n\n In CreateOutput Object:");\r
68   \r
69   //create output objects\r
70   Printf("NSigma: %.1f",fNSigmaPID->GetNSigmaCut());\r
71   Printf("IsMC: %d",fNSigmaPID->GetisMC());\r
72   fOutput = new TList();\r
73   fOutput->SetOwner();\r
74   fOutput->SetName("list");\r
75   \r
76   fHistID = new TH3F("fHistID","fHistID",6,-1.5,4.5,4,-1.5,2.5,38,0.2,4);\r
77   fOutput->Add(fHistID);\r
78   \r
79   if(!fNSigmaPID)AliFatal("PID object should be set in the steering macro");\r
80   fOutput->Add(fNSigmaPID);\r
81   \r
82   //PostData(1, fNSigmaPID  );\r
83   PostData(1, fOutput  );\r
84   \r
85 }\r
86 \r
87 //________________________________________________________________________\r
88 void AliAnalysisTaskContMC::UserExec(Option_t *)\r
89 {\r
90   const Int_t npart=4;\r
91   const Int_t pdgcode[npart+2]={-1,211,321,2212,11,13};\r
92   const Int_t partid[npart+2]={-1,0,1,2,3,4};\r
93   // main event loop\r
94   fAOD = dynamic_cast<AliAODEvent*>(fInputEvent);\r
95   if (!fAOD) {\r
96     AliWarning("ERROR: AliAODEvent not available \n");\r
97     return;\r
98   }\r
99   \r
100   if (strcmp(fAOD->ClassName(), "AliAODEvent"))\r
101     {\r
102       AliFatal("Not processing AODs");\r
103     }\r
104   //MC Loop\r
105   TClonesArray *arrayMC = 0;\r
106   Printf("fIsMC: %d",fIsMC);\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     }\r
114   Double_t centrality = 0;\r
115 \r
116   AliCentrality *centralityObj = fAOD->GetHeader()->GetCentralityP();\r
117   if (centralityObj)\r
118     {\r
119       centrality = centralityObj->GetCentralityPercentileUnchecked("V0M");\r
120       AliInfo(Form("Centrality is %f", centrality));\r
121     }\r
122   else\r
123     {\r
124       Printf("WARNING: Centrality object is 0");\r
125       centrality = -1;\r
126     }\r
127   if (centrality < 0)\r
128     return;\r
129   \r
130 \r
131   //vertex selection\r
132   Int_t nVertex = ((AliAODEvent*)fAOD)->GetNumberOfVertices();\r
133   if( nVertex > 0 ) { \r
134     AliAODVertex* vertex = (AliAODVertex*)((AliAODEvent*)fAOD)->GetPrimaryVertex();\r
135     //Int_t nTracksPrim = vertex->GetNContributors();\r
136     Double_t zVertex = vertex->GetZ();\r
137     //10 cm cut\r
138     if(TMath::Abs(zVertex)>10) return;\r
139     //AliInfo(Form(" Vertex in = %f with %d particles by  %s data ...",zVertex,nTracksPrim,vertex->GetName()));\r
140     // Reject TPC only vertex\r
141     TString name(vertex->GetName());\r
142     if (name.CompareTo("PrimaryVertex") && name.CompareTo("SPDVertex"))return;\r
143   }  \r
144 \r
145   //Int_t count=0;  \r
146   //track loop\r
147   for (Int_t iTracks = 0; iTracks < fAOD->GetNumberOfTracks(); iTracks++) {\r
148     AliAODTrack* track = fAOD->GetTrack(iTracks);\r
149     if(!(track->TestFilterBit(32)))continue;\r
150     if (TMath::Abs(track->Eta()) > .8 || track->Pt() < .2) continue;\r
151     \r
152     Int_t pdg=-999;\r
153     Int_t isph=-999;\r
154     if (fIsMC && arrayMC) {\r
155       AliAODMCParticle *partMC = (AliAODMCParticle*) arrayMC->At(TMath::Abs(track->GetLabel()));\r
156       if (!partMC) { \r
157         AliError("Cannot get MC particle");\r
158         continue; \r
159       }\r
160       pdg=TMath::Abs(partMC->GetPdgCode());\r
161       isph=partMC->IsPhysicalPrimary();\r
162     }\r
163     if(!isph)continue;\r
164     //step 1, TOF Matching\r
165     UInt_t status;\r
166     status=track->GetStatus();\r
167     if((status&AliVTrack::kTOFout)==0 || (status&AliVTrack::kTIME)==0)continue;\r
168     \r
169     //step 2, combined PID\r
170     Int_t IDTPCTOF=fNSigmaPID->GetParticleSpecies(track,1);\r
171     if(IDTPCTOF==999)IDTPCTOF=-1;\r
172     Int_t IDMC=-1;\r
173     for(Int_t ipart=0;ipart<npart+2;ipart++)if(TMath::Abs(pdg)==pdgcode[ipart])IDMC=partid[ipart];  \r
174     fHistID->Fill(IDMC,IDTPCTOF,track->Pt());\r
175   } // end loop on tracks\r
176   \r
177   PostData(1,fOutput);\r
178   \r
179 }\r
180 \r
181 //_________________________________________________________________\r
182 void   AliAnalysisTaskContMC::Terminate(Option_t *)\r
183 {\r
184   // Terminate analysis\r
185   //\r
186   fOutput = dynamic_cast<TList*>(GetOutputData(1));\r
187   if (!fOutput) {\r
188     printf("ERROR: fOutput not available\n");\r
189     return;\r
190   } \r
191   fHistID = dynamic_cast<TH3F*>(fOutput->FindObject("fHistID"));\r
192   fNSigmaPID = dynamic_cast<AliHelperPID*>(fOutput->FindObject("fNSigmaPID"));\r
193 }\r