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