]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGDQ/dielectron/macrosJPSI/ConfigJpsi_mw_pPb_MC.C
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGDQ / dielectron / macrosJPSI / ConfigJpsi_mw_pPb_MC.C
CommitLineData
53cbd1e0 1void SetupTrackCutsDieleData(AliDielectron *diele, Int_t cutDefinition);
2void SetupPairCutsDieleData(AliDielectron *diele, Int_t cutDefinition);
3void InitHistogramsDieleData(AliDielectron *diele, Int_t cutDefinition);
4void InitCFDieleData(AliDielectron *diele, Int_t cutDefinition);
5void SetupV0Cuts(AliDielectron *diele, Int_t cutDefinition);
6void SetSignals(AliDielectron *diele);
7
8
9TString namesDieleData=("basicQ+SPDfirst+pt>1+PID;basicQ+SPDany+pt>1+PID;basicQ+SPDany+pt>1+noexclPIDforcontrolpurpose;nocuts");
10//has to introduce PID cuts....
11//TString namesDieleData=("basicQ+SPDfirst+pt>1+PID");
12
13TObjArray *arrNamesDieleData=namesDieleData.Tokenize(";");
14
15const Int_t nDie=arrNamesDieleData->GetEntries();
16
17AliDielectron* ConfigJpsi_mw_pPb_MC(Int_t cutDefinition)
18{
19 //
20 // Setup the instance of AliDielectron
21 //
22
23 // create the actual framework object
24 TString name=Form("%02d",cutDefinition);
25 if (cutDefinition<arrNamesDieleData->GetEntriesFast()){
26 name=arrNamesDieleData->At(cutDefinition)->GetName();
27 }
28 AliDielectron *diele = new AliDielectron(Form("%s",name.Data()),
29 Form("Track cuts: %s",name.Data()));
30
31
32 // estimators filename
33 //NOTE: what does this mean?: estimator for pp multiplicity, not needed for instance for my pA-purpose(mwinn 16.1.2012)..
34 diele->SetEstimatorFilename("$ALICE_ROOT/PWGDQ/dielectron/files/estimators.root");
35 //diele->SetEstimatorFilename("estimators.root");
36 // cut setup
37 SetupTrackCutsDieleData(diele, cutDefinition);
38 SetupPairCutsDieleData(diele, cutDefinition);
39 SetupV0Cuts(diele, cutDefinition);
40
41 // Set MC signals
42 SetSignals(diele);
43 //
44 // histogram setup
45 // only if an AliDielectronHistos object is attached to the
46 // dielelectron framework histograms will be filled
47 //
48 InitHistogramsDieleData(diele, cutDefinition);
49
50 // the last definition uses no cuts and only the QA histograms should be filled!, now for all cuts
51 if(cutDefinition == 1 ||cutDefinition == 2 || cutDefinition == 3 || cutDefinition == 4) InitCFDieleData(diele, cutDefinition);//first 2 cut sets in 3rd included
52
53
54 return diele;
55}
56
57//______________________________________________________________________________________
58void SetupTrackCutsDieleData(AliDielectron *diele, Int_t cutDefinition)
59{
60 //
61 // Setup the track cuts
62 //
63
64 //NOTE: seems to work, see AliDielectronTrackCuts method IsSelected at the beginning, to be checked with AODs
65 //exclude conversion electrons selected by the tender
66 /* if(!cutDefinition==4){
67 AliDielectronTrackCuts *noconv=new AliDielectronTrackCuts("noConv","noConv");
68 noconv->SetV0DaughterCut(AliPID::kElectron,kTRUE);
69 diele->GetTrackFilter().AddCuts(noconv);
70 }*/
71 // }
72 AliDielectronTrackCuts *trackCuts=new AliDielectronTrackCuts("ITSandgeneral_trackCuts","ITSandgeneral_trackCuts");
73 //ITS related cuts
74 if (cutDefinition==0)
75 trackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,AliESDtrackCuts::kFirst); //does not work in MC???
76 trackCuts->SetRequireTPCRefit(kTRUE);
77 trackCuts->SetRequireITSRefit(kTRUE);
78 diele->GetTrackFilter().AddCuts(trackCuts);
79
80 //Pt cut ----------------------------------------------------------
81 AliDielectronVarCuts *pt = new AliDielectronVarCuts("trackkineandTPCQ","trackkine_and_TPCQ");
82 if ((cutDefinition==0)){
83 pt->AddCut(AliDielectronVarManager::kPt,1.,1e30);
84 }else{ pt->AddCut(AliDielectronVarManager::kP,.8,1e30);}
85 pt->AddCut(AliDielectronVarManager::kKinkIndex0,0.);
86 // SPDany
87 if (cutDefinition==2) pt->AddCut(AliDielectronVarManager::kITSLayerFirstCls,-0.5,1.5);
88 //AOD additions since there are no AliESDtrackCuts -----------------
89 //
90 // TPC #clusteres cut
91 if (!(cutDefinition==1) ) pt->AddCut(AliDielectronVarManager::kNclsTPC,70.,160.);//does not work in MC???
92 if (!(cutDefinition==4) )pt->AddCut(AliDielectronVarManager::kTPCchi2Cl,0.,4.);//to be checked
93 if (!(cutDefinition==1) ) pt->AddCut(AliDielectronVarManager::kEta,-0.9,0.9);//does not work in MC???
94 //TODO: DCA cuts to be investigated!!! NOTE: why?? (mwinn, 15.01.2013)
95 if (!(cutDefinition==4) ) pt->AddCut(AliDielectronVarManager::kImpactParXY,-1.,1.);
96 if (!(cutDefinition==4) ) pt->AddCut(AliDielectronVarManager::kImpactParZ,-3.,3.);
97
98 diele->GetTrackFilter().AddCuts(pt);
99
100 // PID cuts --------------------------------------------------------
101 if(cutDefinition ==0 || cutDefinition ==1 ||cutDefinition ==2){
102 AliDielectronPID *pid = new AliDielectronPID("PID10","TPC nSigma |e|<3 + |Pi|>3.5 + P>3");
103 pid->AddCut(AliDielectronPID::kTPC,AliPID::kElectron,-3.,3.);
104 pid->AddCut(AliDielectronPID::kTPC,AliPID::kPion,-3.5,3.5,0.,0.,kTRUE);//does not work in MC???
105 pid->AddCut(AliDielectronPID::kTPC,AliPID::kProton,-20.,3.,0.,0.,kTRUE);//does not work in MC???
106 diele->GetTrackFilter().AddCuts(pid);
107 }
108 if(cutDefinition ==3){
109 AliDielectronPID *pid = new AliDielectronPID("PID10","TPC nSigma |e|<3");
110 pid->AddCut(AliDielectronPID::kTPC,AliPID::kElectron,-3.,3.);
111 diele->GetTrackFilter().AddCuts(pid);
112 /* //taken out at the 04.02.2013 in order to have some feeling for PID cuts in QAplots...
113 AliDielectronVarCuts *pidsubs = new AliDielectronVarCuts("pidSubs","pidsubs cut");
114 pidsubs->AddCut(AliDielectronVarManager::kP,1.2,1e30);
115 pidsubs->AddCut(AliDielectronVarManager::kTPCsignal,70.,110.);
116 diele->GetTrackFilter().AddCuts(pidsubs);*/
117 }
118}
119
120//______________________________________________________________________________________
121void SetupPairCutsDieleData(AliDielectron *diele, Int_t cutDefinition)
122{
123 //
124 // Setup the pair cuts
125 //
126 // conversion rejection
127 //Double_t gCut = 0.05; // default
128 Double_t gCut = 0.100; // default
129 if(!cutDefinition==4){
130 AliDielectronVarCuts *gammaCut=new AliDielectronVarCuts("gammaCut","gammaCut");
131 gammaCut->AddCut(AliDielectronVarManager::kM,0.,gCut);
132 diele->GetPairPreFilter().AddCuts(gammaCut);
133 diele->SetPreFilterUnlikeOnly();
134 }
135
136 //Invariant mass and rapidity selection
137 if(!cutDefinition==4){
138 AliDielectronVarCuts *pairCut=new AliDielectronVarCuts("|Y|<.9","|Y|<.9");
139 // pairCut->AddCut(AliDielectronVarManager::kM,2.,4.);
140 pairCut->AddCut(AliDielectronVarManager::kY,-0.9,0.9);
141 diele->GetPairFilter().AddCuts(pairCut);
142 }
143}
144//______________________________________________________________________________________
145
146void SetupV0Cuts(AliDielectron *diele, Int_t cutDefinition)
147{
148 //
149 // Setup the V0 cuts
150 //
151 if(!cutDefinition==4){
152
153 AliDielectronV0Cuts *gammaV0Cuts = new AliDielectronV0Cuts("IsGamma","IsGamma");
154 gammaV0Cuts->SetPdgCodes(22,11,11);
155 gammaV0Cuts->SetDefaultPID(13); // TPC+-10, TOF+-3 if available
156 gammaV0Cuts->AddCut(AliDielectronVarManager::kCosPointingAngle, TMath::Cos(0.02), 1.0, kFALSE);
157 gammaV0Cuts->AddCut(AliDielectronVarManager::kChi2NDF, 0.0, 10.0, kFALSE);//to be checked, if properly filled
158 gammaV0Cuts->AddCut(AliDielectronVarManager::kLegDist, 0.0, 0.25, kFALSE);
159 gammaV0Cuts->AddCut(AliDielectronVarManager::kR, 3.0, 90.0, kFALSE);
160 gammaV0Cuts->AddCut(AliDielectronVarManager::kPsiPair, 0.0, 0.05, kFALSE);
161 gammaV0Cuts->AddCut(AliDielectronVarManager::kM, 0.0, 0.05, kFALSE);
162 // gammaV0Cuts->AddCut(AliDielectronVarManager::kOpeningAngle, 0.0, 0.1, kFALSE);
163 gammaV0Cuts->AddCut(AliDielectronVarManager::kArmPt, 0.0, 0.05, kFALSE);
164 // gammaV0Cuts->AddCut(AliDielectronVarManager::kArmAlpha, -0.35, 0.35, kFALSE); // not sure if it works as expected
165 gammaV0Cuts->SetExcludeTracks(kTRUE);
166 gammaV0Cuts->Print();
167
168 // const Double_t |cutAlphaG| < 0.35; && const Double_t cutQTG < 0.05;
169 // const Double_t |cutAlphaG2|[2] = {0.6, 0.8}; && const Double_t cutQTG2 < 0.04;
170
171 diele->GetTrackFilter().AddCuts(gammaV0Cuts);
172 }
173
174}
175
176//______________________________________________________________________________________
177
178void InitHistogramsDieleData(AliDielectron *diele, Int_t cutDefinition)
179{
180 //
181 // Initialise the histograms
182 //
183
184 //Setup histogram Manager
185 AliDielectronHistos *histos=new AliDielectronHistos(diele->GetName(),diele->GetTitle());
186
187 //Initialise histogram classes
188 histos->SetReservedWords("Track;Pair");
189
190 //Track classes
191 //to fill also track info from 2nd event loop until 2
192 for (Int_t i=0; i<2; ++i){
193 histos->AddClass(Form("Track_%s",AliDielectron::TrackClassName(i)));
194 }
195
196 //Pair classes
197 // to fill also mixed event histograms loop until 10
198 for (Int_t i=0; i<3/*for mixing until 10*/; ++i){
199 histos->AddClass(Form("Pair_%s",AliDielectron::PairClassName(i)));
200 }
201
202 //legs from pair
203 for (Int_t i=0; i<3; ++i){
204 histos->AddClass(Form("Track_Legs_%s",AliDielectron::PairClassName(i)));
205 }
206
207
208
209 //add histograms to event class
210 if (cutDefinition==0) {
211 histos->AddClass("Event");
212 histos->AddClass("Event_noCuts");
213 histos->UserHistogram("Event","VtxZ","Vertex Z;Z[cm]",300,-15.,15.,AliDielectronVarManager::kZvPrim);
214 // nAcc
215 histos->UserHistogram("Event","NAccRaw","Accepted raw SPD tracklets, |y|<1; nTrackl; #Entries",301,-0.5,300.5, AliDielectronVarManager::kNaccTrckltsEsd10);
216 // histos->UserHistogram("Event","NAccCorr","Accepted corr SPD tracklets, |y|<1; nTrackl; #Entries",101,-0.5,100.5, AliDielectronVarManager::kNaccTrckltsEsd10Corr);
217 // nAcc vs Zvtx
218 // histos->UserHistogram("Event","NAccRaw_vs_Zvtx","Accepted raw SPD tracklets vs Z vtx, |y|<1; Zvtx[cm]; nTrackl ",300,-15.,15.,101,-0.5,100.5, AliDielectronVarManager::kZvPrim,AliDielectronVarManager::kNaccTrckltsEsd10);
219 // histos->UserHistogram("Event","NAccCorr_vs_Zvtx","Accepted corr SPD tracklets vs Z vtx, |y|<1; Zvtx[cm]; nTrackl ",300,-15.,15.,101,-0.5,100.5, AliDielectronVarManager::kZvPrim,AliDielectronVarManager::kNaccTrckltsEsd10Corr);
220
221
222 // no event cuts
223 histos->UserHistogram("Event_noCuts","VtxZ","Vertex Z;Z[cm]",300,-15.,15.,AliDielectronVarManager::kZvPrim);
224 // nAcc
225 histos->UserHistogram("Event_noCuts","NAccRaw","Accepted raw SPD tracklets, |y|<1; nTrackl; #Entries",301,-0.5,300.5, AliDielectronVarManager::kNaccTrckltsEsd10);
226 histos->UserHistogram("Event_noCuts","NAccRaw |y|<1 eta phi","Eta; Phi; nTrackl; #Entries",100,-1,1,144,0,6.285,301,-0.5,300.5,AliDielectronVarManager::kEta, AliDielectronVarManager::kPhi, AliDielectronVarManager::kNaccTrckltsEsd10);
227 // histos->UserHistogram("Event_noCuts","NAccCorr","Accepted corr SPD tracklets, |y|<1; nTrackl; #Entries",101,-0.5,100.5, AliDielectronVarManager::kNaccTrckltsEsd10Corr);
228 // nAcc vs Zvtx
229 // histos->UserHistogram("Event_noCuts","NAccRaw_vs_Zvtx","Accepted raw SPD tracklets vs Z vtx, |y|<1; Zvtx[cm]; nTrackl ",300,-15.,15.,101,-0.5,100.5, AliDielectronVarManager::kZvPrim,AliDielectronVarManager::kNaccTrckltsEsd10);
230 // histos->UserHistogram("Event_noCuts","NAccCorr_vs_Zvtx","Accepted corr SPD tracklets vs Z vtx, |y|<1; Zvtx[cm]; nTrackl ",300,-15.,15.,101,-0.5,100.5, AliDielectronVarManager::kZvPrim,AliDielectronVarManager::kNaccTrckltsEsd10Corr);
231
232 }
233
234 //add histograms to Track classes
235 histos->UserHistogram("Track","Pt","Pt;Pt [GeV];#tracks",400,0,20.,AliDielectronVarManager::kPt);
236 histos->UserHistogram("Track","TOFPIDBit","TOFPIDBit;TOFPIDBit;#tracks",2,-0.5,1.5,AliDielectronVarManager::kTOFPIDBit);
237 histos->UserHistogram("Track","TPCnCls","Number of Clusters TPC;TPC number clusteres;#tracks",160,-0.5,159.5,AliDielectronVarManager::kNclsTPC);
238 histos->UserHistogram("Track","TPCsignalN","Number of PID Clusters TPC;TPC PID number clusteres;#tracks",160,-0.5,159.5,AliDielectronVarManager::kTPCsignalN);
239 histos->UserHistogram("Track","nClsoverfindablecluster","Number of found Clusters TPC over findably ;TPC number cluster over findable;#tracks",160,0.0,1.1,AliDielectronVarManager::kNFclsTPCrFrac);
240 histos->UserHistogram("Track","dXY","dXY;dXY [cm];#tracks",500,-1.,1.,AliDielectronVarManager::kImpactParXY);
241 histos->UserHistogram("Track","dZ","dZ;dZ [cm];#tracks",600,-3.,3.,AliDielectronVarManager::kImpactParZ);
242 histos->UserHistogram("Track","Eta_Phi","Eta Phi Map; Eta; Phi;#tracks",
243 100,-1,1,144,0,6.285,AliDielectronVarManager::kEta,AliDielectronVarManager::kPhi);
244
245 histos->UserHistogram("Track","dEdx_P","dEdx;P [GeV];TPC signal (arb units);#tracks",
246 200,0.2,20.,100,0.,200.,AliDielectronVarManager::kPIn,AliDielectronVarManager::kTPCsignal,kTRUE);
247 histos->UserHistogram("Track","TPCnSigmaEle_P","TPC number of sigmas Electrons;P [GeV];TPC number of sigmas Electrons;#tracks",
248 200,0.2,20.,100,-10.,10.,AliDielectronVarManager::kPIn,AliDielectronVarManager::kTPCnSigmaEle,kTRUE);
249 histos->UserHistogram("Track","TPCnSigmaPi_P","TPC number of sigmas Kaons;PIN [GeV];TPC number of sigmas Pions;#tracks",
250 200,0.2,20.,100,-10.,10.,AliDielectronVarManager::kPIn,AliDielectronVarManager::kTPCnSigmaPio,kTRUE);
251 histos->UserHistogram("Track","TPCnSigmaPro_P","TPC number of sigmas Protons;PIN [GeV];TPC number of sigmas Protons;#tracks",
252 200,0.2,20.,100,-10.,10.,AliDielectronVarManager::kPIn,AliDielectronVarManager::kTPCnSigmaPro,kTRUE);
253 // histos->UserHistogram("Track","TPCnSigmaKao_P","TPC number of sigmas Kaons;PIN [GeV];TPC number of sigmas Kaons;#tracks",
254 // 200,0.2,20.,100,-10.,10.,AliDielectronVarManager::kPIn,AliDielectronVarManager::kTPCnSigmaKao,kTRUE);
255 histos->UserHistogram("Track","TOFbeta_P","TOF beta;P [GeV];TOF beta;#tracks",
256 200,0.2,20.,100,0.,1.,AliDielectronVarManager::kP,AliDielectronVarManager::kTOFnSigmaEle,kTRUE);
257 // histos->UserHistogram("Track","TOFnSigmaEle_P","TOF number of sigmas Electrons;P [GeV];TOF number of sigmas Electrons;#tracks",
258 // 200,0.2,20.,200,-20.,20.,AliDielectronVarManager::kP,AliDielectronVarManager::kTOFnSigmaEle,kTRUE);
259 histos->UserHistogram("Track","TRDnCls","Number of Clusters TRD;TRD number clusters;#tracks",160,-0.5,159.5,AliDielectronVarManager::kNclsTRD);
260 histos->UserHistogram("Track","TRDntracklets","Number of tracklets TRD;TRD number tracklets;#tracks",7,-0.5,6.5,AliDielectronVarManager::kTRDntracklets);
261 histos->UserHistogram("Track","TRDprobEle_P","TRD electron prob.;P [GeV];TRD electron prob.;#tracks",
262 200,0.2,20.,100,.0,1.,AliDielectronVarManager::kP,AliDielectronVarManager::kTRDprobEle,kTRUE);
263 histos->UserHistogram("Track","TRDprobEle2D_P","TRD electron prob. 2D;P [GeV];TRD electron prob. 2D;#tracks",
264 200,0.2,20.,100,.0,1.,AliDielectronVarManager::kP,AliDielectronVarManager::kTRDprob2DEle,kTRUE);
265
266
267
268
269
270 //add histograms to Pair classes
271 histos->UserHistogram("Pair","InvMass","Inv.Mass;Inv. Mass [GeV];#pairs",
272 125,0.,125*.04,AliDielectronVarManager::kM);
273 histos->UserHistogram("Pair","Rapidity","Rapidity;Rapidity;#pairs",
274 100,-1.,1.,AliDielectronVarManager::kY);
275 histos->UserHistogram("Pair","Pt","Pt;Pt;#pairs",
276 200,0.,20.,AliDielectronVarManager::kPt);
277 histos->UserHistogram("Pair","OpeningAngle","Opening angle;angle",
278 100,0.,3.15,AliDielectronVarManager::kOpeningAngle);
279 histos->UserHistogram("Pair","OpeningAngletransverse","Opening angle transverse;angle",
280 100,0.,3.15,AliDielectronVarManager::kDeltaPhi);
281 histos->UserHistogram("Pair","Chi2NDF","chisquareNDF;chisquare/ndf;#pairs",
282 100,0.,30.,AliDielectronVarManager::kChi2NDF);
283 histos->UserHistogram("Pair","distanceXY","distancelegsXY;distanceXY[cm];#pairs",
284 100,0.,.0001,AliDielectronVarManager::kLegDistXY);
285 histos->UserHistogram("Pair","distance","distancelegs;distance[cm];#pairs",
286 100,0.,.0001,AliDielectronVarManager::kLegDist);
287 histos->UserHistogram("Pair","pseudoproperdecaylength","pseudoproperdecaylength;pseudoproperdecaylength[cm];#pairs",
288 100,0.,.5,AliDielectronVarManager::kPseudoProperTime);
289 histos->UserHistogram("Pair","Armenteros-Podolanski","Armenteros-Podolanski;ArmAlpha;ArmPt[GeV];#tracks",
290 100,-10.0,10.,100,0.,2.,AliDielectronVarManager::kArmAlpha,AliDielectronVarManager::kArmPt,kTRUE);
291
292
293
294 // 3D histos: invMass - Multiplicity - ptJpsi
295 histos->UserHistogram("Pair","InvMass_NaccRaw_PtJpsi","Inv.Mass - NaccRaw - PtJpsi;Inv. Mass [GeV];NaccRaw; pTJpsi[GeV/c]", 125,0.,125*.04,101,-0.5,100.5, 100, 0.,10., AliDielectronVarManager::kM,AliDielectronVarManager::kNaccTrckltsEsd10, AliDielectronVarManager::kPt);
296 // histos->UserHistogram("Pair","InvMass_NaccCorr_PtJpsi","Inv.Mass - NaccCorr - PtJpsi;Inv. Mass [GeV];NaccCor; pTJpsi[GeV/c]", 125,0.,125*.04,101,-0.5,100.5,100,0.,10., AliDielectronVarManager::kM,AliDielectronVarManager::kNaccTrckltsEsd10Corr, AliDielectronVarManager::kPt);
297
298 diele->SetHistogramManager(histos);
299}
300
301
302void InitCFDieleData(AliDielectron *diele, Int_t cutDefinition)
303{//number of dimensions also not the problem
304 //
305 // Setupd the CF Manager if needed
306 //
307
308 AliDielectronCF *cf=new AliDielectronCF(diele->GetName(),diele->GetTitle());
309
310 //pair variables
311 cf->AddVariable(AliDielectronVarManager::kPt,"0.0,1.0,1.3,2.0, 3.0,5.0, 7.0,10.0,100.0");
312 cf->AddVariable(AliDielectronVarManager::kY,"-1,-0.9,-0.8,-0.5,-0.3,0,0.3,0.5,0.8,0.9,1.0");
313 cf->AddVariable(AliDielectronVarManager::kM,125,0.,125*.04); //40Mev Steps
314 cf->AddVariable(AliDielectronVarManager::kPseudoProperTime,150,-0.3,0.3);
315 cf->AddVariable(AliDielectronVarManager::kPseudoProperTimeErr,600,0.,0.3);
316 cf->AddVariable(AliDielectronVarManager::kPairType,11,0,11);
317 // cf->AddVariable(AliDielectronVarManager::kArmAlpha,"-10.,-5.,-2.,-.1.,0.,1.,2.,5.,10.");
318 //cf->AddVariable(AliDielectronVarManager::kArmPt,"0.0,0.4,0.6,0.8,1.0,1.5,2.0,3.0");
319
320
321 //leg variables
322 //TODO: add variables used for cuts, debug-tree??
323 cf->AddVariable(AliDielectronVarManager::kPt,"0.0, 0.8, 1.0, 1.1, 1.2, 1.3, 100.0",kTRUE);
324 cf->AddVariable(AliDielectronVarManager::kNclsTPC,"0, 60, 65, 70, 75, 80, 85, 90, 100, 120, 160",kTRUE);
325 cf->AddVariable(AliDielectronVarManager::kTPCsignalN,"0, 50, 60, 70, 80, 90, 100, 160",kTRUE);
326 cf->AddVariable(AliDielectronVarManager::kTPCchi2Cl,"0, 0.5, 1., 1.5, 2., 2.5, 3., 3.5, 4.",kTRUE);
327 cf->AddVariable(AliDielectronVarManager::kEta, "-1.0, -0.9, -0.8,0.0,0.8, 0.9, 1.0", kTRUE);
328 cf->AddVariable(AliDielectronVarManager::kTPCnSigmaEle,"-2.5,-2,-1.5,-1,-0.5,4.",kTRUE);
329 cf->AddVariable(AliDielectronVarManager::kTPCnSigmaPio,"3.,3.5,4.,100",kTRUE);
330 cf->AddVariable(AliDielectronVarManager::kTPCnSigmaPro,"3.,3.5,4.,100",kTRUE);
331 cf->AddVariable(AliDielectronVarManager::kITSLayerFirstCls,6,-0.5,5.5,kTRUE);
332 cf->AddVariable(AliDielectronVarManager::kP,"0.0, 0.8, 0.9, 0.95, 1.0, 1.05, 1.1, 1.2,1.5,2.0,3.0,4.0,5.0, 100.0",kTRUE);
333 //event variables
334 //cf->AddVariable(AliDielectronVarManager::kNaccTrcklts,"0.0, 9.0, 17.0, 25.0, 36.0, 55.0, 500.0");
335 // cf->AddVariable(AliDielectronVarManager::kNaccTrckltsEsd10,101,-0.5,100.5);
336 // cf->AddVariable(AliDielectronVarManager::kNaccTrckltsEsd10Corr,101,-0.5,100.5);
337 //cf->AddVariable(AliDielectronVarManager::kZvPrim,"-20.,-15.,-10.,-5.,0.,5.,10.,15.,20.");
338 Bool_t hasMC=(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler()!=0x0);
339 if (hasMC){
340 cf->AddVariable(AliDielectronVarManager::kPdgCode,10000,-5000.5,4999.5,kTRUE);
341 cf->AddVariable(AliDielectronVarManager::kPdgCodeMother,10000,-5000.5,4999.5,kTRUE);
342 cf->AddVariable(AliDielectronVarManager::kPdgCodeGrandMother,10000,-5000.5,4999.5,kTRUE);
343
344 //only steps for efficiencies//seems not to help with problem!!!!, not only MC truth step, but also MCreconstructed with Signals, no MC reconstructed in total
345 //cf->SetStepsForMCtruthOnly();
346 }
347
348 //only in this case write MC truth info//NOTE perhaps only wirking for first one????
349 // if (cutDefinition==2){
350 if (hasMC) cf->SetStepForMCtruth();
351 // }
352 cf->SetStepsForSignal();
353
354 diele->SetCFManagerPair(cf);
355
356}
357
358void SetSignals(AliDielectron *diele){
359 AliDielectronSignalMC* promptJpsi = new AliDielectronSignalMC("promptJpsi","Prompt J/psi"); // prompt J/psi (not from beauty decays)
360 promptJpsi->SetLegPDGs(11,-11);
361 promptJpsi->SetMotherPDGs(443,443);//???
362 promptJpsi->SetGrandMotherPDGs(503,503,kTRUE,kTRUE); // not from beauty hadrons
363 promptJpsi->SetMothersRelation(AliDielectronSignalMC::kSame);//options kUndefined, kDifferent
364 promptJpsi->SetFillPureMCStep(kTRUE);
365 promptJpsi->SetLegSources(AliDielectronSignalMC::kFinalState, AliDielectronSignalMC::kFinalState);//kDirect produces empty Histograms...//
366 promptJpsi->SetCheckBothChargesLegs(kTRUE,kTRUE);
367 promptJpsi->SetCheckBothChargesMothers(kTRUE,kTRUE);
368 promptJpsi->SetCheckBothChargesGrandMothers(kTRUE,kTRUE);
369 diele->AddSignalMC(promptJpsi);
370
371
372 AliDielectronSignalMC* beautyJpsi = new AliDielectronSignalMC("producedbeautyJpsi","produced beauty hadron -> J/psi"); // J/psi->e+e- from beauty hadron decays
373 beautyJpsi->SetLegPDGs(11,-11);
374 beautyJpsi->SetMotherPDGs(443,443);
375 beautyJpsi->SetMothersRelation(AliDielectronSignalMC::kSame);//test
376 beautyJpsi->SetMotherSources(AliDielectronSignalMC::kPrimary,AliDielectronSignalMC::kPrimary);//added lines!!!
377 beautyJpsi->SetGrandMotherSources(AliDielectronSignalMC::kPrimary,AliDielectronSignalMC::kPrimary);//added lines!!!!
378 beautyJpsi->SetGrandMotherPDGs(503,503);
379 beautyJpsi->SetFillPureMCStep(kTRUE); //does not help to take it out!!!
380 beautyJpsi->SetLegSources(AliDielectronSignalMC::kFinalState, AliDielectronSignalMC::kFinalState);//possibilities: kFinalState, kDirect, kPrimary, kDontCare, kSecondary, kNoCocktail
381 // beautyJpsi->SetJpsiRadiative(AliDielectronSignalMC::kIsNotRadiative);//tested line, works, but is not solve problem!!!
382 beautyJpsi->SetCheckBothChargesLegs(kTRUE,kTRUE);
383 beautyJpsi->SetCheckBothChargesMothers(kTRUE,kTRUE);
384 beautyJpsi->SetCheckBothChargesGrandMothers(kTRUE,kTRUE);
385 diele->AddSignalMC(beautyJpsi);
386
387
388 AliDielectronSignalMC* background1 = new AliDielectronSignalMC("backgroundnolegjpsi","backgroundnolegjpsi"); //not J/psi->e+e-
389 //background1->SetLegPDGs(11,-11);
390 background1->SetMotherPDGs(443,443, kTRUE,kTRUE );
391 background1->SetMothersRelation(AliDielectronSignalMC::kUndefined);
392 background1->SetGrandMotherPDGs(503,503);
393 background1->SetFillPureMCStep(kFALSE);
394 background1->SetLegSources(AliDielectronSignalMC::kFinalState, AliDielectronSignalMC::kFinalState);
395 background1->SetCheckBothChargesLegs(kTRUE,kTRUE);
396 diele->AddSignalMC(background1);
397
398AliDielectronSignalMC* background2 = new AliDielectronSignalMC("backgroundonelegjpsi","backgroundonelegjpsi"); //not J/psi->e+e-
399// background2->SetLegPDGs(11,-11);
400 background2->SetMotherPDGs(443,443, kTRUE,kFALSE );
401 background2->SetMothersRelation(AliDielectronSignalMC::kUndefined);
402 background2->SetFillPureMCStep(kFALSE);
403 background2->SetLegSources(AliDielectronSignalMC::kFinalState, AliDielectronSignalMC::kFinalState);
404 background2->SetCheckBothChargesLegs(kTRUE,kTRUE);
405 diele->AddSignalMC(background2);
406 /*
407 AliDielectronSignalMC* beautyMesonJpsi = new AliDielectronSignalMC("beautyMesonJpsi","beauty meson -> J/psi"); // J/psi->e+e- from beauty hadron decays
408 beautyMesonJpsi->SetLegPDGs(11,-11);
409 beautyMesonJpsi->SetMotherPDGs(443,443);
410 beautyMesonJpsi->SetMothersRelation(AliDielectronSignalMC::kSame);
411 beautyMesonJpsi->SetGrandMotherPDGs(500,500);
412 beautyMesonJpsi->SetFillPureMCStep(kTRUE);
413 beautyMesonJpsi->SetLegSources(AliDielectronSignalMC::kFinalState, AliDielectronSignalMC::kFinalState);
414 beautyMesonJpsi->SetCheckBothChargesLegs(kTRUE,kTRUE);
415 beautyMesonJpsi->SetCheckBothChargesMothers(kTRUE,kTRUE);
416 beautyMesonJpsi->SetCheckBothChargesGrandMothers(kTRUE,kTRUE);
417 diele->AddSignalMC(beautyMesonJpsi);
418
419
420 // physical backgrounds (electrons from other sources)
421 AliDielectronSignalMC* diEleContinuum = new AliDielectronSignalMC("diEleContinuum","di-electron continuum"); // all di-electrons originating in the collision
422 diEleContinuum->SetLegPDGs(11,-11);
423 diEleContinuum->SetLegSources(AliDielectronSignalMC::kFinalState, AliDielectronSignalMC::kFinalState);
424 diEleContinuum->SetCheckBothChargesLegs(kTRUE,kTRUE);
425 diele->AddSignalMC(diEleContinuum);
426
427 AliDielectronSignalMC* diEleCharm = new AliDielectronSignalMC("diEleCharm","di-electrons from charm"); // dielectrons originating from charm hadrons (not neccessary from same mother)
428 diEleCharm->SetLegPDGs(11,-11);
429 diEleCharm->SetMotherPDGs(403,403);
430 diEleCharm->SetLegSources(AliDielectronSignalMC::kFinalState, AliDielectronSignalMC::kFinalState);
431 diEleCharm->SetCheckBothChargesLegs(kTRUE,kTRUE);
432 diEleCharm->SetCheckBothChargesMothers(kTRUE,kTRUE);
433 diele->AddSignalMC(diEleCharm);
434
435 AliDielectronSignalMC* diEleOpenCharm = new AliDielectronSignalMC("diEleOpenCharm","di-electrons from open charm"); // dielectrons originating from open charm hadrons
436 diEleOpenCharm->SetLegPDGs(11,-11);
437 diEleOpenCharm->SetMotherPDGs(402,402);
438 diEleOpenCharm->SetLegSources(AliDielectronSignalMC::kFinalState, AliDielectronSignalMC::kFinalState);
439 diEleOpenCharm->SetCheckBothChargesLegs(kTRUE,kTRUE);
440 diEleOpenCharm->SetCheckBothChargesMothers(kTRUE,kTRUE);
441 diele->AddSignalMC(diEleOpenCharm);
442
443 AliDielectronSignalMC* diEleOpenCharmJpsi = new AliDielectronSignalMC("diEleOpenCharmJpsi","1 leg from open charm + 1 jpsi leg"); // 1 leg from open charm + 1 leg from jpsi
444 diEleOpenCharmJpsi->SetLegPDGs(11,-11);
445 diEleOpenCharmJpsi->SetMotherPDGs(402,443);
446 diEleOpenCharmJpsi->SetLegSources(AliDielectronSignalMC::kFinalState, AliDielectronSignalMC::kFinalState);
447 diEleOpenCharmJpsi->SetCheckBothChargesLegs(kTRUE,kTRUE);
448 diEleOpenCharmJpsi->SetCheckBothChargesMothers(kTRUE,kTRUE);
449 diele->AddSignalMC(diEleOpenCharmJpsi);
450
451 AliDielectronSignalMC* muonMuonPairs = new AliDielectronSignalMC("muonMuonPairs","muon+muon pairs"); // dimuon pairs, all sources (both from collision and from material interaction)
452 muonMuonPairs->SetLegPDGs(13,13);
453 muonMuonPairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
454 diele->AddSignalMC(muonMuonPairs);
455
456 // background from secondary electrons
457 AliDielectronSignalMC* secondaryElectrons = new AliDielectronSignalMC("secondaryElectrons","Secondary electrons"); // all di-electrons from secondary electrons (interaction with detector)
458 secondaryElectrons->SetLegPDGs(11,-11);
459 secondaryElectrons->SetLegSources(AliDielectronSignalMC::kSecondary, AliDielectronSignalMC::kSecondary);
460 secondaryElectrons->SetCheckBothChargesLegs(kTRUE,kTRUE);
461 diele->AddSignalMC(secondaryElectrons);
462
463 AliDielectronSignalMC* primarySecElePairs = new AliDielectronSignalMC("primarySecElePairs","Primary+Secondary electron pairs"); // primary-secondary pairs
464 primarySecElePairs->SetLegPDGs(11,-11);
465 primarySecElePairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
466 primarySecElePairs->SetLegSources(AliDielectronSignalMC::kFinalState, AliDielectronSignalMC::kSecondary);
467 diele->AddSignalMC(primarySecElePairs);
468
469 AliDielectronSignalMC* conversionElePairs = new AliDielectronSignalMC("conversionElePairs","conversion electron pairs"); // pairs made from conversion (may be also from 2 different conversions)
470 conversionElePairs->SetLegPDGs(11,-11);
471 conversionElePairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
472 conversionElePairs->SetLegSources(AliDielectronSignalMC::kSecondary, AliDielectronSignalMC::kSecondary);
473 conversionElePairs->SetMotherPDGs(22,22);
474 diele->AddSignalMC(conversionElePairs);
475
476 // misidentification
477 AliDielectronSignalMC* allEleMisIdPairs = new AliDielectronSignalMC("allEleMisIdPairs","all electron+misid. pairs"); // one true electron + a mis-id electron (all sources included)
478 allEleMisIdPairs->SetLegPDGs(11,11,kFALSE,kTRUE);
479 allEleMisIdPairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
480 diele->AddSignalMC(allEleMisIdPairs);
481
482 AliDielectronSignalMC* allMisIdMisIdPairs = new AliDielectronSignalMC("allMisIdMisIdPairs","all misid.+misid. pairs"); // mis-id + mis-id
483 allMisIdMisIdPairs->SetLegPDGs(11,11,kTRUE,kTRUE);
484 allMisIdMisIdPairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
485 diele->AddSignalMC(allMisIdMisIdPairs);
486
487 AliDielectronSignalMC* elePionPairs = new AliDielectronSignalMC("elePionPairs","electron+pion pairs"); // true electron + mis-id pion
488 elePionPairs->SetLegPDGs(11,211);
489 elePionPairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
490 diele->AddSignalMC(elePionPairs);
491
492 AliDielectronSignalMC* eleKaonPairs = new AliDielectronSignalMC("eleKaonPairs","electron+kaon pairs"); // true electron + mis-id kaon
493 eleKaonPairs->SetLegPDGs(11,321);
494 eleKaonPairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
495 diele->AddSignalMC(eleKaonPairs);
496
497 AliDielectronSignalMC* eleProtonPairs = new AliDielectronSignalMC("eleProtonPairs","Electron+proton pairs"); // true electron + mis-id proton
498 eleProtonPairs->SetLegPDGs(11,2212);
499 eleProtonPairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
500 diele->AddSignalMC(eleProtonPairs);
501
502 AliDielectronSignalMC* piPiPairs = new AliDielectronSignalMC("piPiPairs","pion+pion pairs"); // mis-id pion + mis-id pion
503 piPiPairs->SetLegPDGs(211,211);
504 piPiPairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
505 diele->AddSignalMC(piPiPairs);
506
507 AliDielectronSignalMC* piKaonPairs = new AliDielectronSignalMC("piKaonPairs","pion+kaon pairs"); // mis-id pion + mis-id kaon
508 piKaonPairs->SetLegPDGs(211,321);
509 piKaonPairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
510 diele->AddSignalMC(piKaonPairs);
511
512 AliDielectronSignalMC* piProtonPairs = new AliDielectronSignalMC("piProtonPairs","pion+proton pairs"); // mis-id pion + mis-id proton
513 piProtonPairs->SetLegPDGs(211,2212);
514 piProtonPairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
515 diele->AddSignalMC(piProtonPairs);
516
517 AliDielectronSignalMC* kaonKaonPairs = new AliDielectronSignalMC("kaonKaonPairs","kaon+kaon pairs"); // mis-id kaon + mis-id kaon
518 kaonKaonPairs->SetLegPDGs(321,321);
519 kaonKaonPairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
520 diele->AddSignalMC(kaonKaonPairs);
521
522 AliDielectronSignalMC* kaonProtonPairs = new AliDielectronSignalMC("kaonProtonPairs","kaon+proton pairs"); // mis-id kaon + mis-id proton
523 kaonProtonPairs->SetLegPDGs(321,2212);
524 kaonProtonPairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
525 diele->AddSignalMC(kaonProtonPairs);
526
527 AliDielectronSignalMC* protonProtonPairs = new AliDielectronSignalMC("protonProtonPairs","proton+proton pairs"); // mis-id proton + mis-id proton
528 protonProtonPairs->SetLegPDGs(2212,2212);
529 protonProtonPairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
530 diele->AddSignalMC(protonProtonPairs);
531
532 AliDielectronSignalMC* muonAllPairs = new AliDielectronSignalMC("muonAllPairs","muon+everything pairs"); // mis-id muon + something else (electron, pion, kaon, proton)
533 muonAllPairs->SetLegPDGs(13,13,kFALSE,kTRUE);
534 muonAllPairs->SetCheckBothChargesLegs(kTRUE,kTRUE);
535 diele->AddSignalMC(muonAllPairs);
536 */
537
538}