]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/totEt/AliAnalysisHadEtMonteCarlo.cxx
- memory leak fixed, created histograms were not deleted
[u/mrichter/AliRoot.git] / PWG4 / totEt / AliAnalysisHadEtMonteCarlo.cxx
1 //_________________________________________________________________________
2 //  Utility Class for transverse energy studies, charged hadrons
3 //  Base class for MC analysis
4 //  - MC output
5 // implementation file
6 //
7 //Created by Christine Nattrass, Rebecca Scott, Irakli Martashvili
8 //University of Tennessee at Knoxville
9 //_________________________________________________________________________
10 #include "AliAnalysisHadEtMonteCarlo.h"
11 #include "AliAnalysisEtCuts.h"
12
13 #include "AliStack.h"
14 #include "AliMCEvent.h"
15 #include "AliESDEvent.h"
16 #include "AliESDtrackCuts.h"
17 #include "AliESDpid.h"
18 #include "AliPID.h"
19 #include "AliESDtrack.h"
20 #include "AliVParticle.h"
21 #include "AliAnalysisTask.h"
22 #include "AliAnalysisManager.h"
23 #include "AliAnalysisHadEtReconstructed.h"
24 #include "AliAnalysisEtCuts.h"
25 #include <iostream>
26 #include "TRandom.h"
27
28 using namespace std;
29
30 ClassImp(AliAnalysisHadEtMonteCarlo);
31
32
33 Int_t AliAnalysisHadEtMonteCarlo::fgNumSmearWidths = 4;
34 Float_t AliAnalysisHadEtMonteCarlo::fgSmearWidths[4] = {0.005,0.006,0.007,0.008};
35
36 AliAnalysisHadEtMonteCarlo::AliAnalysisHadEtMonteCarlo():AliAnalysisHadEt()
37                                                         ,fSimPiKPEt(0)
38                                                         ,fSimHadEt(0)
39                                                         ,fSimTotEt(0) 
40                                                         ,fPtSmearer(0)
41 {
42 //   for(int i=0;i<fgNumSmearWidths;i++){
43 //     //fSimPiKPEtSmeared[i] = 0.0;
44 //   }
45 }
46 AliAnalysisHadEtMonteCarlo::~AliAnalysisHadEtMonteCarlo(){//destructor
47   //if(fSimPiKPEtSmeared) delete [] fSimPiKPEtSmeared;
48   if(fPtSmearer) delete fPtSmearer;
49 }
50
51 void AliAnalysisHadEtMonteCarlo::ResetEventValues(){//resetting event variables
52   AliAnalysisHadEt::ResetEventValues();
53     fSimHadEt=0.0;
54     fSimTotEt=0.0;
55     fSimPiKPEt=0.0;
56 }
57 Int_t AliAnalysisHadEtMonteCarlo::AnalyseEvent(AliVEvent* ev,AliVEvent* ev2)
58 { // analyse MC and real event info
59   FillHisto1D("NEvents",0.5,1);
60
61   AnalyseEvent(ev);
62   AliMCEvent *mcEvent = dynamic_cast<AliMCEvent*>(ev);
63   AliESDEvent *realEvent = dynamic_cast<AliESDEvent*>(ev2);
64   AliStack *stack = mcEvent->Stack();
65
66   //for PID
67   AliESDpid *pID = new AliESDpid();
68   //pID->MakePID(realEvent);
69
70   //This code taken from https://twiki.cern.ch/twiki/bin/view/ALICE/SelectionOfPrimaryTracksForPp2009DataAnalysis
71   //Gets good tracks
72   //=============================================
73   // Primary vertex
74   const AliESDVertex *vertex = realEvent->GetPrimaryVertexTracks();
75   if(vertex->GetNContributors()<1) {
76     // SPD vertex
77     vertex = realEvent->GetPrimaryVertexSPD();
78     if(vertex->GetNContributors()<1) {
79       // NO GOOD VERTEX, SKIP EVENT 
80     }
81   }
82   // apply a cut |zVertex| < CUT, if needed
83
84   
85
86   //fEsdtrackCutsITSTPC->SetEtaRange(-0.8,0.8); // normally, |eta|<0.8
87   //=============================================
88
89   //Roughly following $ALICE_ROOT/PWG0/dNdEta/AlidNdEtaCorrectionTask
90
91   //=============================================TPC&&ITS=============================================
92   TString *strTPC = new TString("TPC");
93   TString *strITS = new TString("ITS");
94   TString *strTPCITS = new TString("TPCITS");
95   for(Int_t cutset=0;cutset<2;cutset++){
96     TString *cutName;
97     TObjArray* list;
98     switch(cutset){
99     case 0:
100       cutName = strTPC;
101       list = fEsdtrackCutsTPC->GetAcceptedTracks(realEvent);
102       break;
103     case 1:
104       cutName = strITS;
105       list = fEsdtrackCutsITS->GetAcceptedTracks(realEvent);
106       break;
107     case 2:
108       cutName = strTPCITS;
109       list = fEsdtrackCutsITSTPC->GetAcceptedTracks(realEvent);
110       break;
111     default:
112       cerr<<"Error:  cannot fill histograms!"<<endl;
113       return -1;
114     }
115     Int_t nGoodTracks = list->GetEntries();
116     for (Int_t iTrack = 0; iTrack < nGoodTracks; iTrack++)
117       {
118         AliESDtrack *track = dynamic_cast<AliESDtrack*> (list->At(iTrack));
119         if (!track)
120           {
121             Printf("ERROR: Could not get track %d", iTrack);
122             continue;
123           }
124         else{
125           Float_t nSigmaPion,nSigmaProton,nSigmaKaon,nSigmaElectron;
126             pID->MakeTPCPID(track);
127             pID->MakeITSPID(track);
128           if(cutset!=1){
129             nSigmaPion = TMath::Abs(pID->NumberOfSigmasTPC(track,AliPID::kPion));
130             nSigmaProton = TMath::Abs(pID->NumberOfSigmasTPC(track,AliPID::kProton));
131             nSigmaKaon = TMath::Abs(pID->NumberOfSigmasTPC(track,AliPID::kKaon));
132             nSigmaElectron = TMath::Abs(pID->NumberOfSigmasTPC(track,AliPID::kElectron));
133           }
134           else{
135             nSigmaPion = TMath::Abs(pID->NumberOfSigmasITS(track,AliPID::kPion));
136             nSigmaProton = TMath::Abs(pID->NumberOfSigmasITS(track,AliPID::kProton));
137             nSigmaKaon = TMath::Abs(pID->NumberOfSigmasITS(track,AliPID::kKaon));
138             nSigmaElectron = TMath::Abs(pID->NumberOfSigmasITS(track,AliPID::kElectron));
139           }
140           bool isPion = (nSigmaPion<3.0 && nSigmaProton>2.0 && nSigmaKaon>2.0);
141           bool isElectron = (nSigmaElectron<2.0 && nSigmaPion>4.0 && nSigmaProton>3.0 && nSigmaKaon>3.0);
142           bool isKaon = (nSigmaPion>3.0 && nSigmaProton>2.0 && nSigmaKaon<2.0);
143           bool isProton = (nSigmaPion>3.0 && nSigmaProton<2.0 && nSigmaKaon>2.0);
144
145           //bool IsElectron = false;
146           bool unidentified = (!isProton && !isKaon && !isElectron);
147           Float_t dEdx = track->GetTPCsignal();
148           if(cutset==1) dEdx = track->GetITSsignal();
149
150           FillHisto2D(Form("dEdxAll%s",cutName->Data()),track->P(),dEdx,1.0);
151           //if(cutset==1) cout<<"filling "<<track->P()<<" "<<dEdx<<endl;
152
153           UInt_t label = (UInt_t)TMath::Abs(track->GetLabel());
154           TParticle  *simPart  = stack->Particle(label);
155           if(!simPart) {
156             Printf("no MC particle\n");                 
157             continue;           
158           }
159           else{//analysis
160             if(stack->IsPhysicalPrimary(label)){
161               if (TMath::Abs(simPart->Eta()) < fCuts->GetCommonEtaCut()){
162               //if (TMath::Abs(simPart->Eta()) < 0.7)       {
163
164                 Int_t pdgCode =  simPart->GetPDG(0)->PdgCode();
165                 Int_t mypid = 0;
166                 if(pdgCode==fgPiPlusCode) mypid = 1;
167                 if(pdgCode==fgProtonCode) mypid = 2;
168                 if(pdgCode==fgKPlusCode) mypid = 3;
169                 if(pdgCode==fgEPlusCode) mypid = 4;
170                 if(pdgCode==fgPiMinusCode) mypid = 1;
171                 if(pdgCode==fgAntiProtonCode) mypid = 2;
172                 if(pdgCode==fgKMinusCode) mypid = 3;
173                 if(pdgCode==fgEMinusCode) mypid = 4;
174                 //cout<<pdgCode->PdgCode()<<" ";
175                 //fPdgDB->GetSimParticle("pi+")->PdgCode();
176                 bool filled = false;      
177                 //============Charged hadrons===================================
178                 //identified...
179                 if(isPion){
180                   if(pdgCode!=fgPiPlusCode && pdgCode!=fgPiMinusCode){
181                     FillHisto2D(Form("MisidentifiedPIDs%s",cutName->Data()),1,mypid,1);
182                     //if(mypid==0)cerr<<"I was misidentified! I'm not a pion! I am a "<<simPart->GetName()<<endl;
183                   }
184                   float myEt = Et(simPart);
185                   if(track->Charge()>0){ FillHisto2D(Form("EtReconstructed%sIdentifiedPiPlus",cutName->Data()),track->Pt(),track->Eta(),myEt);}
186                   else{ FillHisto2D(Form("EtReconstructed%sIdentifiedPiMinus",cutName->Data()),track->Pt(),track->Eta(),myEt);}
187                   FillHisto2D(Form("dEdxPion%s",cutName->Data()),track->P(),dEdx,1.0);
188                 }
189                 if(isProton){
190                   if(pdgCode!=fgProtonCode && pdgCode!=fgAntiProtonCode){
191                     FillHisto2D(Form("MisidentifiedPIDs%s",cutName->Data()),2,mypid,1);
192                     // if(mypid==0)cerr<<"I was misidentified!  I'm not a proton! I am a "<<simPart->GetName()<<endl;
193                   }
194                   float myEt = Et(simPart);
195                   if(track->Charge()>0){ FillHisto2D(Form("EtReconstructed%sIdentifiedProton",cutName->Data()),track->Pt(),track->Eta(),myEt);}
196                   else{ FillHisto2D(Form("EtReconstructed%sIdentifiedAntiProton",cutName->Data()),track->Pt(),track->Eta(),myEt);}
197                   FillHisto2D(Form("dEdxProton%s",cutName->Data()),track->P(),dEdx,1.0);
198                 }
199                 if(isKaon){
200                   if(pdgCode!=fgKMinusCode && pdgCode!=fgKPlusCode){
201                     FillHisto2D(Form("MisidentifiedPIDs%s",cutName->Data()),3,mypid,1);
202                     //if(mypid==0)cerr<<"I was misidentified!  I'm not a kaon! I am a "<<simPart->GetName()<<" p "<<track->P()<<" nSigmaProton "<<nSigmaProton<<" nSigmaPion "<<nSigmaPion<<" nSigmaKaon "<<nSigmaKaon<<endl;
203                   }
204                   float myEt = Et(simPart);
205                   if(track->Charge()>0){ FillHisto2D(Form("EtReconstructed%sIdentifiedKPlus",cutName->Data()),track->Pt(),track->Eta(),myEt);}
206                   else{ FillHisto2D(Form("EtReconstructed%sIdentifiedKMinus",cutName->Data()),track->Pt(),track->Eta(),myEt);}
207                   FillHisto2D(Form("dEdxKaon%s",cutName->Data()),track->P(),dEdx,1.0);
208                 }
209                 if(isElectron){
210                   if(pdgCode!=fgEMinusCode && pdgCode!=fgEPlusCode){
211                     FillHisto2D(Form("MisidentifiedPIDs%s",cutName->Data()),4,mypid,1);
212                     //cerr<<"I was misidentified!  I'm not an electron! I am a "<<simPart->GetName()<<endl;
213                   }
214                   float myEt = Et(simPart);
215                   if(track->Charge()>0){ FillHisto2D(Form("EtReconstructed%sIdentifiedEPlus",cutName->Data()),track->Pt(),track->Eta(),myEt);}
216                   else{ FillHisto2D(Form("EtReconstructed%sIdentifiedEMinus",cutName->Data()),track->Pt(),track->Eta(),myEt);}
217                   FillHisto2D(Form("dEdxElectron%s",cutName->Data()),track->P(),dEdx,1.0);
218                 }
219                 if(unidentified){
220                   if(pdgCode!=fgEMinusCode && pdgCode!=fgEPlusCode){
221                     float myEtPi = Et(simPart,fgPionMass);
222                     float myEt = Et(simPart);
223                     FillHisto2D(Form("EtReconstructed%sUnidentifiedAssumingPion",cutName->Data()),track->Pt(),track->Eta(),myEtPi);
224                     FillHisto2D(Form("EtReconstructed%sUnidentified",cutName->Data()),track->Pt(),track->Eta(),myEt);
225                     FillHisto2D(Form("EtNReconstructed%sUnidentified",cutName->Data()),track->Pt(),track->Eta(),1.0);
226                   }
227                   FillHisto2D(Form("dEdxUnidentified%s",cutName->Data()),track->P(),dEdx,1.0);
228                   //cout<<"I was not identified.  I am a "<<simPart->GetName()<<" PID "<<pdgCode<<endl;
229                   //track what was not identified successfully
230                   FillHisto1D(Form("UnidentifiedPIDs%s",cutName->Data()),mypid,1);
231                 }
232                 //...simulated
233                 if(pdgCode == fgPiPlusCode){
234                   //cout<<"I'm a real primary "<<simPart->GetName()<<"! "<<"my label is "<<simPart->GetFirstMother()<<" track no "<<iTrack<<"/"<<realEvent->GetNumberOfTracks()<<endl;//<<" "<<label<<" "<<pdgCode<<endl;
235                 
236                   float myEt = Et(simPart);
237                   FillHisto2D(Form("EtReconstructed%sPiPlus",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
238                   FillHisto2D(Form("EtReconstructed%sChargedHadron",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
239                   FillHisto2D(Form("EtNReconstructed%sPiPlus",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
240                   FillHisto2D(Form("EtNReconstructed%sChargedHadron",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
241                   FillHisto2D(Form("EtReconstructed%sChargedHadronAssumingPion",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
242                   filled = true;
243                 }
244                 if(pdgCode == fgPiMinusCode){
245                   float myEt = Et(simPart);
246                   FillHisto2D(Form("EtReconstructed%sPiMinus",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
247                   FillHisto2D(Form("EtReconstructed%sChargedHadron",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
248                   FillHisto2D(Form("EtNReconstructed%sPiMinus",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
249                   FillHisto2D(Form("EtNReconstructed%sChargedHadron",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
250                   FillHisto2D(Form("EtReconstructed%sChargedHadronAssumingPion",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
251                   filled = true;
252                 }
253                 if(pdgCode == fgKPlusCode){
254                   float myEt = Et(simPart);
255                   float myEtPi = Et(simPart,fgPionMass);
256                   FillHisto2D(Form("EtReconstructed%sKPlus",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
257                   FillHisto2D(Form("EtReconstructed%sChargedHadron",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
258                   FillHisto2D(Form("EtNReconstructed%sKPlus",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
259                   FillHisto2D(Form("EtNReconstructed%sChargedHadron",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
260                   FillHisto2D(Form("EtReconstructed%sChargedHadronAssumingPion",cutName->Data()),simPart->Pt(),simPart->Eta(),myEtPi);
261                   FillHisto2D(Form("EtReconstructed%sKPlusAssumingPion",cutName->Data()),simPart->Pt(),simPart->Eta(),myEtPi);
262                   filled = true;
263                 }
264                 if(pdgCode == fgKMinusCode){
265                   float myEt = Et(simPart);
266                   float myEtPi = Et(simPart,fgPionMass);
267                   FillHisto2D(Form("EtReconstructed%sKMinus",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
268                   FillHisto2D(Form("EtReconstructed%sChargedHadron",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
269                   FillHisto2D(Form("EtNReconstructed%sKMinus",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
270                   FillHisto2D(Form("EtNReconstructed%sChargedHadron",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
271                   FillHisto2D(Form("EtReconstructed%sChargedHadronAssumingPion",cutName->Data()),simPart->Pt(),simPart->Eta(),myEtPi);
272                   FillHisto2D(Form("EtReconstructed%sKMinusAssumingPion",cutName->Data()),simPart->Pt(),simPart->Eta(),myEtPi);
273                   filled = true;
274                 }
275                 if(pdgCode == fgProtonCode){
276                   float myEt = Et(simPart);
277                   float myEtPi = Et(simPart,fgPionMass);
278                   FillHisto2D(Form("EtReconstructed%sProton",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
279                   FillHisto2D(Form("EtReconstructed%sChargedHadron",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
280                   FillHisto2D(Form("EtNReconstructed%sProton",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
281                   FillHisto2D(Form("EtNReconstructed%sChargedHadron",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
282                   FillHisto2D(Form("EtReconstructed%sChargedHadronAssumingPion",cutName->Data()),simPart->Pt(),simPart->Eta(),myEtPi);
283                   FillHisto2D(Form("EtReconstructed%sProtonAssumingPion",cutName->Data()),simPart->Pt(),simPart->Eta(),myEtPi);
284                   filled = true;
285                 }
286                 if(pdgCode == fgAntiProtonCode){
287                   float myEt = Et(simPart);
288                   float myEtPi = Et(simPart,fgPionMass);
289                   FillHisto2D(Form("EtReconstructed%sAntiProton",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
290                   FillHisto2D(Form("EtReconstructed%sChargedHadron",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
291                   FillHisto2D(Form("EtNReconstructed%sAntiProton",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
292                   FillHisto2D(Form("EtNReconstructed%sChargedHadron",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
293                   FillHisto2D(Form("EtReconstructed%sChargedHadronAssumingPion",cutName->Data()),simPart->Pt(),simPart->Eta(),myEtPi);
294                   FillHisto2D(Form("EtReconstructed%sAntiProtonAssumingPion",cutName->Data()),simPart->Pt(),simPart->Eta(),myEtPi);
295                   filled = true;
296                 }
297                 if(pdgCode == fgEPlusCode){
298                   float myEt = Et(simPart);
299                   FillHisto2D(Form("EtReconstructed%sEPlus",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
300                   if(!isElectron || unidentified){
301                     float myEtPi = Et(simPart,fgPionMass);
302                     FillHisto2D(Form("EtReconstructed%sMisidentifiedElectrons",cutName->Data()),track->Pt(),track->Eta(),myEtPi);
303                   }
304                   filled = true;
305                 }
306                 if(pdgCode == fgEMinusCode){
307                   if(!isElectron || unidentified){
308                     float myEtPi = Et(simPart,fgPionMass);
309                     FillHisto2D(Form("EtReconstructed%sMisidentifiedElectrons",cutName->Data()),track->Pt(),track->Eta(),myEtPi);
310                   }
311                   float myEt = Et(simPart);
312                     FillHisto2D(Form("EtReconstructed%sEMinus",cutName->Data()),simPart->Pt(),simPart->Eta(),myEt);
313                     filled = true;
314                 }
315                 //if(!filled){
316                   //TParticlePDG *pc = simPart->GetPDG(0);
317                   //if( strcmp(pc->ParticleClass(),"Baryon")==0 || strcmp(pc->ParticleClass(),"Meson")==0 ){
318                   //cout<<"Did not find a place for "<<simPart->GetName()<<" "<<pdgCode<<" which is a "<<pc->ParticleClass()<<endl;
319                   //}
320                   //}
321               }
322               
323             }
324             else{//not a primary - we're after V0 daughters!
325               //cout<<"I'm a secondary "<<simPart->GetName()<<"!";//<<endl;
326               if (TMath::Abs(simPart->Eta()) < fCuts->GetCommonEtaCut()){
327                 TParticle *mom = stack->Particle(simPart->GetFirstMother());
328                 if(mom){
329                   TParticlePDG *pc = mom->GetPDG(0);
330                   if(pc){
331                     Int_t pdgCode =  mom->GetPDG(0)->PdgCode();
332                     if(pdgCode == fgLambdaCode){
333                       float myEt = Et(simPart);
334                       FillHisto2D(Form("EtReconstructed%sLambdaDaughters",cutName->Data()),track->Pt(),track->Eta(),myEt);
335                     }
336                     if(pdgCode == fgAntiLambdaCode){
337                       float myEt = Et(simPart);
338                       FillHisto2D(Form("EtReconstructed%sAntiLambdaDaughters",cutName->Data()),track->Pt(),track->Eta(),myEt);
339                     }
340                     if(pdgCode == fgK0SCode){
341                       float myEt = Et(simPart);
342                       FillHisto2D(Form("EtReconstructed%sK0SDaughters",cutName->Data()),track->Pt(),track->Eta(),myEt);
343                     }
344                     if(pdgCode == fgXiCode){
345                       float myEt = Et(simPart);
346                       FillHisto2D(Form("EtReconstructed%sXiDaughters",cutName->Data()),track->Pt(),track->Eta(),myEt);
347                     }
348                     if(pdgCode == fgAntiXiCode){
349                       float myEt = Et(simPart);
350                       FillHisto2D(Form("EtReconstructed%sAntiXiDaughters",cutName->Data()),track->Pt(),track->Eta(),myEt);
351                     }
352                     if(pdgCode == fgOmegaCode){
353                       float myEt = Et(simPart);
354                       FillHisto2D(Form("EtReconstructed%sOmegaDaughters",cutName->Data()),track->Pt(),track->Eta(),myEt);
355                     }
356                     if(pdgCode == fgXiCode){
357                       float myEt = Et(simPart);
358                       FillHisto2D(Form("EtReconstructed%sAntiOmegaDaughters",cutName->Data()),track->Pt(),track->Eta(),myEt);
359                     }
360
361                     if(mom->GetFirstMother()>0){
362                       TParticle *grandma = stack->Particle(mom->GetFirstMother());
363                       if(grandma){
364                         Int_t pdgCodeMom =  mom->GetPDG(0)->PdgCode();
365                         if(pdgCodeMom==fgPiPlusCode || pdgCodeMom==fgPiMinusCode || pdgCodeMom==fgProtonCode ||pdgCodeMom==fgAntiProtonCode || pdgCodeMom==fgKPlusCode || pdgCode==fgKMinusCode){
366                           //cout<<" my grandmother is "<<grandma->GetName()<<" "<<endl;
367                           Int_t pdgCodeGrandma =  grandma->GetPDG(0)->PdgCode();
368                       
369                           if(pdgCodeGrandma == fgXiCode){
370                             float myEt = Et(simPart);
371                             FillHisto2D(Form("EtReconstructed%sXiDaughters",cutName->Data()),track->Pt(),track->Eta(),myEt);
372                           }
373                           if(pdgCodeGrandma == fgAntiXiCode){
374                             float myEt = Et(simPart);
375                             FillHisto2D(Form("EtReconstructed%sAntiXiDaughters",cutName->Data()),track->Pt(),track->Eta(),myEt);
376                           }
377                           if(pdgCodeGrandma == fgOmegaCode){
378                             float myEt = Et(simPart);
379                             FillHisto2D(Form("EtReconstructed%sOmegaDaughters",cutName->Data()),track->Pt(),track->Eta(),myEt);
380                           }
381                           if(pdgCodeGrandma == fgXiCode){
382                             float myEt = Et(simPart);
383                             FillHisto2D(Form("EtReconstructed%sAntiOmegaDaughters",cutName->Data()),track->Pt(),track->Eta(),myEt);
384                           }
385
386                         }
387                       }
388                     }
389                   }
390                 }
391               }
392             }
393           }
394
395         }
396       }
397     delete list;
398   }
399   delete pID;
400   delete strTPC;
401   delete strITS;
402   delete strTPCITS;
403   //delete AliESDpid;
404   return 1;
405 }
406 Int_t AliAnalysisHadEtMonteCarlo::AnalyseEvent(AliVEvent* ev)
407 { // analyse MC event
408      ResetEventValues();
409      
410     // Get us an mc event
411     AliMCEvent *mcEvent = dynamic_cast<AliMCEvent*>(ev);
412
413     // Let's play with the stack!
414     AliStack *stack = mcEvent->Stack();
415
416     Int_t nPrim = stack->GetNtrack();
417
418     Float_t fSimPiKPEtPtSmeared = 0;
419     Float_t fSimPiKPEtEfficiencySmeared = 0;
420     Float_t fSimPiKPEtPtCutSmearedTPC = 0;
421     Float_t fSimPiKPEtPtCutSmearedITS = 0;
422     Float_t fSimPiKPEtPIDSmeared = 0;
423     Float_t fSimPiKPEtPIDSmearedNoID = 0;
424     //=================Tracks which may or may not have been reconstructed=================
425
426     for (Int_t iPart = 0; iPart < nPrim; iPart++)
427     {
428
429         TParticle *part = stack->Particle(iPart);
430
431         if (!part)
432           {
433             Printf("ERROR: Could not get particle %d", iPart);
434             continue;
435           }
436
437         TParticlePDG *pc = part->GetPDG(0);
438
439         // Check if it is a primary particle
440         if (stack->IsPhysicalPrimary(iPart)){//primaries
441
442           if (TMath::Abs(part->Eta()) < fCuts->GetCommonEtaCut())           {
443
444             Int_t pdgCode =  part->GetPDG(0)->PdgCode();
445             //cout<<pdgCode->PdgCode()<<" ";
446             //fPdgDB->GetParticle("pi+")->PdgCode();
447             bool filled = false;
448             //Investigating smearing...
449             //Numbers are realistic correction factors from previous studies
450             if(pdgCode==fgPiPlusCode ||pdgCode==fgPiMinusCode ||pdgCode==fgKPlusCode ||pdgCode==fgKMinusCode ||pdgCode==fgProtonCode ||pdgCode==fgAntiProtonCode){
451               //To investigate Smearing...
452               Float_t myet = Et(part);
453               fSimPiKPEt += myet;
454               Float_t theta = part->Theta();
455               Short_t charge = 1;
456               Float_t momentum = part->P();
457               //pt smearing
458               Float_t pSmeared = momentum *  fPtSmearer->Gaus(1,0.005);//Gaussian centered around 1
459               fSimPiKPEtPtSmeared += Et(pSmeared,theta,pdgCode,charge);
460               //Efficiency smearing
461               float efficiency = 2.26545*TMath::Exp(-TMath::Power(9.99977e-01/part->Pt(),7.85488e-02));//simple rough efficiency from fitting curve
462               if(fPtSmearer->Binomial(1,efficiency) ==1){
463                 fSimPiKPEtEfficiencySmeared += (1.0/efficiency)*myet;
464               }
465               //pT cut smeared
466               if(part->Pt()>0.10){fSimPiKPEtPtCutSmearedITS +=1.00645645*myet;}
467               if(part->Pt()>0.15){fSimPiKPEtPtCutSmearedTPC +=1.02000723*myet;}
468               //PID smearing
469               fSimPiKPEtPIDSmearedNoID += 1.02679314*Et(momentum,theta,fgPiPlusCode,charge);
470               if(part->P()<1.0){//then the particle would have been ID'd
471                 fSimPiKPEtPIDSmeared += 1.0085942*myet;
472               }
473               else{//Then it would have been assumed to be a pion
474                 fSimPiKPEtPIDSmeared += 1.0085942*Et(momentum,theta,fgPiPlusCode,charge);
475               }
476             }
477
478             //============Charged hadrons===================================
479             if(pdgCode == fgPiPlusCode){
480               //cout<<"I'm a simulated primary "<<part->GetName()<<"! "<<"my label is "<<part->GetFirstMother()<<" pt "<<part->Pt()<<endl;
481               float myEt = Et(part);
482
483               fSimHadEt += myEt;
484               fSimTotEt += myEt;
485
486               FillHisto2D("EtSimulatedPiPlus",part->Pt(),part->Eta(),myEt);
487               FillHisto2D("EtNSimulatedPiPlus",part->Pt(),part->Eta(),1.0);
488               FillHisto2D("EtSimulatedChargedHadron",part->Pt(),part->Eta(),myEt);
489               FillHisto2D("EtNSimulatedChargedHadron",part->Pt(),part->Eta(),1.0);
490               FillHisto2D("EtSimulatedChargedHadronAssumingPion",part->Pt(),part->Eta(),myEt);
491               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
492               Short_t charge = 1;
493               Float_t myEtLow = Et(0.0,part->Theta(),pdgCode,charge);
494               Float_t myEtITS = Et(0.10/TMath::Sin(part->Theta()),part->Theta(),pdgCode,charge);
495               Float_t myEtTPC = Et(0.15/TMath::Sin(part->Theta()),part->Theta(),pdgCode,charge);
496               FillHisto2D("EtSimulatedChargedHadronAssumingNoPt",part->Pt(),part->Eta(),myEtLow);
497               FillHisto2D("EtSimulatedChargedHadronAssumingPtTPCCut",part->Pt(),part->Eta(),myEtTPC);
498               FillHisto2D("EtSimulatedChargedHadronAssumingPtITSCut",part->Pt(),part->Eta(),myEtITS);
499               filled = true;
500             }
501             if(pdgCode == fgPiMinusCode){
502               float myEt = Et(part);
503               fSimHadEt += myEt;
504               fSimTotEt += myEt;
505               FillHisto2D("EtSimulatedPiMinus",part->Pt(),part->Eta(),myEt);
506               FillHisto2D("EtNSimulatedPiMinus",part->Pt(),part->Eta(),1.0);
507               FillHisto2D("EtSimulatedChargedHadron",part->Pt(),part->Eta(),myEt);
508               FillHisto2D("EtNSimulatedChargedHadron",part->Pt(),part->Eta(),1.0);
509               FillHisto2D("EtSimulatedChargedHadronAssumingPion",part->Pt(),part->Eta(),myEt);
510               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
511               Short_t charge = -1;
512               Float_t myEtLow = Et(0.0,part->Theta(),pdgCode,charge);
513               Float_t myEtITS = Et(0.10/TMath::Sin(part->Theta()),part->Theta(),pdgCode,charge);
514               Float_t myEtTPC = Et(0.15/TMath::Sin(part->Theta()),part->Theta(),pdgCode,charge);
515               FillHisto2D("EtSimulatedChargedHadronAssumingNoPt",part->Pt(),part->Eta(),myEtLow);
516               FillHisto2D("EtSimulatedChargedHadronAssumingPtTPCCut",part->Pt(),part->Eta(),myEtTPC);
517               FillHisto2D("EtSimulatedChargedHadronAssumingPtITSCut",part->Pt(),part->Eta(),myEtITS);
518               filled = true;
519             }
520             if(pdgCode == fgKPlusCode){
521               float myEt = Et(part);
522               float myEtPi = Et(part,fgPionMass);
523               fSimHadEt += myEt;
524               fSimTotEt += myEt;
525               FillHisto2D("EtSimulatedKPlus",part->Pt(),part->Eta(),myEt);
526               FillHisto2D("EtNSimulatedKPlus",part->Pt(),part->Eta(),1.0);
527               FillHisto2D("EtSimulatedChargedHadron",part->Pt(),part->Eta(),myEt);
528               FillHisto2D("EtNSimulatedChargedHadron",part->Pt(),part->Eta(),1.0);
529               FillHisto2D("EtSimulatedChargedHadronAssumingPion",part->Pt(),part->Eta(),myEtPi);
530               FillHisto2D("EtSimulatedKPlusAssumingPion",part->Pt(),part->Eta(),myEtPi);
531               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
532               Short_t charge = 1;
533               Float_t myEtLow = Et(0.0,part->Theta(),pdgCode,charge);
534               Float_t myEtITS = Et(0.10/TMath::Sin(part->Theta()),part->Theta(),pdgCode,charge);
535               Float_t myEtTPC = Et(0.15/TMath::Sin(part->Theta()),part->Theta(),pdgCode,charge);
536               FillHisto2D("EtSimulatedChargedHadronAssumingNoPt",part->Pt(),part->Eta(),myEtLow);
537               FillHisto2D("EtSimulatedChargedHadronAssumingPtTPCCut",part->Pt(),part->Eta(),myEtTPC);
538               FillHisto2D("EtSimulatedChargedHadronAssumingPtITSCut",part->Pt(),part->Eta(),myEtITS);
539               filled = true;
540             }
541             if(pdgCode == fgKMinusCode){
542               float myEt = Et(part);
543               float myEtPi = Et(part,fgPionMass);
544               fSimHadEt += myEt;
545               fSimTotEt += myEt;
546               FillHisto2D("EtSimulatedKMinus",part->Pt(),part->Eta(),myEt);
547               FillHisto2D("EtNSimulatedKMinus",part->Pt(),part->Eta(),1.0);
548               FillHisto2D("EtSimulatedChargedHadron",part->Pt(),part->Eta(),myEt);
549               FillHisto2D("EtNSimulatedChargedHadron",part->Pt(),part->Eta(),1.0);
550               FillHisto2D("EtSimulatedChargedHadronAssumingPion",part->Pt(),part->Eta(),myEtPi);
551               FillHisto2D("EtSimulatedKMinusAssumingPion",part->Pt(),part->Eta(),myEtPi);
552               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
553               Short_t charge = -1;
554               Float_t myEtLow = Et(0.0,part->Theta(),pdgCode,charge);
555               Float_t myEtITS = Et(0.10/TMath::Sin(part->Theta()),part->Theta(),pdgCode,charge);
556               Float_t myEtTPC = Et(0.15/TMath::Sin(part->Theta()),part->Theta(),pdgCode,charge);
557               FillHisto2D("EtSimulatedChargedHadronAssumingNoPt",part->Pt(),part->Eta(),myEtLow);
558               FillHisto2D("EtSimulatedChargedHadronAssumingPtTPCCut",part->Pt(),part->Eta(),myEtTPC);
559               FillHisto2D("EtSimulatedChargedHadronAssumingPtITSCut",part->Pt(),part->Eta(),myEtITS);
560               filled = true;
561             }
562             if(pdgCode == fgProtonCode){
563               float myEt = Et(part);
564               float myEtPi = Et(part,fgPionMass);
565               fSimHadEt += myEt;
566               fSimTotEt += myEt;
567               FillHisto2D("EtSimulatedProton",part->Pt(),part->Eta(),myEt);
568               FillHisto2D("EtNSimulatedProton",part->Pt(),part->Eta(),1.0);
569               FillHisto2D("EtSimulatedChargedHadron",part->Pt(),part->Eta(),myEt);
570               FillHisto2D("EtNSimulatedChargedHadron",part->Pt(),part->Eta(),1.0);
571               FillHisto2D("EtSimulatedChargedHadronAssumingPion",part->Pt(),part->Eta(),myEtPi);
572               FillHisto2D("EtSimulatedProtonAssumingPion",part->Pt(),part->Eta(),myEtPi);
573               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
574               Short_t charge = 1;
575               Float_t myEtLow = Et(0.0,part->Theta(),pdgCode,charge);
576               Float_t myEtITS = Et(0.10/TMath::Sin(part->Theta()),part->Theta(),pdgCode,charge);
577               Float_t myEtTPC = Et(0.15/TMath::Sin(part->Theta()),part->Theta(),pdgCode,charge);
578               FillHisto2D("EtSimulatedChargedHadronAssumingNoPt",part->Pt(),part->Eta(),myEtLow);
579               FillHisto2D("EtSimulatedChargedHadronAssumingPtTPCCut",part->Pt(),part->Eta(),myEtTPC);
580               FillHisto2D("EtSimulatedChargedHadronAssumingPtITSCut",part->Pt(),part->Eta(),myEtITS);
581               filled = true;
582             }
583             if(pdgCode == fgAntiProtonCode){
584               float myEt = Et(part);
585               float myEtPi = Et(part,fgPionMass);
586               fSimHadEt += myEt;
587               fSimTotEt += myEt;
588               FillHisto2D("EtSimulatedAntiProton",part->Pt(),part->Eta(),myEt);
589               FillHisto2D("EtNSimulatedAntiProton",part->Pt(),part->Eta(),1.0);
590               FillHisto2D("EtSimulatedChargedHadron",part->Pt(),part->Eta(),myEt);
591               FillHisto2D("EtNSimulatedChargedHadron",part->Pt(),part->Eta(),1.0);
592               FillHisto2D("EtSimulatedChargedHadronAssumingPion",part->Pt(),part->Eta(),myEtPi);
593               FillHisto2D("EtSimulatedAntiProtonAssumingPion",part->Pt(),part->Eta(),myEtPi);
594               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
595               Short_t charge = -1;
596               Float_t myEtLow = Et(0.0,part->Theta(),pdgCode,charge);
597               Float_t myEtITS = Et(0.10/TMath::Sin(part->Theta()),part->Theta(),pdgCode,charge);
598               Float_t myEtTPC = Et(0.15/TMath::Sin(part->Theta()),part->Theta(),pdgCode,charge);
599               FillHisto2D("EtSimulatedChargedHadronAssumingNoPt",part->Pt(),part->Eta(),myEtLow);
600               FillHisto2D("EtSimulatedChargedHadronAssumingPtTPCCut",part->Pt(),part->Eta(),myEtTPC);
601               FillHisto2D("EtSimulatedChargedHadronAssumingPtITSCut",part->Pt(),part->Eta(),myEtITS);
602               filled = true;
603             }
604             //============Other hadrons===================================
605
606             if(pdgCode == fgNeutronCode){
607               float myEt = Et(part);
608               fSimHadEt += myEt;
609               fSimTotEt += myEt;
610               FillHisto2D("EtSimulatedNeutron",part->Pt(),part->Eta(),myEt);
611               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
612               filled = true;
613             }
614             if(pdgCode == fgAntiNeutronCode){
615               float myEt = Et(part);
616               fSimHadEt += myEt;
617               fSimTotEt += myEt;
618               FillHisto2D("EtSimulatedAntiNeutron",part->Pt(),part->Eta(),myEt);
619               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
620               filled = true;
621             }
622             if(pdgCode == fgLambdaCode){
623               float myEt = Et(part);
624               fSimHadEt += myEt;
625               fSimTotEt += myEt;
626               //cout<<"I am a simulated lambda! pt "<<part->Pt()<<" eta "<<part->Eta()<<endl;
627               FillHisto2D("EtSimulatedLambda",part->Pt(),part->Eta(),myEt);
628               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
629               Int_t ndaughters = part->GetNDaughters();
630               for(Int_t idaughter = 0;idaughter<ndaughters;idaughter++){
631                 Int_t daughterindex = part->GetDaughter(idaughter);
632                 if(daughterindex<0 || daughterindex>1e5) continue;
633                 TParticle *daughter = stack->ParticleFromTreeK(daughterindex);
634                 if(daughter){
635                   if(daughter->GetPDG(0)){
636                     Int_t daughtercode = daughter->GetPDG(0)->PdgCode();
637                     if(daughtercode==fgPiMinusCode || daughtercode==fgProtonCode){
638                       myEt = Et(daughter);
639                       FillHisto2D("EtSimulatedLambdaDaughters",daughter->Pt(),daughter->Eta(),myEt);
640                       //cout<<"Lambda daughter is a "<<daughter->GetName()<<endl;
641                     }
642                   }
643                   else{
644                     //cout<<"Lambda daughter is a "<<daughter->GetName()<<endl;
645                   }
646                 }
647               }
648               filled = true;
649             }
650             if(pdgCode == fgAntiLambdaCode){
651               float myEt = Et(part);
652               fSimHadEt += myEt;
653               fSimTotEt += myEt;
654               FillHisto2D("EtSimulatedAntiLambda",part->Pt(),part->Eta(),myEt);
655               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
656               Int_t ndaughters = part->GetNDaughters();
657               for(Int_t idaughter = 0;idaughter<ndaughters;idaughter++){
658                 Int_t daughterindex = part->GetDaughter(idaughter);
659                 if(daughterindex<0 || daughterindex>1e5) continue;
660                 TParticle *daughter = stack->ParticleFromTreeK(daughterindex);
661                 if(daughter){
662                   if(daughter->GetPDG(0)){
663                     Int_t daughtercode = daughter->GetPDG(0)->PdgCode();
664                     if(daughtercode==fgPiPlusCode || daughtercode==fgAntiProtonCode){
665                       myEt = Et(daughter);
666                       FillHisto2D("EtSimulatedAntiLambdaDaughters",daughter->Pt(),daughter->Eta(),myEt);
667                       //cout<<"AntiLambda daughter is a "<<daughter->GetName()<<endl;
668                     }
669                   }
670                   else{
671                     //cout<<"AntiLambda daughter is a "<<daughter->GetName()<<endl;
672                   }
673                 }
674               }
675               filled = true;
676             }
677             if(pdgCode == fgK0SCode){
678               float myEt = Et(part);
679               fSimHadEt += myEt;
680               fSimTotEt += myEt;
681               FillHisto2D("EtSimulatedK0S",part->Pt(),part->Eta(),myEt);
682               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
683               Int_t ndaughters = part->GetNDaughters();
684               for(Int_t idaughter = 0;idaughter<ndaughters;idaughter++){
685                 Int_t daughterindex = part->GetDaughter(idaughter);
686                 if(daughterindex<0 || daughterindex>1e5) continue;
687                 TParticle *daughter = stack->ParticleFromTreeK(daughterindex);
688                 if(daughter){
689                   if(daughter->GetPDG(0)){
690
691                     Int_t daughtercode = daughter->GetPDG(0)->PdgCode();
692                     if(daughtercode==fgPiMinusCode || daughtercode==fgPiPlusCode){
693                       myEt = Et(daughter);
694                       FillHisto2D("EtSimulatedK0SDaughters",daughter->Pt(),daughter->Eta(),myEt);
695                       //cout<<"K0S daughter is a "<<daughter->GetName()<<endl;
696                     }
697                   }
698                   else{
699                     //cout<<"K0S daughter is a "<<daughter->GetName()<<endl;
700                   }
701                 }
702               }
703               filled = true;
704             }
705             if(pdgCode == fgK0LCode){
706               float myEt = Et(part);
707               fSimHadEt += myEt;
708               fSimTotEt += myEt;
709               FillHisto2D("EtSimulatedK0L",part->Pt(),part->Eta(),myEt);
710               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
711               filled = true;
712             }
713             if(pdgCode == fgOmegaCode){
714               float myEt = Et(part);
715               fSimHadEt += myEt;
716               fSimTotEt += myEt;
717               FillHisto2D("EtSimulatedOmega",part->Pt(),part->Eta(),myEt);
718               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
719               Int_t ndaughters = part->GetNDaughters();
720               for(Int_t idaughter = 0;idaughter<ndaughters;idaughter++){
721                 Int_t daughterindex = part->GetDaughter(idaughter);
722                 if(daughterindex<0 || daughterindex>1e5) continue;
723                 TParticle *daughter = stack->Particle(daughterindex);
724                 if(daughter){
725                   if(daughter->GetPDG(0)){
726
727                     Int_t daughtercode = daughter->GetPDG(0)->PdgCode();
728                     if(daughtercode==fgPiPlusCode || daughtercode==fgProtonCode || daughtercode==fgKMinusCode){
729                       myEt = Et(daughter);
730                       FillHisto2D("EtSimulatedOmegaDaughters",daughter->Pt(),daughter->Eta(),myEt);
731                     //cout<<"Omega daughter is a "<<daughter->GetName()<<endl;
732                     }
733                   }
734                   else{
735                     //cout<<"Omega daughter is a "<<daughter->GetName()<<endl;
736                   }
737                 }
738               }
739               filled = true;
740             }
741             if(pdgCode == fgAntiOmegaCode){
742               float myEt = Et(part);
743               fSimHadEt += myEt;
744               fSimTotEt += myEt;
745               FillHisto2D("EtSimulatedOmega",part->Pt(),part->Eta(),myEt);
746               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
747               Int_t ndaughters = part->GetNDaughters();
748               for(Int_t idaughter = 0;idaughter<ndaughters;idaughter++){
749                 Int_t daughterindex = part->GetDaughter(idaughter);
750                 if(daughterindex<0 || daughterindex>1e5) continue;
751                 TParticle *daughter = stack->ParticleFromTreeK(daughterindex);
752                 if(daughter){
753                   if(daughter->GetPDG(0)){
754                     Int_t daughtercode = daughter->GetPDG(0)->PdgCode();
755                     if(daughtercode==fgPiMinusCode || daughtercode==fgAntiProtonCode || daughtercode==fgKPlusCode){
756                       myEt = Et(daughter);
757                       FillHisto2D("EtSimulatedAntiOmegaDaughters",daughter->Pt(),daughter->Eta(),myEt);
758                       //cout<<"AntiOmega daughter is a "<<daughter->GetName()<<endl;
759                     }
760                   }
761                   else{
762                     //cout<<"AntiOmega daughter is a "<<daughter->GetName()<<endl;
763                   }
764                 }
765               }
766               filled = true;
767             }
768             //There are two codes for Sigmas
769             if(pdgCode == fgSigmaCode || pdgCode == -3222){
770               float myEt = Et(part);
771               fSimHadEt += myEt;
772               fSimTotEt += myEt;
773               FillHisto2D("EtSimulatedSigma",part->Pt(),part->Eta(),myEt);
774               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
775               filled = true;
776             }
777             if(pdgCode == fgAntiSigmaCode || pdgCode == 3222){
778               float myEt = Et(part);
779               fSimHadEt += myEt;
780               fSimTotEt += myEt;
781               FillHisto2D("EtSimulatedAntiSigma",part->Pt(),part->Eta(),myEt);
782               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
783               filled = true;
784             }
785             if(pdgCode == fgXiCode){
786               float myEt = Et(part);
787               fSimHadEt += myEt;
788               fSimTotEt += myEt;
789               FillHisto2D("EtSimulatedXi",part->Pt(),part->Eta(),myEt);
790               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
791               Int_t ndaughters = part->GetNDaughters();
792               for(Int_t idaughter = 0;idaughter<ndaughters;idaughter++){
793                 Int_t daughterindex = part->GetDaughter(idaughter);
794                 if(daughterindex<0 || daughterindex>1e5 || daughterindex>1e5) continue;
795                 //cerr<<"Daughter index "<<daughterindex<<" npart "<<nPrim<<endl;
796                 TParticle *daughter = stack->ParticleFromTreeK(daughterindex);
797                 if(daughter){
798                   if(daughter->GetPDG(0)){
799
800                     Int_t daughtercode = daughter->GetPDG(0)->PdgCode();
801                     if(daughtercode==fgPiPlusCode || daughtercode==fgProtonCode || daughtercode==fgPiMinusCode){
802                       myEt = Et(daughter);
803                       FillHisto2D("EtSimulatedXiDaughters",daughter->Pt(),daughter->Eta(),myEt);
804                     //cout<<"Xi daughter is a "<<daughter->GetName()<<endl;
805                     }
806                   }
807                   else{
808                     //cout<<"Xi daughter is a "<<daughter->GetName()<<endl;
809                   }
810                 }
811               }
812               filled = true;
813             }
814             if(pdgCode == fgAntiXiCode){
815               float myEt = Et(part);
816               fSimHadEt += myEt;
817               fSimTotEt += myEt;
818               FillHisto2D("EtSimulatedAntiXi",part->Pt(),part->Eta(),myEt);
819               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
820               Int_t ndaughters = part->GetNDaughters();
821               for(Int_t idaughter = 0;idaughter<ndaughters;idaughter++){
822                 Int_t daughterindex = part->GetDaughter(idaughter);
823                 if(daughterindex<0 || daughterindex>1e5) continue;
824                 TParticle *daughter = stack->ParticleFromTreeK(daughterindex);
825                 if(daughter){
826                   if(daughter->GetPDG(0)){
827                     Int_t daughtercode = daughter->GetPDG(0)->PdgCode();
828                     if(daughtercode==fgPiPlusCode || daughtercode==fgAntiProtonCode || daughtercode==fgPiMinusCode){
829                       myEt = Et(daughter);
830                       FillHisto2D("EtSimulatedAntiXiDaughters",daughter->Pt(),daughter->Eta(),myEt);
831                       //cout<<"AntiXi daughter is a "<<daughter->GetName()<<endl;
832                     }
833                   }
834                   else{
835                     //cout<<"AntiXi daughter is a "<<daughter->GetName()<<endl;
836                   }
837                 }
838               }
839               filled = true;
840             }
841             if(pdgCode == fgXi0Code){
842               float myEt = Et(part);
843               fSimHadEt += myEt;
844               fSimTotEt += myEt;
845               FillHisto2D("EtSimulatedXi0",part->Pt(),part->Eta(),myEt);
846               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
847               filled = true;
848             }
849             if(pdgCode == fgAntiXi0Code){
850               float myEt = Et(part);
851               fSimHadEt += myEt;
852               fSimTotEt += myEt;
853               FillHisto2D("EtSimulatedAntiXi0",part->Pt(),part->Eta(),myEt);
854               FillHisto2D("EtSimulatedAllHadron",part->Pt(),part->Eta(),myEt);
855               filled = true;
856             }
857             //============electrons===================================
858
859             if(pdgCode == fgEPlusCode){
860               float myEt = Et(part);
861               fSimTotEt += myEt;
862               FillHisto2D("EtSimulatedEPlus",part->Pt(),part->Eta(),myEt);
863               filled = true;
864             }
865             if(pdgCode == fgEMinusCode){
866               float myEt = Et(part);
867               fSimTotEt += myEt;
868               FillHisto2D("EtSimulatedEMinus",part->Pt(),part->Eta(),myEt);
869               filled = true;
870             }
871             //============neutrals===================================
872             if(pdgCode == fgGammaCode){
873               TParticle *mom = stack->Particle(part->GetFirstMother());
874               Int_t pdgCodeMom =  mom->GetPDG(0)->PdgCode();
875               //cout<<"I am a gamma and my mom is "<<mom->GetName()<<endl;
876               //We want to separate the gammas by pi0, eta, omega0 but we don't want to double count energy so we get the et from the gamma daughter
877               if(pdgCodeMom == fgEtaCode){
878                 float myEt = Et(part);
879                 fSimTotEt += myEt;
880                 FillHisto2D("EtSimulatedEta",mom->Pt(),mom->Eta(),myEt);
881                 filled = true;
882               }
883               if(pdgCodeMom == fgPi0Code){
884                 float myEt = Et(part);
885                 fSimTotEt += myEt;
886                 FillHisto2D("EtSimulatedPi0",mom->Pt(),mom->Eta(),myEt);
887                 filled = true;
888               }
889               if(pdgCodeMom == fgOmega0Code){
890                 float myEt = Et(part);
891                 fSimTotEt += myEt;
892                 FillHisto2D("EtSimulatedOmega0",mom->Pt(),mom->Eta(),myEt);
893                 filled = true;
894               }
895               if(!filled){
896                 float myEt = Et(part);
897                 fSimTotEt += myEt;
898                 FillHisto2D("EtSimulatedGamma",part->Pt(),part->Eta(),myEt);
899                 filled = true;
900               }
901             }
902             if(pdgCode == fgEtaCode){
903               float myEt = Et(part);
904               fSimTotEt += myEt;
905               FillHisto2D("EtSimulatedEta",part->Pt(),part->Eta(),myEt);
906               filled = true;
907             }
908             if(pdgCode == fgPi0Code){
909               float myEt = Et(part);
910               fSimTotEt += myEt;
911               FillHisto2D("EtSimulatedPi0",part->Pt(),part->Eta(),myEt);
912               filled = true;
913             }
914             if(pdgCode == fgOmega0Code){
915               float myEt = Et(part);
916               fSimTotEt += myEt;
917               FillHisto2D("EtSimulatedOmega0",part->Pt(),part->Eta(),myEt);
918               filled = true;
919             }
920             if(!filled){
921               //if( strcmp(pc->ParticleClass(),"Baryon")==0 || strcmp(pc->ParticleClass(),"Meson")==0 ){
922                 cout<<"Did not find a place for "<<part->GetName()<<" "<<pdgCode<<" which is a "<<pc->ParticleClass()<<endl;
923                 //}
924             }
925           }
926         }
927     }
928
929     if(fSimTotEt>0.0)FillHisto1D("SimTotEt",fSimTotEt,1.0);
930     if(fSimHadEt>0.0)FillHisto1D("SimHadEt",fSimHadEt,1.0);
931     if(fSimPiKPEt>0.0)FillHisto1D("SimPiKPEt",fSimPiKPEt,1.0);
932
933     //Smearing histograms
934     if(fSimPiKPEt>0.0) FillHisto2D("SimPiKPEtMinusSimPtSmeared",fSimPiKPEt,(fSimPiKPEt-fSimPiKPEtPtSmeared)/fSimPiKPEt,1.0);
935     FillHisto1D("SimPiKPEtPtSmeared",fSimPiKPEtPtSmeared,1.0);
936     if(fSimPiKPEt>0.0) FillHisto2D("SimPiKPEtMinusSimEfficiencySmeared",fSimPiKPEt,(fSimPiKPEt-fSimPiKPEtEfficiencySmeared)/fSimPiKPEt,1.0);
937     FillHisto1D("SimPiKPEtEfficiencySmeared",fSimPiKPEtEfficiencySmeared,1.0);
938     if(fSimPiKPEt>0.0) FillHisto2D("SimPiKPEtMinusSimPtCutSmearedTPC",fSimPiKPEt,(fSimPiKPEt-fSimPiKPEtPtCutSmearedTPC)/fSimPiKPEt,1.0);
939     FillHisto1D("SimPiKPEtPtCutSmearedTPC",fSimPiKPEtPtCutSmearedTPC,1.0);
940     if(fSimPiKPEt>0.0) FillHisto2D("SimPiKPEtMinusSimPtCutSmearedITS",fSimPiKPEt,(fSimPiKPEt-fSimPiKPEtPtCutSmearedITS)/fSimPiKPEt,1.0);
941     FillHisto1D("SimPiKPEtPtCutSmearedITS",fSimPiKPEtPtCutSmearedTPC,1.0);
942     if(fSimPiKPEt>0.0) FillHisto2D("SimPiKPEtMinusSimPIDSmeared",fSimPiKPEt,(fSimPiKPEt-fSimPiKPEtPIDSmeared)/fSimPiKPEt,1.0);
943     //if(fSimPiKPEt>0.0)cout<<"Filling SimPiKPEtMinusSimPIDSmeared with "<<fSimPiKPEt<<","<<(fSimPiKPEt-fSimPiKPEtPIDSmeared)/fSimPiKPEt<<endl;
944     FillHisto1D("SimPiKPEtPIDSmeared",fSimPiKPEtPIDSmeared,1.0);
945     if(fSimPiKPEt>0.0) FillHisto2D("SimPiKPEtMinusSimPIDSmearedNoID",fSimPiKPEt,(fSimPiKPEt-fSimPiKPEtPIDSmearedNoID)/fSimPiKPEt,1.0);
946     //if(fSimPiKPEt>0.0)cout<<"Filling SimPiKPEtMinusSimPIDSmearedNoID with "<<fSimPiKPEt<<","<<(fSimPiKPEt-fSimPiKPEtPIDSmearedNoID)/fSimPiKPEt<<endl;
947     FillHisto1D("SimPiKPEtPIDSmearedNoID",fSimPiKPEtPIDSmearedNoID,1.0);
948
949     return 1;
950     
951 }
952
953 void AliAnalysisHadEtMonteCarlo::Init()
954 { // Init
955     AliAnalysisHadEt::Init();
956     if(!fPtSmearer) fPtSmearer = new TRandom();
957 }
958
959 void AliAnalysisHadEtMonteCarlo::CreateHistograms(){
960   //for simulated Et only (no reconstruction)
961   CreateEtaPtHisto2D(TString("EtSimulatedPiPlus"),TString("Simulated E_{T} from #pi^{+}"));
962   CreateEtaPtHisto2D("EtSimulatedPiMinus","Simulated E_{T} from #pi^{-}");
963   CreateEtaPtHisto2D("EtSimulatedKPlus","Simulated E_{T} from K^{+}");
964   CreateEtaPtHisto2D("EtSimulatedKMinus","Simulated E_{T} from K^{-}");
965   CreateEtaPtHisto2D("EtSimulatedProton","Simulated E_{T} from p");
966   CreateEtaPtHisto2D("EtSimulatedAntiProton","Simulated E_{T} from #bar{p}");
967   CreateEtaPtHisto2D("EtSimulatedChargedHadron","Simulated E_{T} from charged hadrons");
968   CreateEtaPtHisto2D(TString("EtNSimulatedPiPlus"),TString("Number of Simulated #pi^{+}"));
969   CreateEtaPtHisto2D("EtNSimulatedPiMinus","Number of simulated #pi^{-}");
970   CreateEtaPtHisto2D("EtNSimulatedKPlus","Number of simulated K^{+}");
971   CreateEtaPtHisto2D("EtNSimulatedKMinus","Number of simulated K^{-}");
972   CreateEtaPtHisto2D("EtNSimulatedProton","Number of simulated p");
973   CreateEtaPtHisto2D("EtNSimulatedAntiProton","Number of simulated #bar{p}");
974   CreateEtaPtHisto2D("EtNSimulatedChargedHadron","Number of simulated charged hadrons");
975   CreateEtaPtHisto2D("EtSimulatedChargedHadronAssumingNoPt","Simulated E_{T} from charged hadrons assuming p_{T}=0");
976   CreateEtaPtHisto2D("EtSimulatedChargedHadronAssumingPtTPCCut","Simulated E_{T} from charged hadrons assuming p_{T}=0.15");
977   CreateEtaPtHisto2D("EtSimulatedChargedHadronAssumingPtITSCut","Simulated E_{T} from charged hadrons assuming p_{T}=0.10");
978
979   CreateEtaPtHisto2D("EtSimulatedChargedHadronAssumingPion","Simulated E_{T} from charged hadrons assuming they are all pions");
980   CreateEtaPtHisto2D("EtSimulatedKPlusAssumingPion","Simulated E_{T} from K^{+} assuming #pi mass");
981   CreateEtaPtHisto2D("EtSimulatedKMinusAssumingPion","Simulated E_{T} from K^{-} assuming #pi mass");
982   CreateEtaPtHisto2D("EtSimulatedProtonAssumingPion","Simulated E_{T} from p assuming #pi mass");
983   CreateEtaPtHisto2D("EtSimulatedAntiProtonAssumingPion","Simulated E_{T} from #bar{p} assuming #pi mass");
984
985   CreateEtaPtHisto2D("EtSimulatedLambda","Simulated E_{T} from #Lambda");
986   CreateEtaPtHisto2D("EtSimulatedAntiLambda","Simulated E_{T} from #bar{#Lambda}");
987   CreateEtaPtHisto2D("EtSimulatedK0S","Simulated E_{T} from K^{0}_{S}");
988   CreateEtaPtHisto2D("EtSimulatedK0L","Simulated E_{T} from K^{0}_{L}");
989   CreateEtaPtHisto2D("EtSimulatedNeutron","Simulated E_{T} from neutrons");
990   CreateEtaPtHisto2D("EtSimulatedAntiNeutron","Simulated E_{T} from #bar{n}");
991   CreateEtaPtHisto2D("EtSimulatedEPlus","Simulated E_{T} from e^{+}");
992   CreateEtaPtHisto2D("EtSimulatedEMinus","Simulated E_{T} from e^{-}");
993   CreateEtaPtHisto2D("EtSimulatedOmega","Simulated E_{T} from #Omega^{-}");
994   CreateEtaPtHisto2D("EtSimulatedAntiOmega","Simulated E_{T} from #Omega^{+}");
995   CreateEtaPtHisto2D("EtSimulatedXi","Simulated E_{T} from #Xi^{-}");
996   CreateEtaPtHisto2D("EtSimulatedAntiXi","Simulated E_{T} from #Xi^{+}");
997   CreateEtaPtHisto2D("EtSimulatedSigma","Simulated E_{T} from #Xi^{-}");
998   CreateEtaPtHisto2D("EtSimulatedAntiSigma","Simulated E_{T} from #Xi^{+}");
999   CreateEtaPtHisto2D("EtSimulatedXi0","Simulated E_{T} from #Xi^{0}");
1000   CreateEtaPtHisto2D("EtSimulatedAntiXi0","Simulated E_{T} from #Xi^{0}");
1001   CreateEtaPtHisto2D("EtSimulatedAllHadron","Simulated E_{T} from all hadrons");
1002
1003
1004   CreateEtaPtHisto2D("EtSimulatedLambdaDaughters","Simulated E_{T} from #Lambda Daughters");
1005   CreateEtaPtHisto2D("EtSimulatedAntiLambdaDaughters","Simulated E_{T} from #bar{#Lambda} Daughters");
1006   CreateEtaPtHisto2D("EtSimulatedK0SDaughters","Simulated E_{T} from K^{0}_{S} Daughters");
1007   CreateEtaPtHisto2D("EtSimulatedOmegaDaughters","Simulated E_{T} from #Omega^{-} Daughters");
1008   CreateEtaPtHisto2D("EtSimulatedAntiOmegaDaughters","Simulated E_{T} from #Omega^{+} Daughters");
1009   CreateEtaPtHisto2D("EtSimulatedXiDaughters","Simulated E_{T} from #Xi^{-} Daughters");
1010   CreateEtaPtHisto2D("EtSimulatedAntiXiDaughters","Simulated E_{T} from #Xi^{+} Daughters");
1011
1012
1013   CreateEtaPtHisto2D("EtSimulatedGamma","Simulated E_{T} from #gamma");
1014   CreateEtaPtHisto2D("EtSimulatedEta","Simulated E_{T} from #eta");
1015   CreateEtaPtHisto2D("EtSimulatedPi0","Simulated E_{T} from #pi^{0}");
1016   CreateEtaPtHisto2D("EtSimulatedOmega0","Simulated E_{T} from #omega");
1017
1018   TString *strTPC = new TString("TPC");
1019   TString *strITS = new TString("ITS");
1020   TString *strTPCITS = new TString("TPCITS");
1021   for(Int_t i=0;i<2;i++){
1022     TString *cutName;
1023     Float_t maxPtdEdx = 10;
1024     Float_t mindEdx = 35;
1025     Float_t maxdEdx = 150.0;
1026     switch(i){
1027     case 0:
1028       cutName = strTPC;
1029       break;
1030     case 1:
1031       cutName = strITS;
1032       maxPtdEdx = 5;
1033       maxdEdx = 500.0;
1034       break;
1035     case 2:
1036       cutName = strTPCITS;
1037       break;
1038     default:
1039       cerr<<"Error:  cannot make histograms!"<<endl;
1040       return;
1041     }
1042
1043     CreateEtaPtHisto2D(Form("EtReconstructed%sIdentifiedPiPlus",cutName->Data()),"Reconstructed E_{T} from identified #pi^{+}");
1044     CreateEtaPtHisto2D(Form("EtReconstructed%sIdentifiedPiMinus",cutName->Data()),"Reconstructed E_{T} from identified #pi^{-}");
1045     CreateEtaPtHisto2D(Form("EtReconstructed%sIdentifiedKPlus",cutName->Data()),"Reconstructed E_{T} from identified K^{+}");
1046     CreateEtaPtHisto2D(Form("EtReconstructed%sIdentifiedEMinus",cutName->Data()),"Reconstructed E_{T} from identified e^{-}");
1047     CreateEtaPtHisto2D(Form("EtReconstructed%sIdentifiedEPlus",cutName->Data()),"Reconstructed E_{T} from identified e^{+}");
1048     CreateEtaPtHisto2D(Form("EtReconstructed%sIdentifiedKMinus",cutName->Data()),"Reconstructed E_{T} from identified K^{-}");
1049     CreateEtaPtHisto2D(Form("EtReconstructed%sIdentifiedProton",cutName->Data()),"Reconstructed E_{T} from identified p");
1050     CreateEtaPtHisto2D(Form("EtReconstructed%sIdentifiedAntiProton",cutName->Data()),"Reconstructed E_{T} from identified #bar{p}");
1051     CreateEtaPtHisto2D(Form("EtNReconstructed%sUnidentified",cutName->Data()),"Number of Reconstructed unidentified particles");
1052     CreateEtaPtHisto2D(Form("EtReconstructed%sUnidentifiedAssumingPion",cutName->Data()),"Reconstructed E_{T} from unidentified particles assuming pion mass");
1053     CreateEtaPtHisto2D(Form("EtReconstructed%sUnidentified",cutName->Data()),"Reconstructed E_{T} from unidentified particles using real mass");
1054     CreateEtaPtHisto2D(Form("EtReconstructed%sMisidentifiedElectrons",cutName->Data()),"Reconstructed E_{T} from misidentified electrons");
1055
1056
1057     CreateEtaPtHisto2D(Form("EtReconstructed%sPiPlus",cutName->Data()),"Reconstructed E_{T} from #pi^{+}");
1058     CreateEtaPtHisto2D(Form("EtReconstructed%sPiMinus",cutName->Data()),"Reconstructed E_{T} from #pi^{-}");
1059     CreateEtaPtHisto2D(Form("EtReconstructed%sKPlus",cutName->Data()),"Reconstructed E_{T} from K^{+}");
1060     CreateEtaPtHisto2D(Form("EtReconstructed%sKMinus",cutName->Data()),"Reconstructed E_{T} from K^{-}");
1061     CreateEtaPtHisto2D(Form("EtReconstructed%sProton",cutName->Data()),"Reconstructed E_{T} from p");
1062     CreateEtaPtHisto2D(Form("EtReconstructed%sAntiProton",cutName->Data()),"Reconstructed E_{T} from #bar{p}");
1063     CreateEtaPtHisto2D(Form("EtReconstructed%sChargedHadron",cutName->Data()),"Reconstructed E_{T} from charged hadrons");
1064     CreateEtaPtHisto2D(Form("EtNReconstructed%sPiPlus",cutName->Data()),"Reconstructed E_{T} from #pi^{+}");
1065     CreateEtaPtHisto2D(Form("EtNReconstructed%sPiMinus",cutName->Data()),"Reconstructed E_{T} from #pi^{-}");
1066     CreateEtaPtHisto2D(Form("EtNReconstructed%sKPlus",cutName->Data()),"Reconstructed E_{T} from K^{+}");
1067     CreateEtaPtHisto2D(Form("EtNReconstructed%sKMinus",cutName->Data()),"Reconstructed E_{T} from K^{-}");
1068     CreateEtaPtHisto2D(Form("EtNReconstructed%sProton",cutName->Data()),"Reconstructed E_{T} from p");
1069     CreateEtaPtHisto2D(Form("EtNReconstructed%sAntiProton",cutName->Data()),"Reconstructed E_{T} from #bar{p}");
1070     CreateEtaPtHisto2D(Form("EtNReconstructed%sChargedHadron",cutName->Data()),"Reconstructed E_{T} from charged hadrons");
1071
1072     CreateEtaPtHisto2D(Form("EtReconstructed%sChargedHadronAssumingPion",cutName->Data()),"Reconstructed E_{T} from charged hadrons assuming they are all pions");
1073     CreateEtaPtHisto2D(Form("EtReconstructed%sKPlusAssumingPion",cutName->Data()),"Reconstructed E_{T} from K^{+} assuming #pi mass");
1074     CreateEtaPtHisto2D(Form("EtReconstructed%sKMinusAssumingPion",cutName->Data()),"Reconstructed E_{T} from K^{-} assuming #pi mass");
1075     CreateEtaPtHisto2D(Form("EtReconstructed%sProtonAssumingPion",cutName->Data()),"Reconstructed E_{T} from p assuming #pi mass");
1076     CreateEtaPtHisto2D(Form("EtReconstructed%sAntiProtonAssumingPion",cutName->Data()),"Reconstructed E_{T} from #bar{p} assuming #pi mass");
1077
1078     CreateEtaPtHisto2D(Form("EtReconstructed%sEPlus",cutName->Data()),"Reconstructed E_{T} from e^{+}");
1079     CreateEtaPtHisto2D(Form("EtReconstructed%sEMinus",cutName->Data()),"Reconstructed E_{T} from e^{-}");
1080
1081
1082
1083     CreateEtaPtHisto2D(Form("EtReconstructed%sLambdaDaughters",cutName->Data()),"Reconstructed E_{T} from #Lambda Daughters");
1084   CreateEtaPtHisto2D(Form("EtReconstructed%sAntiLambdaDaughters",cutName->Data()),"Reconstructed E_{T} from #bar{#Lambda} Daughters");
1085   CreateEtaPtHisto2D(Form("EtReconstructed%sK0SDaughters",cutName->Data()),"Reconstructed E_{T} from K^{0}_{S} Daughters");
1086   CreateEtaPtHisto2D(Form("EtReconstructed%sOmegaDaughters",cutName->Data()),"Reconstructed E_{T} from #Omega^{-} Daughters");
1087   CreateEtaPtHisto2D(Form("EtReconstructed%sAntiOmegaDaughters",cutName->Data()),"Reconstructed E_{T} from #Omega^{+} Daughters");
1088   CreateEtaPtHisto2D(Form("EtReconstructed%sXiDaughters",cutName->Data()),"Reconstructed E_{T} from #Xi^{-} Daughters");
1089   CreateEtaPtHisto2D(Form("EtReconstructed%sAntiXiDaughters",cutName->Data()),"Reconstructed E_{T} from #Xi^{+} Daughters");
1090
1091     CreateIntHisto1D(Form("UnidentifiedPIDs%s",cutName->Data()),"PIDs of unidentified particles", "PID", "Number of particles",9, -4,4);
1092     CreateHisto2D(Form("MisidentifiedPIDs%s",cutName->Data()),"PIDs of misidentified particles", "PID real","PID identified",5, -.5,4.5,5, -.5,4.5);
1093     CreateHisto2D(Form("dEdxAll%s",cutName->Data()),"dE/dx for all particles","momentum (GeV/c)","dE/dx",400,0.0,maxPtdEdx,200,mindEdx,maxdEdx);
1094     CreateHisto2D(Form("dEdxPion%s",cutName->Data()),"dE/dx for #pi^{#pm}","momentum (GeV/c)","dE/dx",400,0.0,maxPtdEdx,200,mindEdx,maxdEdx);
1095     CreateHisto2D(Form("dEdxKaon%s",cutName->Data()),"dE/dx for K^{#pm}","momentum (GeV/c)","dE/dx",400,0.0,maxPtdEdx,200,mindEdx,maxdEdx);
1096     CreateHisto2D(Form("dEdxProton%s",cutName->Data()),"dE/dx for p(#bar{p})","momentum (GeV/c)","dE/dx",400,0.0,maxPtdEdx,200,mindEdx,maxdEdx);
1097     CreateHisto2D(Form("dEdxElectron%s",cutName->Data()),"dE/dx for e^{#pm}","momentum (GeV/c)","dE/dx",400,0.0,maxPtdEdx,200,mindEdx,maxdEdx);
1098     CreateHisto2D(Form("dEdxUnidentified%s",cutName->Data()),"dE/dx for unidentified particles","momentum (GeV/c)","dE/dx",400,0.0,maxPtdEdx,200,mindEdx,maxdEdx);
1099   }
1100   delete strTPC;
1101   delete strITS;
1102   delete strTPCITS;
1103
1104   Float_t minEt = 0.0;
1105   Float_t maxEt = 100.0;
1106   Int_t nbinsEt = 100;
1107   char histoname[200];
1108   char histotitle[200];
1109   char xtitle[50];
1110   char ytitle[50];
1111   TString *sTPC = new TString("TPC");
1112   TString *sITS = new TString("ITS");
1113   TString *sTPCpt = new TString("0.15");
1114   TString *sITSpt = new TString("0.10");
1115   TString *sPID = new TString("");
1116   TString *sNoPID = new TString("NoPID");
1117   TString *sNoPIDString = new TString(", No PID");
1118   TString *sHadEt = new TString("HadEt");
1119   TString *sTotEt = new TString("TotEt");
1120   TString *sTotEtString = new TString("total E_{T}");
1121   TString *sHadEtString = new TString("hadronic E_{T}");
1122   TString *sFull = new TString("Full");
1123   TString *sEMCAL = new TString("EMCAL");
1124   TString *sPHOS = new TString("PHOS");
1125   float etDiff = 1.5;
1126   
1127   for(int tpc = 0;tpc<2;tpc++){
1128     TString *detector;
1129     TString *ptstring;
1130     if(tpc==1) {detector = sTPC; ptstring = sTPCpt;}
1131     else{detector = sITS; ptstring = sITSpt;}
1132     for(int hadet = 0;hadet<2;hadet++){
1133       TString *et;
1134       TString *etstring;
1135       if(hadet==1) {et = sHadEt; etstring = sHadEtString;}
1136       else{et = sTotEt; etstring = sTotEtString;}
1137       for(int type = 0;type<1;type++){
1138         TString *acceptance;
1139         switch(type){
1140         case 0:
1141           acceptance = sFull;
1142           break;
1143         case 1:
1144           acceptance = sEMCAL;
1145           break;
1146         case 2:
1147           acceptance = sPHOS;
1148           break;
1149         default:
1150           acceptance = sFull;
1151         }
1152         sprintf(histoname,"Sim%sMinusRawEt%sAcceptance%s",et->Data(),acceptance->Data(),detector->Data());
1153         sprintf(histotitle,"(Simulated %s - raw reconstructed)/(Simulated %s) with %s acceptance for p_{T}>%s GeV/c",etstring->Data(),etstring->Data(),acceptance->Data(),ptstring->Data());
1154         sprintf(ytitle,"(Simulated %s - raw reconstructed)/(Simulated %s)",etstring->Data(),etstring->Data());
1155         sprintf(xtitle,"Simulated %s",etstring->Data());
1156         CreateHisto2D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt,nbinsEt,-etDiff,etDiff);
1157
1158         for(int pid = 0;pid<2;pid++){
1159           TString *partid;
1160           TString *partidstring;
1161           if(pid==1){partid = sPID; partidstring = sPID;}
1162           else{partid = sNoPID; partidstring = sNoPIDString;}
1163
1164           sprintf(histoname,"Sim%sVsReco%s%sAcceptance%s%s",et->Data(),et->Data(),acceptance->Data(),detector->Data(),partid->Data());
1165           sprintf(histotitle,"Simulated %s vs reconstructed %s with %s acceptance for p_{T}>%s GeV/c%s",etstring->Data(),etstring->Data(),acceptance->Data(),ptstring->Data(),partidstring->Data());
1166           sprintf(xtitle,"Simulated %s",etstring->Data());
1167           sprintf(ytitle,"Reconstructed %s (%s acc., p_{T}>%s GeV/c,%s)",etstring->Data(),acceptance->Data(),ptstring->Data(),partidstring->Data());
1168           CreateHisto2D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt,nbinsEt,minEt,maxEt);
1169
1170           sprintf(histoname,"Sim%sMinusReco%s%sAcceptance%s%s",et->Data(),et->Data(),acceptance->Data(),detector->Data(),partid->Data());
1171           sprintf(histotitle,"(Simulated %s - reconstructed %s)/(Simulated %s) with %s acceptance for p_{T}>%s GeV/c%s",etstring->Data(),etstring->Data(),etstring->Data(),acceptance->Data(),ptstring->Data(),partidstring->Data());
1172           sprintf(ytitle,"(Simulated %s - reconstructed %s)/(Simulated %s)",etstring->Data(),etstring->Data(),etstring->Data());
1173           sprintf(xtitle,"Simulated %s",etstring->Data());
1174           CreateHisto2D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt,nbinsEt,-etDiff,etDiff);
1175
1176           if(hadet==0){//we only want to do this once...  not the most elegant way of coding but hey...
1177             sprintf(histoname,"SimPiKPMinusRecoPiKP%sAcceptance%s%s",acceptance->Data(),detector->Data(),partid->Data());
1178             sprintf(histotitle,"(Sim PiKP - reco PiKP)/(Sim PiKP) with %s acceptance for p_{T}>%s GeV/c%s",acceptance->Data(),ptstring->Data(),partidstring->Data());
1179             sprintf(ytitle,"(Sim PiKP - reco PiKP)/(Sim PiKP)");
1180             sprintf(xtitle,"Simulated E_{T}^{#pi,K,p}");
1181             CreateHisto2D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt,nbinsEt,-etDiff,etDiff);
1182             //cout<<"Creating "<<histoname<<endl;
1183           }
1184         }
1185       }
1186     }
1187   }
1188    CreateHisto1D("SimPiKPEt","Simulated #pi,K,p E_{T}","Simulated #pi,K,p E_{T}","Number of events",nbinsEt,minEt,maxEt);
1189   CreateHisto1D("SimTotEt","Simulated Total E_{T}","Simulated Total E_{T}","Number of events",nbinsEt*4,minEt,maxEt);
1190   CreateHisto1D("SimHadEt","Simulated Hadronic E_{T}","Simulated Hadronic E_{T}","Number of events",nbinsEt*4,minEt,maxEt);
1191
1192   etDiff = 0.15;
1193
1194   //======================================================================
1195
1196   sprintf(histoname,"SimPiKPEtMinusSimPtSmeared");
1197   sprintf(histotitle,"Simulated (true-smeared)/true for 0.5 percent momentum smearing");
1198   sprintf(ytitle,"(true-smeared)/true");
1199   sprintf(xtitle,"true p, K, p E_{T}");
1200   CreateHisto2D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt,nbinsEt,-etDiff/10.0,etDiff/10.0);
1201   sprintf(histoname,"SimPiKPEtPtSmeared");
1202   sprintf(histotitle,"Simulated E_{T} for 0.5 percent momentum smearing");
1203   sprintf(ytitle,"Number of events");
1204   sprintf(xtitle,"p, K, p E_{T}");
1205   CreateHisto1D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt);
1206
1207   //======================================================================
1208
1209   sprintf(histoname,"SimPiKPEtMinusSimEfficiencySmeared");
1210   sprintf(histotitle,"Simulated (true-smeared)/true for efficiency smearing");
1211   sprintf(ytitle,"(true-smeared)/true");
1212   sprintf(xtitle,"true p, K, p E_{T}");
1213   CreateHisto2D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt,nbinsEt,-etDiff*5,etDiff*5);
1214   sprintf(histoname,"SimPiKPEtEfficiencySmeared");
1215   sprintf(histotitle,"Simulated E_{T} for efficiency smearing");
1216   sprintf(ytitle,"Number of events");
1217   sprintf(xtitle,"p, K, p E_{T}");
1218   CreateHisto1D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt);
1219
1220   //======================================================================
1221
1222   sprintf(histoname,"SimPiKPEtMinusSimPtCutSmearedTPC");
1223   sprintf(histotitle,"Simulated (true-smeared)/true for p_{T}>0.15 GeV/c smearing");
1224   sprintf(ytitle,"(true-smeared)/true");
1225   sprintf(xtitle,"true p, K, p E_{T}");
1226   CreateHisto2D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt,nbinsEt,-etDiff,etDiff);
1227   sprintf(histoname,"SimPiKPEtPtCutSmearedTPC");
1228   sprintf(histotitle,"Simulated E_{T} for p_{T}>0.15 GeV/c smearing");
1229   sprintf(ytitle,"Number of events");
1230   sprintf(xtitle,"p, K, p E_{T}");
1231   CreateHisto1D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt);
1232
1233
1234   //======================================================================
1235
1236   sprintf(histoname,"SimPiKPEtMinusSimPtCutSmearedITS");
1237   sprintf(histotitle,"Simulated (true-smeared)/true for p_{T}>0.10 GeV/c smearing");
1238   sprintf(ytitle,"(true-smeared)/true");
1239   sprintf(xtitle,"true p, K, p E_{T}");
1240   CreateHisto2D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt,nbinsEt,-etDiff,etDiff);
1241   sprintf(histoname,"SimPiKPEtPtCutSmearedITS");
1242   sprintf(histotitle,"Simulated E_{T} for p_{T}>0.10 GeV/c smearing");
1243   sprintf(ytitle,"Number of events");
1244   sprintf(xtitle,"p, K, p E_{T}");
1245   CreateHisto1D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt);
1246
1247   //======================================================================
1248
1249   sprintf(histoname,"SimPiKPEtMinusSimPIDSmeared");
1250   sprintf(histotitle,"Simulated (true-smeared)/true for PID smearing");
1251   sprintf(ytitle,"(true-smeared)/true");
1252   sprintf(xtitle,"true p, K, p E_{T}");
1253   CreateHisto2D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt,nbinsEt,-etDiff,etDiff);
1254   sprintf(histoname,"SimPiKPEtPIDSmeared");
1255   sprintf(histotitle,"Simulated E_{T} for PID smearing");
1256   sprintf(ytitle,"Number of events");
1257   sprintf(xtitle,"p, K, p E_{T}");
1258   CreateHisto1D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt);
1259
1260   //======================================================================
1261
1262   sprintf(histoname,"SimPiKPEtMinusSimPIDSmearedNoID");
1263   sprintf(histotitle,"Simulated (true-smeared)/true for PID smearing No ID");
1264   sprintf(ytitle,"(true-smeared)/true");
1265   sprintf(xtitle,"true p, K, p E_{T}");
1266   CreateHisto2D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt,nbinsEt,-etDiff,etDiff);
1267   sprintf(histoname,"SimPiKPEtPIDSmearedNoID");
1268   sprintf(histotitle,"Simulated E_{T} for PID smearing No ID");
1269   sprintf(ytitle,"Number of events");
1270   sprintf(xtitle,"p, K, p E_{T}");
1271   CreateHisto1D(histoname,histotitle,xtitle,ytitle,nbinsEt,minEt,maxEt);
1272
1273   delete sTPC;
1274   delete sITS;
1275   delete sTPCpt;
1276   delete sITSpt;
1277   delete sPID;
1278   delete sNoPID;
1279   delete sNoPIDString;
1280   delete sHadEt;
1281   delete sTotEt;
1282   delete sTotEtString;
1283   delete sHadEtString;
1284   delete sFull;
1285   delete sEMCAL;
1286   delete sPHOS;
1287   CreateIntHisto1D("NEvents","Number of events","number of events","Number of events",1,0,1);
1288
1289 }
1290