]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FEMTOSCOPY/macros/Train/TwoTrackQA/ConfigFemtoAnalysisCentral.C
8ad15bae25b9f5ebccc45a692b8b942c895f0cf2
[u/mrichter/AliRoot.git] / PWG2 / FEMTOSCOPY / macros / Train / TwoTrackQA / ConfigFemtoAnalysisCentral.C
1 /*********************************************************************
2  *                                                                   *
3  * ConfigFemtoAnalysis.C - configuration macro for the femtoscopic   *
4  * analysis, meant as a QA process for two-particle effects          *
5  *                                                                   *
6  * Author: Adam Kisiel (Adam.Kisiel@cern.ch)                         *
7  *                                                                   *
8  *********************************************************************/
9
10 #if !defined(__CINT__) || defined(__MAKECINT_)
11 #include "AliFemtoManager.h"
12 #include "AliFemtoEventReaderESDChain.h"
13 #include "AliFemtoEventReaderESDChainKine.h"
14 #include "AliFemtoSimpleAnalysis.h"
15 #include "AliFemtoBasicEventCut.h"
16 #include "AliFemtoESDTrackCut.h"
17 #include "AliFemtoCorrFctn.h"
18 #include "AliFemtoCutMonitorParticleYPt.h"
19 #include "AliFemtoCutMonitorParticleVertPos.h"
20 #include "AliFemtoCutMonitorParticleMomRes.h"
21 #include "AliFemtoCutMonitorEventMult.h"
22 #include "AliFemtoCutMonitorEventVertex.h"
23 #include "AliFemtoShareQualityTPCEntranceSepPairCut.h"
24 #include "AliFemtoQinvCorrFctn.h"
25 #include "AliFemtoShareQualityCorrFctn.h"
26 #include "AliFemtoTPCInnerCorrFctn.h"
27 #include "AliFemtoVertexMultAnalysis.h"
28 #include "AliFemtoCorrFctn3DSpherical.h"
29 #include "AliFemtoChi2CorrFctn.h"
30 #include "AliFemtoCorrFctnTPCNcls.h"
31 #include "AliFemtoModelBPLCMSCorrFctn.h"
32 #include "AliFemtoModelCorrFctn3DLCMSSym.h"
33 #include "AliFemtoModelCorrFctn3DSpherical.h"
34 #include "AliFemtoModelGausLCMSFreezeOutGenerator.h"
35 #include "AliFemtoModelGausRinvFreezeOutGenerator.h"
36 #include "AliFemtoModelManager.h"
37 #include "AliFemtoModelWeightGeneratorBasic.h"
38 #include "AliFemtoModelWeightGeneratorLednicky.h"
39 #include "AliFemtoCorrFctnDirectYlm.h"
40 #include "AliFemtoModelCorrFctnDirectYlm.h"
41 #include "AliFemtoModelCorrFctnSource.h"
42 #include "AliFemtoCutMonitorParticlePtPDG.h"
43 #include "AliFemtoKTPairCut.h"
44 #endif
45
46 //________________________________________________________________________
47 AliFemtoManager* ConfigFemtoAnalysis(int doESD=1, int iCollision=1, int iCent=0) {
48
49   double PionMass = 0.13956995;
50   double KaonMass = 0.493677;
51   
52 //   AliFemtoEventReaderESDChainKine* Reader=new AliFemtoEventReaderESDChainKine();
53 //   Reader->SetConstrained(true);
54 //   Reader->SetUseTPCOnly(false);
55
56   int multmin=0, multmax=0;
57   double zmin, zmax;
58
59   int ppranges[9] = { 1, 12, 17, 23, 30, 37, 45, 58, 150 };
60   int aaranges[11] = { 0, 50, 100, 200, 300, 400, 500, 600, 700, 800, 900 };
61
62   if (iCollision == 0) { // pp setup
63     if ((iCent >=0) && (iCent<8)) {
64       multmin = ppranges[iCent];
65       multmax = ppranges[iCent+1];
66     }
67     zmin = -12.0;
68     zmax = 12.0;
69   }
70   if (iCollision == 1) { // PbPb setup
71     if ((iCent >=0) && (iCent<10)) {
72       multmin = aaranges[iCent];
73       multmax = aaranges[iCent+1];
74     }
75     zmin = -8.0;
76     zmax = 8.0;
77   }
78
79   AliFemtoEventReaderESDChain* Reader=new AliFemtoEventReaderESDChain();
80   Reader->SetConstrained(true);
81   Reader->SetUseTPCOnly(false);
82
83   if (iCollision == 0) 
84     Reader->SetUseMultiplicity(AliFemtoEventReaderESDChain::kGlobalCount);
85   if (iCollision == 1) 
86     Reader->SetUseMultiplicity(AliFemtoEventReaderESDChain::kV0Centrality);
87
88   AliFemtoManager* Manager=new AliFemtoManager();
89   Manager->SetEventReader(Reader);
90         
91   // Switches for QA analyses
92   int runPiPlusStandard = 1;
93   int runPiMinusStandard = 1;
94   int runPositiveQA = 1;
95   int runNegativeQA = 1;
96   int runPositiveTPCQA = 1;
97   int runNegativeTPCQA = 0;
98
99   // *** First QA task - standard HBT analysis with all the cut on ***
100
101   // *** Begin pion-pion (positive) analysis ***
102   if (runPiPlusStandard) {
103     AliFemtoVertexMultAnalysis *anpipstd = new AliFemtoVertexMultAnalysis(1, zmin, zmax, 5, multmin, multmax);
104     anpipstd->SetNumEventsToMix(iCollision?3:10);
105     anpipstd->SetMinSizePartCollection(2);
106
107     AliFemtoBasicEventCut* mecpipstd = new AliFemtoBasicEventCut();
108     mecpipstd->SetEventMult(0,100000);
109     mecpipstd->SetVertZPos(-1000,1000);
110         
111     AliFemtoESDTrackCut* dtcpipstd = new AliFemtoESDTrackCut();
112     dtcpipstd->SetPidProbPion(0.2,1.001);
113     dtcpipstd->SetPidProbMuon(0.0,1.0);
114     dtcpipstd->SetPidProbKaon(0.0,1.0);
115     dtcpipstd->SetPidProbProton(0.0,1.0);
116     dtcpipstd->SetCharge(1.0);
117     dtcpipstd->SetPt(0.15,0.5);
118     dtcpipstd->SetMass(PionMass);
119     // Track quality cuts
120     dtcpipstd->SetStatus(AliESDtrack::kTPCrefit|AliESDtrack::kITSrefit);
121     //  dtcpipstd->SetStatus(AliESDtrack::kTPCrefit);
122     dtcpipstd->SetminTPCncls(50);
123     dtcpipstd->SetRemoveKinks(kTRUE);
124     dtcpipstd->SetLabel(kFALSE);
125     dtcpipstd->SetMaxITSChiNdof(2.5);
126     dtcpipstd->SetMaxTPCChiNdof(3.0);
127     dtcpipstd->SetMaxImpactXY(3.0);
128     dtcpipstd->SetMaxImpactZ(3.0);
129
130     AliFemtoCutMonitorParticleYPt *cutPassYPtpipstd = new AliFemtoCutMonitorParticleYPt("cutPasspipstd", 0.13957);
131     AliFemtoCutMonitorParticleYPt *cutFailYPtpipstd = new AliFemtoCutMonitorParticleYPt("cutFailpipstd", 0.13957);
132     dtcpipstd->AddCutMonitor(cutPassYPtpipstd, cutFailYPtpipstd);
133
134     AliFemtoShareQualityTPCEntranceSepPairCut *sqpcpipstd = new AliFemtoShareQualityTPCEntranceSepPairCut();
135     sqpcpipstd->SetShareQualityMax(0.0);
136     sqpcpipstd->SetShareFractionMax(0.02);
137     sqpcpipstd->SetRemoveSameLabel(kFALSE);
138     sqpcpipstd->SetTPCEntranceSepMinimum(2.0);
139
140     anpipstd->SetEventCut(mecpipstd);
141     anpipstd->SetFirstParticleCut(dtcpipstd);
142     anpipstd->SetSecondParticleCut(dtcpipstd);
143     anpipstd->SetPairCut(sqpcpipstd);
144     
145     AliFemtoShareQualityCorrFctn *csqqinvpipstd= new AliFemtoShareQualityCorrFctn("sqqinvcfpipstd",40,0.0,0.4);
146     anpipstd->AddCorrFctn(csqqinvpipstd);
147
148     AliFemtoCorrFctnDirectYlm *cylmpipstd = new AliFemtoCorrFctnDirectYlm("cylmpipstd",3,60,0.0,0.3,1);
149     anpipstd->AddCorrFctn(cylmpipstd);
150     
151     AliFemtoQinvCorrFctn *cqinvpipstd = new AliFemtoQinvCorrFctn("qinvcfpipstd", 100,0.0,1.0);
152     anpipstd->AddCorrFctn(cqinvpipstd);
153
154     AliFemtoChi2CorrFctn *cchiqinvpipstd= new AliFemtoChi2CorrFctn("chicfpipstd",40,0.0,0.4);
155     anpipstd->AddCorrFctn(cchiqinvpipstd);
156
157     AliFemtoCorrFctnTPCNcls *cqtpcnclspipstd = new AliFemtoCorrFctnTPCNcls("cqtpcnclspipstd",40,0.0,0.4);
158     anpipstd->AddCorrFctn(cqtpcnclspipstd);
159
160     AliFemtoCorrFctn3DLCMSSym *cq3dlcmspipstd = new AliFemtoCorrFctn3DLCMSSym("cq3dlcmspipstd",60,0.3);
161     anpipstd->AddCorrFctn(cq3dlcmspipstd);
162
163     Manager->AddAnalysis(anpipstd);     
164   }
165   // *** End pion-pion (positive) analysis
166
167   // *** Begin pion-pion (negative) analysis ***
168   if (runPiMinusStandard) {
169     AliFemtoVertexMultAnalysis *anpimstd = new AliFemtoVertexMultAnalysis(1, zmin, zmax, 5, multmin, multmax);
170     anpimstd->SetNumEventsToMix(iCollision?3:10);
171     anpimstd->SetMinSizePartCollection(2);
172
173     AliFemtoBasicEventCut* mecpimstd = new AliFemtoBasicEventCut();
174     mecpimstd->SetEventMult(0,100000);
175     mecpimstd->SetVertZPos(-1000,1000);
176         
177     AliFemtoESDTrackCut* dtcpimstd = new AliFemtoESDTrackCut();
178     dtcpimstd->SetPidProbPion(0.2,1.001);
179     dtcpimstd->SetPidProbMuon(0.0,1.0);
180     dtcpimstd->SetPidProbKaon(0.0,1.0);
181     dtcpimstd->SetPidProbProton(0.0,1.0);
182     dtcpimstd->SetCharge(-1.0);
183     dtcpimstd->SetPt(0.15,0.5);
184     dtcpimstd->SetMass(PionMass);
185     // Track quality cuts
186     dtcpimstd->SetStatus(AliESDtrack::kTPCrefit|AliESDtrack::kITSrefit);
187     //  dtcpimstd->SetStatus(AliESDtrack::kTPCrefit);
188     dtcpimstd->SetminTPCncls(50);
189     dtcpimstd->SetRemoveKinks(kTRUE);
190     dtcpimstd->SetLabel(kFALSE);
191     dtcpimstd->SetMaxITSChiNdof(2.5);
192     dtcpimstd->SetMaxTPCChiNdof(3.0);
193     dtcpimstd->SetMaxImpactXY(3.0);
194     dtcpimstd->SetMaxImpactZ(3.0);
195
196     AliFemtoCutMonitorParticleYPt *cutPassYPtpimstd = new AliFemtoCutMonitorParticleYPt("cutPasspimstd", 0.13957);
197     AliFemtoCutMonitorParticleYPt *cutFailYPtpimstd = new AliFemtoCutMonitorParticleYPt("cutFailpimstd", 0.13957);
198     dtcpimstd->AddCutMonitor(cutPassYPtpimstd, cutFailYPtpimstd);
199
200     AliFemtoShareQualityTPCEntranceSepPairCut *sqpcpimstd = new AliFemtoShareQualityTPCEntranceSepPairCut();
201     sqpcpimstd->SetShareQualityMax(0.0);
202     sqpcpimstd->SetShareFractionMax(0.02);
203     sqpcpimstd->SetRemoveSameLabel(kFALSE);
204     sqpcpimstd->SetTPCEntranceSepMinimum(2.0);
205
206     anpimstd->SetEventCut(mecpimstd);
207     anpimstd->SetFirstParticleCut(dtcpimstd);
208     anpimstd->SetSecondParticleCut(dtcpimstd);
209     anpimstd->SetPairCut(sqpcpimstd);
210     
211     AliFemtoShareQualityCorrFctn *csqqinvpimstd= new AliFemtoShareQualityCorrFctn("sqqinvcfpimstd",40,0.0,0.4);
212     anpimstd->AddCorrFctn(csqqinvpimstd);
213
214     AliFemtoCorrFctnDirectYlm *cylmpimstd = new AliFemtoCorrFctnDirectYlm("cylmpimstd",3,60,0.0,0.3,1);
215     anpimstd->AddCorrFctn(cylmpimstd);
216     
217     AliFemtoQinvCorrFctn *cqinvpimstd = new AliFemtoQinvCorrFctn("qinvcfpimstd", 100,0.0,1.0);
218     anpimstd->AddCorrFctn(cqinvpimstd);
219
220     AliFemtoChi2CorrFctn *cchiqinvpimstd= new AliFemtoChi2CorrFctn("chicfpimstd",40,0.0,0.4);
221     anpimstd->AddCorrFctn(cchiqinvpimstd);
222
223     AliFemtoCorrFctnTPCNcls *cqtpcnclspimstd = new AliFemtoCorrFctnTPCNcls("cqtpcnclspimstd",40,0.0,0.4);
224     anpimstd->AddCorrFctn(cqtpcnclspimstd);
225
226     AliFemtoCorrFctn3DLCMSSym *cq3dlcmspimstd = new AliFemtoCorrFctn3DLCMSSym("cq3dlcmspimstd",60,0.3);
227     anpimstd->AddCorrFctn(cq3dlcmspimstd);
228
229     Manager->AddAnalysis(anpimstd);     
230   }
231   // *** End pion-pion (negative) analysis
232
233   // *** Second QA task - HBT analysis with all pair cuts off ***
234   // *** Begin pion-pion (positive) analysis ***
235   if (runPositiveQA) {
236     AliFemtoVertexMultAnalysis *anpipnct = new AliFemtoVertexMultAnalysis(1, zmin, zmax, 5, multmin, multmax);
237     anpipnct->SetNumEventsToMix(iCollision?3:10);
238     anpipnct->SetMinSizePartCollection(2);
239
240     AliFemtoBasicEventCut* mecpipnct = new AliFemtoBasicEventCut();
241     mecpipnct->SetEventMult(0,100000);
242     mecpipnct->SetVertZPos(-1000,1000);
243         
244     AliFemtoESDTrackCut* dtcpipnct = new AliFemtoESDTrackCut();
245     dtcpipnct->SetPidProbPion(0.0,1.001);
246     dtcpipnct->SetPidProbMuon(0.0,1.0);
247     dtcpipnct->SetPidProbKaon(0.0,1.0);
248     dtcpipnct->SetPidProbProton(0.0,1.0);
249     dtcpipnct->SetCharge(1.0);
250     dtcpipnct->SetPt(0.1,1.0);
251     dtcpipnct->SetMass(PionMass);
252     // Track quality cuts
253     dtcpipnct->SetStatus(AliESDtrack::kTPCrefit|AliESDtrack::kITSrefit);
254     //  dtcpipnct->SetStatus(AliESDtrack::kTPCrefit);
255     dtcpipnct->SetminTPCncls(50);
256     dtcpipnct->SetRemoveKinks(kTRUE);
257     dtcpipnct->SetLabel(kFALSE);
258     dtcpipnct->SetMaxITSChiNdof(6.0);
259     dtcpipnct->SetMaxTPCChiNdof(6.0);
260     dtcpipnct->SetMaxImpactXY(3.0);
261     dtcpipnct->SetMaxImpactZ(3.0);
262
263     AliFemtoCutMonitorParticleYPt *cutPassYPtpipnct = new AliFemtoCutMonitorParticleYPt("cutPasspipnct", 0.13957);
264     AliFemtoCutMonitorParticleYPt *cutFailYPtpipnct = new AliFemtoCutMonitorParticleYPt("cutFailpipnct", 0.13957);
265     dtcpipnct->AddCutMonitor(cutPassYPtpipnct, cutFailYPtpipnct);
266
267     AliFemtoShareQualityTPCEntranceSepPairCut *sqpcpipnct = new AliFemtoShareQualityTPCEntranceSepPairCut();
268     sqpcpipnct->SetShareQualityMax(1.0);
269     sqpcpipnct->SetShareFractionMax(1.0);
270     sqpcpipnct->SetRemoveSameLabel(kFALSE);
271     sqpcpipnct->SetTPCEntranceSepMinimum(0.0);
272
273     anpipnct->SetEventCut(mecpipnct);
274     anpipnct->SetFirstParticleCut(dtcpipnct);
275     anpipnct->SetSecondParticleCut(dtcpipnct);
276     anpipnct->SetPairCut(sqpcpipnct);
277     
278     AliFemtoShareQualityCorrFctn *csqqinvpipnct= new AliFemtoShareQualityCorrFctn("sqqinvcfpipnct",40,0.0,0.4);
279     anpipnct->AddCorrFctn(csqqinvpipnct);
280
281     AliFemtoCorrFctnDirectYlm *cylmpipnct = new AliFemtoCorrFctnDirectYlm("cylmpipnct",3,60,0.0,0.3,1);
282     anpipnct->AddCorrFctn(cylmpipnct);
283     
284     AliFemtoQinvCorrFctn *cqinvpipnct = new AliFemtoQinvCorrFctn("qinvcfpipnct", 100,0.0,1.0);
285     anpipnct->AddCorrFctn(cqinvpipnct);
286
287     AliFemtoChi2CorrFctn *cchiqinvpipnct= new AliFemtoChi2CorrFctn("chicfpipnct",40,0.0,0.4);
288     anpipnct->AddCorrFctn(cchiqinvpipnct);
289
290     AliFemtoCorrFctnTPCNcls *cqtpcnclspipnct = new AliFemtoCorrFctnTPCNcls("cqtpcnclspipnct",40,0.0,0.4);
291     anpipnct->AddCorrFctn(cqtpcnclspipnct);
292
293     Manager->AddAnalysis(anpipnct);     
294   }
295   // *** End pion-pion (positive) analysis
296
297   // *** Begin pion-pion (negative) analysis ***
298   if (runNegativeQA) {
299     AliFemtoVertexMultAnalysis *anpimnct = new AliFemtoVertexMultAnalysis(1, zmin, zmax, 5, multmin, multmax);
300     anpimnct->SetNumEventsToMix(iCollision?3:10);
301     anpimnct->SetMinSizePartCollection(2);
302
303     AliFemtoBasicEventCut* mecpimnct = new AliFemtoBasicEventCut();
304     mecpimnct->SetEventMult(0,100000);
305     mecpimnct->SetVertZPos(-1000,1000);
306         
307     AliFemtoESDTrackCut* dtcpimnct = new AliFemtoESDTrackCut();
308     dtcpimnct->SetPidProbPion(0.0,1.001);
309     dtcpimnct->SetPidProbMuon(0.0,1.0);
310     dtcpimnct->SetPidProbKaon(0.0,1.0);
311     dtcpimnct->SetPidProbProton(0.0,1.0);
312     dtcpimnct->SetCharge(-1.0);
313     dtcpimnct->SetPt(0.1,1.0);
314     dtcpimnct->SetMass(PionMass);
315     // Track quality cuts
316     dtcpimnct->SetStatus(AliESDtrack::kTPCrefit|AliESDtrack::kITSrefit);
317     //  dtcpimnct->SetStatus(AliESDtrack::kTPCrefit);
318     dtcpimnct->SetminTPCncls(50);
319     dtcpimnct->SetRemoveKinks(kTRUE);
320     dtcpimnct->SetLabel(kFALSE);
321     dtcpimnct->SetMaxITSChiNdof(6.0);
322     dtcpimnct->SetMaxTPCChiNdof(6.0);
323     dtcpimnct->SetMaxImpactXY(3.0);
324     dtcpimnct->SetMaxImpactZ(3.0);
325
326     AliFemtoCutMonitorParticleYPt *cutPassYPtpimnct = new AliFemtoCutMonitorParticleYPt("cutPasspimnct", 0.13957);
327     AliFemtoCutMonitorParticleYPt *cutFailYPtpimnct = new AliFemtoCutMonitorParticleYPt("cutFailpimnct", 0.13957);
328     dtcpimnct->AddCutMonitor(cutPassYPtpimnct, cutFailYPtpimnct);
329
330     AliFemtoShareQualityTPCEntranceSepPairCut *sqpcpimnct = new AliFemtoShareQualityTPCEntranceSepPairCut();
331     sqpcpimnct->SetShareQualityMax(1.0);
332     sqpcpimnct->SetShareFractionMax(1.0);
333     sqpcpimnct->SetRemoveSameLabel(kFALSE);
334     sqpcpimnct->SetTPCEntranceSepMinimum(0.0);
335
336     anpimnct->SetEventCut(mecpimnct);
337     anpimnct->SetFirstParticleCut(dtcpimnct);
338     anpimnct->SetSecondParticleCut(dtcpimnct);
339     anpimnct->SetPairCut(sqpcpimnct);
340     
341     AliFemtoShareQualityCorrFctn *csqqinvpimnct= new AliFemtoShareQualityCorrFctn("sqqinvcfpimnct",40,0.0,0.4);
342     anpimnct->AddCorrFctn(csqqinvpimnct);
343
344     AliFemtoCorrFctnDirectYlm *cylmpimnct = new AliFemtoCorrFctnDirectYlm("cylmpimnct",3,60,0.0,0.3,1);
345     anpimnct->AddCorrFctn(cylmpimnct);
346     
347     AliFemtoQinvCorrFctn *cqinvpimnct = new AliFemtoQinvCorrFctn("qinvcfpimnct", 100,0.0,1.0);
348     anpimnct->AddCorrFctn(cqinvpimnct);
349
350     AliFemtoChi2CorrFctn *cchiqinvpimnct= new AliFemtoChi2CorrFctn("chicfpimnct",40,0.0,0.4);
351     anpimnct->AddCorrFctn(cchiqinvpimnct);
352
353     AliFemtoCorrFctnTPCNcls *cqtpcnclspimnct = new AliFemtoCorrFctnTPCNcls("cqtpcnclspimnct",40,0.0,0.4);
354     anpimnct->AddCorrFctn(cqtpcnclspimnct);
355
356     Manager->AddAnalysis(anpimnct);     
357   }
358   // *** End pion-pion (negative) analysis
359
360   // *** Third QA task - HBT analysis with all pair cuts off, TPC only ***
361   // *** Begin pion-pion (positive) analysis ***
362   if (runPositiveTPCQA) {
363     AliFemtoVertexMultAnalysis *anpiptpc = new AliFemtoVertexMultAnalysis(1, zmin, zmax, 5, multmin, multmax);
364     anpiptpc->SetNumEventsToMix(iCollision?3:10);
365     anpiptpc->SetMinSizePartCollection(2);
366
367     AliFemtoBasicEventCut* mecpiptpc = new AliFemtoBasicEventCut();
368     mecpiptpc->SetEventMult(0,100000);
369     mecpiptpc->SetVertZPos(-1000,1000);
370
371     AliFemtoCutMonitorEventMult *cutPassEvMpiptpc = new AliFemtoCutMonitorEventMult("cutPasspiptpc");
372     AliFemtoCutMonitorEventMult *cutFailEvMpiptpc = new AliFemtoCutMonitorEventMult("cutFailpiptpc");
373     mecpiptpc->AddCutMonitor(cutPassEvMpiptpc, cutFailEvMpiptpc);
374
375     AliFemtoCutMonitorEventVertex *cutPassEvVpiptpc = new AliFemtoCutMonitorEventVertex("cutPasspiptpc");
376     AliFemtoCutMonitorEventVertex *cutFailEvVpiptpc = new AliFemtoCutMonitorEventVertex("cutFailpiptpc");
377     mecpiptpc->AddCutMonitor(cutPassEvVpiptpc, cutFailEvVpiptpc);
378         
379     AliFemtoESDTrackCut* dtcpiptpc = new AliFemtoESDTrackCut();
380 //     dtcpiptpc->SetPidProbPion(0.0,1.001);
381 //     dtcpiptpc->SetPidProbMuon(0.0,1.0);
382 //     dtcpiptpc->SetPidProbKaon(0.0,1.0);
383 //     dtcpiptpc->SetPidProbProton(0.0,1.0);
384     dtcpiptpc->SetCharge(0.0);
385     dtcpiptpc->SetPt(0.01,1000.0);
386     dtcpiptpc->SetMass(PionMass);
387     // Track quality cuts
388     //dtcpiptpc->SetStatus(AliESDtrack::kTPCrefit|AliESDtrack::kITSrefit);
389     dtcpiptpc->SetStatus(AliESDtrack::kTPCrefit);
390     dtcpiptpc->SetminTPCncls(50);
391     dtcpiptpc->SetRemoveKinks(kTRUE);
392     dtcpiptpc->SetLabel(kFALSE);
393 //     dtcpiptpc->SetMaxITSChiNdof(6.0);
394     dtcpiptpc->SetMaxTPCChiNdof(4.0);
395 //     dtcpiptpc->SetMaxImpactXY(3.0);
396 //     dtcpiptpc->SetMaxImpactZ(3.0);
397
398     AliFemtoCutMonitorParticleYPt *cutPassYPtpiptpc = new AliFemtoCutMonitorParticleYPt("cutPasspiptpc", 0.13957);
399     AliFemtoCutMonitorParticleYPt *cutFailYPtpiptpc = new AliFemtoCutMonitorParticleYPt("cutFailpiptpc", 0.13957);
400     dtcpiptpc->AddCutMonitor(cutPassYPtpiptpc, cutFailYPtpiptpc);
401
402     AliFemtoShareQualityTPCEntranceSepPairCut *sqpcpiptpc = new AliFemtoShareQualityTPCEntranceSepPairCut();
403     sqpcpiptpc->SetShareQualityMax(1.0);
404     sqpcpiptpc->SetShareFractionMax(1.0);
405     sqpcpiptpc->SetRemoveSameLabel(kFALSE);
406     sqpcpiptpc->SetTPCEntranceSepMinimum(0.0);
407
408     anpiptpc->SetEventCut(mecpiptpc);
409     anpiptpc->SetFirstParticleCut(dtcpiptpc);
410     anpiptpc->SetSecondParticleCut(dtcpiptpc);
411     anpiptpc->SetPairCut(sqpcpiptpc);
412     
413     AliFemtoShareQualityCorrFctn *csqqinvpiptpc= new AliFemtoShareQualityCorrFctn("sqqinvcfpiptpc",40,0.0,0.4);
414     anpiptpc->AddCorrFctn(csqqinvpiptpc);
415
416 //     AliFemtoCorrFctnDirectYlm *cylmpiptpc = new AliFemtoCorrFctnDirectYlm("cylmpiptpc",3,60,0.0,0.3,1);
417 //     anpiptpc->AddCorrFctn(cylmpiptpc);
418     
419     AliFemtoQinvCorrFctn *cqinvpiptpc = new AliFemtoQinvCorrFctn("qinvcfpiptpc", 100,0.0,1.0);
420     anpiptpc->AddCorrFctn(cqinvpiptpc);
421
422     AliFemtoChi2CorrFctn *cchiqinvpiptpc= new AliFemtoChi2CorrFctn("chicfpiptpc",40,0.0,0.4);
423     anpiptpc->AddCorrFctn(cchiqinvpiptpc);
424
425     AliFemtoCorrFctnTPCNcls *cqtpcnclspiptpc = new AliFemtoCorrFctnTPCNcls("cqtpcnclspiptpc",40,0.0,0.4);
426     anpiptpc->AddCorrFctn(cqtpcnclspiptpc);
427
428     Manager->AddAnalysis(anpiptpc);     
429   }
430   // *** End pion-pion (positive) analysis
431
432   // *** Begin pion-pion (negative) analysis ***
433   if (runNegativeTPCQA) {
434     AliFemtoVertexMultAnalysis *anpimtpc = new AliFemtoVertexMultAnalysis(1, zmin, zmax, 5, multmin, multmax);
435     anpimtpc->SetNumEventsToMix(iCollision?3:10);
436     anpimtpc->SetMinSizePartCollection(2);
437
438     AliFemtoBasicEventCut* mecpimtpc = new AliFemtoBasicEventCut();
439     mecpimtpc->SetEventMult(0,100000);
440     mecpimtpc->SetVertZPos(-1000,1000);
441         
442     AliFemtoESDTrackCut* dtcpimtpc = new AliFemtoESDTrackCut();
443 //     dtcpimtpc->SetPidProbPion(0.0,1.001);
444 //     dtcpimtpc->SetPidProbMuon(0.0,1.0);
445 //     dtcpimtpc->SetPidProbKaon(0.0,1.0);
446 //     dtcpimtpc->SetPidProbProton(0.0,1.0);
447     dtcpimtpc->SetCharge(0.0);
448     dtcpimtpc->SetPt(0.001,100.0);
449     dtcpimtpc->SetMass(PionMass);
450     // Track quality cuts
451     //dtcpimtpc->SetStatus(AliESDtrack::kTPCrefit|AliESDtrack::kITSrefit);
452     //    dtcpimtpc->SetStatus(AliESDtrack::kTPCrefit);
453 //     dtcpimtpc->SetminTPCncls(50);
454     dtcpimtpc->SetRemoveKinks(kTRUE);
455     dtcpimtpc->SetLabel(kFALSE);
456 //     dtcpimtpc->SetMaxITSChiNdof(6.0);
457 //     dtcpimtpc->SetMaxTPCChiNdof(6.0);
458     dtcpimtpc->SetMaxImpactXY(1.0);
459     dtcpimtpc->SetMaxImpactZ(1.0);
460
461     AliFemtoCutMonitorParticleYPt *cutPassYPtpimtpc = new AliFemtoCutMonitorParticleYPt("cutPasspimtpc", 0.13957);
462     AliFemtoCutMonitorParticleYPt *cutFailYPtpimtpc = new AliFemtoCutMonitorParticleYPt("cutFailpimtpc", 0.13957);
463     dtcpimtpc->AddCutMonitor(cutPassYPtpimtpc, cutFailYPtpimtpc);
464
465     AliFemtoShareQualityTPCEntranceSepPairCut *sqpcpimtpc = new AliFemtoShareQualityTPCEntranceSepPairCut();
466     sqpcpimtpc->SetShareQualityMax(1.0);
467     sqpcpimtpc->SetShareFractionMax(1.0);
468     sqpcpimtpc->SetRemoveSameLabel(kFALSE);
469     sqpcpimtpc->SetTPCEntranceSepMinimum(0.0);
470
471     anpimtpc->SetEventCut(mecpimtpc);
472     anpimtpc->SetFirstParticleCut(dtcpimtpc);
473     anpimtpc->SetSecondParticleCut(dtcpimtpc);
474     anpimtpc->SetPairCut(sqpcpimtpc);
475     
476     AliFemtoShareQualityCorrFctn *csqqinvpimtpc= new AliFemtoShareQualityCorrFctn("sqqinvcfpimtpc",40,0.0,0.4);
477     anpimtpc->AddCorrFctn(csqqinvpimtpc);
478
479     AliFemtoCorrFctnDirectYlm *cylmpimtpc = new AliFemtoCorrFctnDirectYlm("cylmpimtpc",3,60,0.0,0.3,1);
480     anpimtpc->AddCorrFctn(cylmpimtpc);
481     
482     AliFemtoQinvCorrFctn *cqinvpimtpc = new AliFemtoQinvCorrFctn("qinvcfpimtpc", 100,0.0,1.0);
483     anpimtpc->AddCorrFctn(cqinvpimtpc);
484
485     AliFemtoChi2CorrFctn *cchiqinvpimtpc= new AliFemtoChi2CorrFctn("chicfpimtpc",40,0.0,0.4);
486     anpimtpc->AddCorrFctn(cchiqinvpimtpc);
487
488     AliFemtoCorrFctnTPCNcls *cqtpcnclspimtpc = new AliFemtoCorrFctnTPCNcls("cqtpcnclspimtpc",40,0.0,0.4);
489     anpimtpc->AddCorrFctn(cqtpcnclspimtpc);
490
491     Manager->AddAnalysis(anpimtpc);     
492   }
493   // *** End pion-pion (negative) analysis
494
495
496
497   return Manager;
498 }                         
499