]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/STRANGENESS/Hypernuclei/AliAnalysisTaskHelium3Pi.cxx
Why the h*ll do we make a remote commit when pulling?
[u/mrichter/AliRoot.git] / PWGLF / STRANGENESS / Hypernuclei / AliAnalysisTaskHelium3Pi.cxx
1 /**************************************************************************\r
2  * Contributors are not mentioned at all.                                 *\r
3  *                                                                        *\r
4  * Permission to use, copy, modify and distribute this software and its   *\r
5  * documentation strictly for non-commercial purposes is hereby granted   *\r
6  * without fee, provided that the above copyright notice appears in all   *\r
7  * copies and that both the copyright notice and this permission notice   *\r
8  * appear in the supporting documentation. The authors make no claims     *\r
9  * about the suitability of this software for any purpose. It is          *\r
10  * provided "as is" without express or implied warranty.                  *\r
11  **************************************************************************/\r
12 //\r
13 //-----------------------------------------------------------------\r
14 //                 AliAnalysisTaskHelium3Pi class\r
15 //-----------------------------------------------------------------\r
16 \r
17 \r
18 class TTree;\r
19 class TParticle;\r
20 class TVector3;\r
21 \r
22 #include "AliAnalysisManager.h"\r
23 #include <AliMCEventHandler.h>\r
24 #include <AliMCEvent.h>\r
25 #include <AliStack.h>\r
26 \r
27 class AliESDVertex;\r
28 class AliAODVertex;\r
29 class AliESDv0;\r
30 class AliAODv0; \r
31 class AliCascadeVertexer;\r
32 \r
33 #include <iostream>\r
34 #include "AliAnalysisTaskSE.h"\r
35 #include "TList.h"\r
36 #include "TH1.h"\r
37 #include "TH2.h"\r
38 #include "TH3.h"\r
39 #include "TNtuple.h"\r
40 #include "TGraph.h"\r
41 #include "TCutG.h"\r
42 #include "TF1.h"\r
43 #include "TCanvas.h"\r
44 #include "TMath.h"\r
45 #include "TChain.h"\r
46 #include "Riostream.h"\r
47 #include "AliLog.h"\r
48 #include "AliCascadeVertexer.h"\r
49 #include "AliESDEvent.h"\r
50 #include "AliESDtrack.h"\r
51 #include "AliExternalTrackParam.h"\r
52 #include "AliAODEvent.h"\r
53 #include "AliInputEventHandler.h"\r
54 #include "AliESDcascade.h"\r
55 #include "AliAODcascade.h"\r
56 #include "AliAnalysisTaskHelium3Pi.h"\r
57 #include "AliESDtrackCuts.h"\r
58 #include "AliCentrality.h"\r
59 #include "TString.h"\r
60 #include <TDatime.h>\r
61 #include <TRandom3.h>\r
62 #include <TLorentzVector.h>\r
63 \r
64 const Int_t AliAnalysisTaskHelium3Pi::fgNrot = 15;\r
65 \r
66 ClassImp(AliAnalysisTaskHelium3Pi)\r
67 \r
68 //________________________________________________________________________\r
69 AliAnalysisTaskHelium3Pi::AliAnalysisTaskHelium3Pi() \r
70 : AliAnalysisTaskSE(),\r
71   fAnalysisType("ESD"), \r
72   fCollidingSystems(0), \r
73   fDataType("REAL"),\r
74   fListHistCascade(0), \r
75   fHistEventMultiplicity(0),         \r
76   fHistTrackMultiplicity(0),      \r
77   fHistTrackMultiplicityCent(0),      \r
78   fHistTrackMultiplicitySemiCent(0),  \r
79   fHistTrackMultiplicityMB(0),        \r
80   fHistTrackMultiplicityPVCent(0),      \r
81   fHistTrackMultiplicityPVSemiCent(0),  \r
82   fHistTrackMultiplicityPVMB(0),        \r
83   fHistMult(0),\r
84   fhBB(0),    \r
85   fhTOF(0),   \r
86   fhMassTOF(0),\r
87   fhBBPions(0),\r
88   fhBBHe(0),   \r
89   fhNaPos(0),  \r
90   fhNaNeg(0),  \r
91   fBetavsTPCsignalPos(0),  \r
92   fBetavsTPCsignalNeg(0),  \r
93   fHelium3TOF(0),   \r
94   fNtuple1(0),\r
95   fNtuple4(0)\r
96 \r
97 {\r
98   // Dummy Constructor \r
99 }\r
100 \r
101 //________________________________________________________________________\r
102 AliAnalysisTaskHelium3Pi::AliAnalysisTaskHelium3Pi(const char *name) \r
103   : AliAnalysisTaskSE(name), \r
104     fAnalysisType("ESD"), \r
105     fCollidingSystems(0), \r
106     fDataType("REAL"),\r
107     fListHistCascade(0), \r
108     fHistEventMultiplicity(0),    \r
109     fHistTrackMultiplicity(0),            \r
110     fHistTrackMultiplicityCent(0),      \r
111     fHistTrackMultiplicitySemiCent(0),  \r
112     fHistTrackMultiplicityMB(0),        \r
113     fHistTrackMultiplicityPVCent(0),      \r
114     fHistTrackMultiplicityPVSemiCent(0),  \r
115     fHistTrackMultiplicityPVMB(0),        \r
116     fHistMult(0),\r
117     fhBB(0),    \r
118     fhTOF(0),   \r
119     fhMassTOF(0),\r
120     fhBBPions(0),\r
121     fhBBHe(0),   \r
122     fhNaPos(0),  \r
123     fhNaNeg(0),  \r
124     fBetavsTPCsignalPos(0),  \r
125     fBetavsTPCsignalNeg(0),  \r
126     fHelium3TOF(0),                       \r
127     fNtuple1(0),\r
128     fNtuple4(0)  \r
129   \r
130 {\r
131   // Define input and output slots here\r
132   // Input slot #0 works with a TChain\r
133   //DefineInput(0, TChain::Class());\r
134   // Output slot #0 writes into a TList container (Cascade)\r
135 \r
136   DefineOutput(1, TList::Class());\r
137 }\r
138 //_______________________________________________________\r
139 AliAnalysisTaskHelium3Pi::~AliAnalysisTaskHelium3Pi() \r
140\r
141   // Destructor\r
142   if (fListHistCascade) {\r
143     delete fListHistCascade;\r
144     fListHistCascade = 0;\r
145   }\r
146 \r
147 }\r
148 //=================DEFINITION BETHE BLOCH==============================\r
149 \r
150 Double_t AliAnalysisTaskHelium3Pi::BetheBloch(Double_t betaGamma,Double_t charge,Bool_t isPbPb) {\r
151 \r
152   Double_t kp1, kp2, kp3, kp4, kp5;\r
153   \r
154   if(isPbPb){\r
155 \r
156     //    pass1 2011\r
157     kp1 = 4.7*charge*charge;\r
158     kp2 = 8.98482806165147636e+00;\r
159     kp3 = 1.54000000000000005e-05;\r
160     kp4 = 2.30445734159456084e+00;\r
161     kp5 = 2.25624744086878559e+00;\r
162 \r
163   }\r
164   \r
165   else{\r
166 \r
167     // to be defined ...\r
168     //pass1 2011\r
169     kp1 = 4.7*charge*charge;\r
170     kp2 = 8.98482806165147636e+00;\r
171     kp3 = 1.54000000000000005e-05;\r
172     kp4 = 2.30445734159456084e+00;\r
173     kp5 = 2.25624744086878559e+00;\r
174 \r
175   }\r
176 \r
177   Double_t beta = betaGamma / TMath::Sqrt(1.0 + betaGamma * betaGamma);\r
178   \r
179   Double_t aa = TMath::Power(beta, kp4);\r
180   Double_t bb = TMath::Power(1.0 / betaGamma, kp5);\r
181   \r
182   bb = TMath::Log(kp3 + bb);\r
183   \r
184   Double_t out = (kp2 - aa - bb) * kp1 / aa;\r
185 \r
186   return out;\r
187  \r
188 }\r
189 \r
190 //==================DEFINITION OF OUTPUT OBJECTS==============================\r
191 \r
192 void AliAnalysisTaskHelium3Pi::UserCreateOutputObjects()\r
193 {\r
194 \r
195   fListHistCascade = new TList();\r
196   fListHistCascade->SetOwner();  // IMPORTANT!\r
197 \r
198   if(! fHistEventMultiplicity ){\r
199     fHistEventMultiplicity   = new TH1F( "fHistEventMultiplicity" , "Nb of Events" , 10 , 0, 10);\r
200     fHistEventMultiplicity->GetXaxis()->SetBinLabel(1,"All Events");\r
201     fHistEventMultiplicity->GetXaxis()->SetBinLabel(2,"Events w/PV");\r
202     fHistEventMultiplicity->GetXaxis()->SetBinLabel(3,"Events w/|Vz|<10cm");\r
203     fHistEventMultiplicity->GetXaxis()->SetBinLabel(4,"Central Events");\r
204     fHistEventMultiplicity->GetXaxis()->SetBinLabel(5,"SemiCentral Events");\r
205     fHistEventMultiplicity->GetXaxis()->SetBinLabel(6,"MB Events");\r
206     fHistEventMultiplicity->GetXaxis()->SetBinLabel(7,"Central Events  w/|Vz|<10cm");\r
207     fHistEventMultiplicity->GetXaxis()->SetBinLabel(8,"SemiCentral Events  w/|Vz|<10cm");\r
208     fHistEventMultiplicity->GetXaxis()->SetBinLabel(9,"MB Events   w/|Vz|<10cm");\r
209 \r
210     fListHistCascade->Add(fHistEventMultiplicity);\r
211   }\r
212 \r
213   if(! fHistTrackMultiplicity ){\r
214     fHistTrackMultiplicity   = new TH2F( "fHistTrackMultiplicity" , "Nb of Tracks", 25000,0, 25000,105,-1,104);\r
215     fHistTrackMultiplicity->GetXaxis()->SetTitle("Number of tracks");\r
216     fHistTrackMultiplicity->GetYaxis()->SetTitle("Percentile");\r
217     fListHistCascade->Add(fHistTrackMultiplicity);\r
218   } \r
219 \r
220   if(! fHistTrackMultiplicityCent ){\r
221     fHistTrackMultiplicityCent   = new TH2F( "fHistTrackMultiplicityCent", "Nb of Tracks Central Events", 25000,0, 25000,105,-1,104 );\r
222     fHistTrackMultiplicityCent->GetXaxis()->SetTitle("Number of tracks");\r
223     fHistTrackMultiplicityCent->GetYaxis()->SetTitle("Percentile");\r
224     fListHistCascade->Add(fHistTrackMultiplicityCent);\r
225   } \r
226 \r
227   if(! fHistTrackMultiplicitySemiCent ){\r
228     fHistTrackMultiplicitySemiCent   = new TH2F( "fHistTrackMultiplicitySemiCent" , "Nb of Tracks SemiCentral Events", 25000,0, 25000 ,105,-1,104);\r
229     fHistTrackMultiplicitySemiCent->GetXaxis()->SetTitle("Number of tracks");\r
230     fHistTrackMultiplicitySemiCent->GetYaxis()->SetTitle("Percentile");\r
231     fListHistCascade->Add(fHistTrackMultiplicitySemiCent);\r
232   } \r
233  \r
234   if(! fHistTrackMultiplicityMB ){\r
235     fHistTrackMultiplicityMB   = new TH2F( "fHistTrackMultiplicityMB" , "Nb of Tracks MBral Events", 25000,0, 25000,105,-1,104 );\r
236     fHistTrackMultiplicityMB->GetXaxis()->SetTitle("Number of tracks");\r
237     fHistTrackMultiplicityMB->GetYaxis()->SetTitle("Percentile");\r
238     fListHistCascade->Add(fHistTrackMultiplicityMB);\r
239   } \r
240 \r
241   if(! fHistTrackMultiplicityPVCent ){\r
242     fHistTrackMultiplicityPVCent   = new TH2F( "fHistTrackMultiplicityPVCent" , "Nb of Tracks Central Events", 25000,0, 25000,105,-1,104 );\r
243     fHistTrackMultiplicityPVCent->GetXaxis()->SetTitle("Number of tracks");\r
244     fHistTrackMultiplicityPVCent->GetYaxis()->SetTitle("Percentile");\r
245     fListHistCascade->Add(fHistTrackMultiplicityPVCent);\r
246   } \r
247 \r
248   if(! fHistTrackMultiplicityPVSemiCent ){\r
249     fHistTrackMultiplicityPVSemiCent   = new TH2F( "fHistTrackMultiplicityPVSemiCent" , "Nb of Tracks SemiCentral Events", 25000,0, 25000 ,105,-1,104);\r
250     fHistTrackMultiplicityPVSemiCent->GetXaxis()->SetTitle("Number of tracks");\r
251     fHistTrackMultiplicityPVSemiCent->GetYaxis()->SetTitle("Percentile");\r
252     fListHistCascade->Add(fHistTrackMultiplicityPVSemiCent);\r
253   } \r
254  \r
255   if(! fHistTrackMultiplicityPVMB ){\r
256     fHistTrackMultiplicityPVMB   = new TH2F( "fHistTrackMultiplicityPVMB" , "Nb of Tracks MBral Events", 25000,0, 25000,105,-1,104 );\r
257     fHistTrackMultiplicityPVMB->GetXaxis()->SetTitle("Number of tracks");\r
258     fHistTrackMultiplicityPVMB->GetYaxis()->SetTitle("Percentile");\r
259     fListHistCascade->Add(fHistTrackMultiplicityPVMB);\r
260   } \r
261 \r
262   if(! fHistMult){\r
263     fHistMult=new TH1F ("fHistMult","Number neg-pos", 10, -1,9);\r
264     fHistMult->GetXaxis()->SetTitle("Type of tracks");\r
265     fListHistCascade->Add(fHistMult);\r
266   }\r
267  \r
268   if(! fhBB ){\r
269     fhBB = new TH2F( "fhBB" , "BetheBlochTPC" , 1000,-6,6,1500,0,5000);\r
270     fhBB->GetXaxis()->SetTitle("p/z (GeV/#it{c})");\r
271     fhBB->GetYaxis()->SetTitle("TPC Signal");\r
272     fListHistCascade->Add(fhBB);\r
273   }\r
274 \r
275   if(! fhTOF ){\r
276     fhTOF = new TH2F( "fhTOF" , "Scatter Plot TOF" , 1000,-6,6,1000,0,1.2);\r
277     fhTOF->GetXaxis()->SetTitle("p/z (GeV/#it{c})");\r
278     fhTOF->GetYaxis()->SetTitle("#beta");\r
279     fListHistCascade->Add(fhTOF);\r
280   }\r
281 \r
282   if(! fhMassTOF){\r
283     fhMassTOF=new TH1F ("fhMassTOF","Particle Mass - TOF", 300,0 ,5);\r
284     fhMassTOF->GetXaxis()->SetTitle("Mass (GeV/#it{c}^{2})");\r
285     fListHistCascade->Add(fhMassTOF);\r
286   }\r
287 \r
288   if(! fhBBPions ){\r
289     fhBBPions = new TH2F( "fhBBPions" , "Bethe-Bloch TPC Pions" , 1000,-6,6,1500,0,5000);\r
290     fhBBPions->GetXaxis()->SetTitle("p/z (GeV/#it{c})");\r
291     fhBBPions->GetYaxis()->SetTitle("TPC Signal");\r
292     fListHistCascade->Add(fhBBPions);\r
293   }\r
294   \r
295   if(! fhBBHe ){\r
296     fhBBHe = new TH2F( "fhBBHe" , "Bethe-Bloch TPC He" , 1000,-6,6,1500,0,5000);\r
297     fhBBHe->GetXaxis()->SetTitle("p/z (GeV/#it{c})");\r
298     fhBBHe->GetYaxis()->SetTitle("TPC Signal");\r
299     fListHistCascade->Add(fhBBHe);\r
300   }\r
301   \r
302   if(! fhNaPos ){\r
303     fhNaPos = new TH2F( "fhNaPos" , "Distribution Pos" , 500,0,5,500,-2,2);\r
304     fhNaPos->GetXaxis()->SetTitle("p/z (GeV/#it{c})");\r
305     fhNaPos->GetYaxis()->SetTitle("(TPCSignal-bbtheo)/bbtheo (He)");\r
306     fListHistCascade->Add(fhNaPos);\r
307   }\r
308   \r
309   if(! fhNaNeg ){\r
310     fhNaNeg = new TH2F( "fhNaNeg" , "Distribution Neg" , 500,0,5,500,-2,2);\r
311     fhNaNeg->GetXaxis()->SetTitle("p/z (GeV/#it{c})");\r
312     fhNaNeg->GetYaxis()->SetTitle("(TPCSignal-bbtheo)/bbtheo (He)");\r
313     fListHistCascade->Add(fhNaNeg);\r
314   }\r
315 \r
316   if(! fBetavsTPCsignalPos ){\r
317     fBetavsTPCsignalPos = new TH2F("fBetavsTPCsignalPos","fBetavsTPCsignalPos",1000,0,1.2,1500,0,5000);\r
318     fBetavsTPCsignalPos->GetXaxis()->SetTitle("#beta");\r
319     fBetavsTPCsignalPos->GetYaxis()->SetTitle("TPC Signal");\r
320     fListHistCascade->Add(fBetavsTPCsignalPos);\r
321   }\r
322   \r
323   if(! fBetavsTPCsignalNeg ){\r
324     fBetavsTPCsignalNeg = new TH2F("fBetavsTPCsignalNeg","fBetavsTPCsignalNeg",1000,0,1.2,1500,0,5000);\r
325     fBetavsTPCsignalNeg->GetXaxis()->SetTitle("#beta");\r
326     fBetavsTPCsignalNeg->GetYaxis()->SetTitle("TPC Signal");\r
327     fListHistCascade->Add(fBetavsTPCsignalNeg);\r
328   }\r
329   \r
330   if(! fHelium3TOF){\r
331     fHelium3TOF = new TH2F("fHelium3massTOF","Helium3 beta vs p/z",1000,0,6,1000,0,1.2);\r
332     fHelium3TOF->GetXaxis()->SetTitle("p/z (GeV/#it{c})");\r
333     fHelium3TOF->GetYaxis()->SetTitle("#beta");\r
334     fListHistCascade->Add(fHelium3TOF);\r
335   }\r
336   \r
337   if(! fNtuple1 ) {\r
338     fNtuple1 = new TNtuple("fNtuple1","Ntuple1","runNumber:bunchcross:orbit:period:eventtype:TrackNumber:percentile:xPrimaryVertex:yPrimaryVertex:zPrimaryVertex:xSecondaryVertex:ySecondaryVertex:zSecondaryVertex:dcaTracks:CosPointingAngle:DCAV0toPrimaryVertex:HeSign:HepInTPC:HeTPCsignal:DcaHeToPrimVertex:HeEta:momHex:momHey:momHez:momHeAtSVx:momHeAtSVy:momHeAtSVz:HeTPCNcls:HeimpactXY:HeimpactZ:HeITSClusterMap:IsHeITSRefit:PionSign:PionpInTPC:PionTPCsignal:DcaPionToPrimVertex:PionEta:momPionx:momPiony:momPionz:momNegPionAtSVx:momNegPionAtSVy:momNegPionAtSVz:PionTPCNcls:PionimpactXY:PionimpactZ:PionITSClusterMap:IsPiITSRefit:xn:xp");\r
339     \r
340     fListHistCascade->Add(fNtuple1);\r
341   }\r
342   \r
343   if(! fNtuple4 ) {\r
344     fNtuple4 = new TNtuple("fNtuple4","Ntuple4","runNumber:BCNumber:OrbitNumber:PeriodNumber:eventtype:isHeITSrefit:percentile:Sign:pinTPC:GetTPCsignal:Px:Py:Pz:Eta:isTOF:poutTPC:timeTOF:trackLenghtTOF:impactXY:impactZ:mapITS:TPCNcls:TRDsignal:xPrimaryVertex:yPrimaryVertex:zPrimaryVertex:chi2PerClusterTPC");\r
345     fListHistCascade->Add(fNtuple4);\r
346   } \r
347 \r
348   PostData(1,  fListHistCascade);\r
349 \r
350 }// end UserCreateOutputObjects\r
351 \r
352 \r
353 \r
354 //====================== USER EXEC ========================\r
355 \r
356 void AliAnalysisTaskHelium3Pi::UserExec(Option_t *) \r
357 {\r
358   //_______________________________________________________________________\r
359   \r
360   //!*********************!//\r
361   //!  Define variables   !//\r
362   //!*********************!//\r
363   Float_t vett1[50];\r
364   for(Int_t i=0;i<50;i++) vett1[i]=0;\r
365 \r
366   Float_t vett4[40];\r
367   for(Int_t i=0;i<40;i++) vett4[i]=0;\r
368   \r
369   Double_t ITSsample[4];\r
370   for(Int_t i=0;i<4;i++)ITSsample[i]=0;\r
371 \r
372   Double_t  pinTPC=0.,poutTPC=0.,TPCSignal=0.;\r
373   Double_t xPrimaryVertex=0.,yPrimaryVertex=0.,zPrimaryVertex=0.;\r
374   Double_t massTOF=0.,timeTOF=0.,trackLenghtTOF=0.,betaTOF=0.;\r
375 \r
376   ULong_t  status=0;\r
377   ULong_t  statusT=0;\r
378   ULong_t  statusPi=0;\r
379 \r
380   Bool_t   isTPC=kFALSE,isTOF=kFALSE,isTOFHe=kFALSE,IsHeITSRefit=kFALSE,isTOFPi=kFALSE,IsPiITSRefit=kFALSE ;\r
381 \r
382   Float_t nSigmaNegPion=0.;\r
383   Float_t nSigmaNegPion1=0.;\r
384   Float_t nSigmaNegPion2=0.;\r
385 \r
386   Double_t cutNSigma = 3;\r
387   Double_t bbtheoM=0.,bbtheoP=0.,bbtheo=0.;\r
388   Double_t zNathashaNeg=0;\r
389   Double_t zNathashaPos=0;\r
390   Double_t fPos[3]={0.,0.,0.};\r
391   Double_t runNumber=0.;\r
392   Double_t evNumber=0.;\r
393 \r
394   Double_t BCNumber=0.;\r
395   Double_t OrbitNumber=0.;\r
396   Double_t PeriodNumber=0.;\r
397 \r
398   Double_t        Helium3Mass = 2.80839; \r
399   Double_t        PionMass    = 0.13957; \r
400   // TLORENTZ vectors\r
401   \r
402   TLorentzVector  vPion,vHelium,vSum;\r
403 \r
404   //!----------------------------------------------------------------\r
405 \r
406   //! A set of very loose parameters for cuts \r
407   \r
408   Double_t fgChi2max=33.;     //! max chi2\r
409   Double_t fgDNmin=0.05;      //! min imp parameter for the 1st daughter = 500um\r
410   Double_t fgDCAmax=1.0;      //! max DCA between the daughter tracks in cm\r
411   Double_t fgCPAmin=0.99;     //! min cosine of V0's pointing angle  \r
412   //  Double_t fgRmin=0.2;    //! min radius of the fiducial volume //original\r
413   Double_t fgRmin=0.1;        //! min radius of the fiducial volume = 1 mm \r
414   Double_t fgRmax=200.;       //! max radius of the fiducial volume = 2 m\r
415 \r
416   //------------------------------------------\r
417 \r
418   // Main loop\r
419   // Called for EACH event\r
420 \r
421   AliVEvent *event = InputEvent();\r
422   if (!event) { Printf("ERROR: Could not retrieve event"); return; }\r
423     \r
424   Info("AliAnalysisTaskHelium3Pi","Starting UserExec");  \r
425 \r
426   SetDataType("REAL");\r
427   \r
428   // create pointer to event\r
429   AliESDEvent* lESDevent = dynamic_cast<AliESDEvent*>(event);\r
430   if (!lESDevent) {\r
431     AliError("Cannot get the ESD event");\r
432     return;\r
433   }  \r
434 \r
435   fHistEventMultiplicity->Fill(0);\r
436 \r
437   Double_t lMagneticField=lESDevent->GetMagneticField();\r
438   Int_t TrackNumber = -1;\r
439 \r
440    \r
441   //*****************//  \r
442   //*   Centrality  *//\r
443   //*****************//\r
444  \r
445   AliCentrality *centrality = lESDevent->GetCentrality();\r
446   Float_t percentile=centrality->GetCentralityPercentile("V0M");\r
447 \r
448   TrackNumber = lESDevent->GetNumberOfTracks();\r
449   if (TrackNumber<2) return;  \r
450 \r
451   fHistTrackMultiplicity->Fill(TrackNumber,percentile); //tracce per evento\r
452 \r
453   //****************************************\r
454 \r
455   Int_t eventtype=-99;\r
456   \r
457   Bool_t isSelectedCentral     = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & AliVEvent::kCentral);\r
458   Bool_t isSelectedSemiCentral = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & AliVEvent::kSemiCentral);\r
459   Bool_t isSelectedMB          = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & AliVEvent::kMB);\r
460  \r
461   if(isSelectedCentral){\r
462     fHistEventMultiplicity->Fill(3);\r
463     fHistTrackMultiplicityCent->Fill(TrackNumber,percentile); \r
464     eventtype=1;\r
465   }\r
466 \r
467   if(isSelectedSemiCentral){\r
468     fHistEventMultiplicity->Fill(4);\r
469     fHistTrackMultiplicitySemiCent->Fill(TrackNumber,percentile); \r
470     eventtype=2;\r
471   }\r
472 \r
473   if(isSelectedMB){\r
474     fHistEventMultiplicity->Fill(5);\r
475     fHistTrackMultiplicityMB->Fill(TrackNumber,percentile); \r
476     eventtype=3;\r
477   }\r
478   \r
479   if(isSelectedCentral || isSelectedSemiCentral || isSelectedMB){\r
480     \r
481     // ANALISYS\r
482     \r
483     // Primary vertex cut\r
484     \r
485     const AliESDVertex *vtx = lESDevent->GetPrimaryVertexTracks();\r
486     \r
487     if(vtx->GetNContributors()<1) {\r
488       \r
489       // SPD vertex cut\r
490       vtx = lESDevent->GetPrimaryVertexSPD();\r
491       \r
492       if(vtx->GetNContributors()<1) {\r
493         Info("AliAnalysisTaskHelium3Pi","No good vertex, skip event");\r
494         return; // NO GOOD VERTEX, SKIP EVENT \r
495       }\r
496     }\r
497     \r
498     fHistEventMultiplicity->Fill(1); // analyzed events with PV\r
499  \r
500     xPrimaryVertex=vtx->GetXv();\r
501     yPrimaryVertex=vtx->GetYv();\r
502     zPrimaryVertex=vtx->GetZv();  \r
503     \r
504     if(TMath::Abs(zPrimaryVertex)>10) return;\r
505     \r
506     if(eventtype==1){\r
507       fHistTrackMultiplicityPVCent->Fill(TrackNumber,percentile); \r
508       fHistEventMultiplicity->Fill(6); \r
509     }\r
510     \r
511     if(eventtype==2){\r
512       fHistTrackMultiplicityPVSemiCent->Fill(TrackNumber,percentile); \r
513       fHistEventMultiplicity->Fill(7); \r
514     }\r
515     \r
516     if(eventtype==3){\r
517       fHistTrackMultiplicityPVMB->Fill(TrackNumber,percentile); \r
518       fHistEventMultiplicity->Fill(8); \r
519     }\r
520     \r
521     \r
522     fHistEventMultiplicity->Fill(2);\r
523     \r
524     //Find Pair candidates\r
525     \r
526     TArrayI PionsTPC(TrackNumber);        //Neg pions\r
527     Int_t nPionsTPC=0;\r
528     \r
529     TArrayI HeTPC(TrackNumber);        //helium3\r
530     Int_t nHeTPC=0;\r
531     \r
532     const Double_t speedOfLight =  TMath::C()*1E2*1E-12; // cm/ps\r
533     \r
534     Float_t impactXY=-999, impactZ=-999;\r
535     Float_t impactXYpi=-999, impactZpi=-999;\r
536     \r
537     //!   SELECTIONS:\r
538     //! - No ITSpureSA\r
539     //! - ITSrefit\r
540     //! - TPCrefit\r
541     //! - ITSmap !=0\r
542     \r
543     // ******************* Track Cuts Definitions ********************//\r
544   \r
545     //! ITS\r
546 \r
547     AliESDtrackCuts* esdtrackCutsITS = new AliESDtrackCuts("esdtrackCutsITS");\r
548     esdtrackCutsITS->SetRequireITSStandAlone(kFALSE);\r
549     esdtrackCutsITS->SetRequireITSPureStandAlone(kFALSE);\r
550     \r
551     //! TPC\r
552     \r
553     Int_t    minclsTPC=60;\r
554     Double_t maxchi2perTPCcl=5.;\r
555     \r
556     AliESDtrackCuts* esdtrackCutsTPC = new AliESDtrackCuts("esdtrackCutsTPC");\r
557     esdtrackCutsTPC->SetRequireTPCRefit(kTRUE);\r
558     esdtrackCutsTPC->SetAcceptKinkDaughters(kFALSE);\r
559     esdtrackCutsTPC->SetMinNClustersTPC(minclsTPC);\r
560     esdtrackCutsTPC->SetMaxChi2PerClusterTPC(maxchi2perTPCcl);\r
561     \r
562     //*********************************************+\r
563     \r
564     runNumber = lESDevent->GetRunNumber();\r
565     evNumber  = lESDevent->GetEventNumberInFile();\r
566     \r
567     BCNumber    = lESDevent->GetBunchCrossNumber();\r
568     OrbitNumber = lESDevent->GetOrbitNumber();\r
569     PeriodNumber= lESDevent->GetPeriodNumber();\r
570     \r
571     //*************************************************************\r
572     \r
573     TF1 *foPion=new TF1("foPion", "[0]*([1]*TMath::Power(TMath::Sqrt(1 + (x/0.13957)*(x/0.13957))/(x/0.13957) , [3]) - 1 - TMath::Power(TMath::Sqrt(1 + (x/0.13957)*(x/0.13957))/(x/0.13957) , [3])*TMath::Log([2] + 1/TMath::Power((x/0.13957), [4])))",0.01,20);\r
574     foPion->SetParameters(4.1,8.98482806165147636e+00,1.54000000000000005e-05,2.30445734159456084e+00,2.25624744086878559e+00);\r
575     \r
576     //--------------------------------------------\r
577     \r
578     for (Int_t j=0; j<TrackNumber; j++) { //loop on tracks\r
579       \r
580       AliESDtrack *esdtrack=lESDevent->GetTrack(j);\r
581       \r
582       if(!esdtrack) { \r
583         AliError(Form("ERROR: Could not retrieve esdtrack %d",j)); \r
584         continue; \r
585       }\r
586       \r
587       // ************** Track cuts ****************\r
588       \r
589       status  = (ULong_t)esdtrack->GetStatus();\r
590       \r
591       isTPC   = (((status) & AliESDtrack::kTPCin)  != 0);\r
592       isTOF   = ((((status) & AliESDtrack::kTOFout) != 0) && (((status) & AliESDtrack::kTIME) != 0));\r
593       \r
594       Bool_t IsTrackAcceptedTPC =  esdtrackCutsTPC->AcceptTrack(esdtrack);\r
595       Bool_t IsTrackAcceptedITS =  esdtrackCutsITS->AcceptTrack(esdtrack);\r
596       \r
597       if (!(IsTrackAcceptedTPC && IsTrackAcceptedITS)) continue;\r
598       \r
599       UInt_t mapITS=esdtrack->GetITSClusterMap();\r
600       //    if(mapITS==0) continue;\r
601       \r
602       //----------------------------------------------\r
603       \r
604       //****** Cuts from  AliV0Vertex.cxx *************\r
605       \r
606       Double_t d=esdtrack->GetD(xPrimaryVertex,yPrimaryVertex,lMagneticField);\r
607       //    if (TMath::Abs(d)<fgDPmin) continue;\r
608       if (TMath::Abs(d)>fgRmax) continue;\r
609       \r
610       //---- (Usefull) Stuff\r
611       \r
612       TPCSignal=esdtrack->GetTPCsignal(); \r
613       \r
614       if (TPCSignal<10)continue;\r
615       \r
616       if(!isTPC)continue;\r
617       \r
618       if(!esdtrack->GetTPCInnerParam())continue;\r
619       \r
620       AliExternalTrackParam trackIn(*esdtrack->GetInnerParam()); \r
621       pinTPC = trackIn.GetP(); \r
622       \r
623       poutTPC=pinTPC;\r
624       \r
625       fHistMult->Fill(0);\r
626       \r
627       if((status) & (AliESDtrack::kITSrefit!=0)){\r
628         fHistMult->Fill(1);\r
629         fhBB->Fill(pinTPC*esdtrack->GetSign(),TPCSignal);\r
630       }\r
631       \r
632       timeTOF=esdtrack->GetTOFsignal();                 // ps\r
633       trackLenghtTOF= esdtrack->GetIntegratedLength();  // cm\r
634       \r
635       if(isTOF){\r
636         \r
637         if(!esdtrack->GetOuterParam())continue;    \r
638         \r
639         AliExternalTrackParam trackOut(*esdtrack->GetOuterParam()); \r
640         \r
641         poutTPC = trackOut.GetP(); \r
642         \r
643         betaTOF= (trackLenghtTOF/timeTOF)/2.99792458e-2;\r
644         \r
645         fhTOF->Fill(poutTPC*esdtrack->GetSign(),betaTOF);\r
646         \r
647         Double_t mass2=(poutTPC*poutTPC)*((((speedOfLight*speedOfLight)*(timeTOF*timeTOF))-(trackLenghtTOF*trackLenghtTOF))/(trackLenghtTOF*trackLenghtTOF));\r
648         if(mass2>0) massTOF=TMath::Sqrt(mass2);\r
649         fhMassTOF->Fill(massTOF);\r
650         \r
651         if(esdtrack->GetSign() < 0.)fBetavsTPCsignalNeg->Fill(betaTOF,TPCSignal);\r
652         if(esdtrack->GetSign() > 0.)fBetavsTPCsignalPos->Fill(betaTOF,TPCSignal);\r
653         \r
654       }\r
655       \r
656       //pass2\r
657       \r
658       bbtheo =BetheBloch((2*pinTPC)/3.,2,kTRUE);    //! OK\r
659       bbtheoM=(1 - 0.08*5)*bbtheo;                  //! OK \r
660       bbtheoP=(1 + 0.08*5)*bbtheo;                  //! OK\r
661       \r
662       fHistMult->Fill(2);\r
663       \r
664       if(esdtrack->GetSign()<0){\r
665         zNathashaNeg=(TPCSignal-bbtheo)/bbtheo;\r
666         fhNaNeg->Fill(pinTPC,zNathashaNeg); \r
667       }\r
668       \r
669       if(esdtrack->GetSign() > 0.){\r
670         zNathashaPos=(TPCSignal-bbtheo)/bbtheo;\r
671         fhNaPos->Fill(pinTPC,zNathashaPos); \r
672       }\r
673       \r
674       nSigmaNegPion1 = TMath::Abs((TPCSignal - foPion->Eval(pinTPC))/foPion->Eval(pinTPC))/0.07; \r
675       nSigmaNegPion2 = TMath::Abs((TPCSignal - foPion->Eval(pinTPC))/foPion->Eval(pinTPC))/0.04; \r
676       \r
677       if(pinTPC<0.6)\r
678         nSigmaNegPion=nSigmaNegPion1;\r
679       if(pinTPC>=0.6)\r
680         nSigmaNegPion=nSigmaNegPion2;\r
681       \r
682       if ( (nSigmaNegPion < cutNSigma)){ \r
683         \r
684         fhBBPions->Fill(pinTPC*esdtrack->GetSign(),TPCSignal);\r
685         \r
686         if(pinTPC<3){\r
687           PionsTPC[nPionsTPC++]=j;\r
688         }\r
689       }\r
690       \r
691       if( TPCSignal > bbtheoM ) {\r
692         \r
693         if(pinTPC>0.6){\r
694           \r
695           fhBBHe->Fill(pinTPC*esdtrack->GetSign(),TPCSignal);\r
696           HeTPC[nHeTPC++]=j;\r
697           \r
698           \r
699           Bool_t isHeITSrefit=((status) & (AliESDtrack::kITSrefit));\r
700           \r
701           esdtrack->GetImpactParameters(impactXY, impactZ);\r
702           esdtrack->GetITSdEdxSamples(ITSsample);\r
703           \r
704           \r
705           Int_t  fIdxInt[200]; //dummy array\r
706           Int_t nClustersTPC = esdtrack->GetTPCclusters(fIdxInt);\r
707           \r
708           Float_t chi2PerClusterTPC = esdtrack->GetTPCchi2()/(Float_t)(nClustersTPC);\r
709           \r
710           vett4[0] =(Float_t)runNumber;\r
711           vett4[1] =(Float_t)BCNumber;\r
712           vett4[2] =(Float_t)OrbitNumber;\r
713           vett4[3] =(Float_t)PeriodNumber;\r
714           vett4[4] =(Float_t)eventtype;\r
715           vett4[5] =(Float_t)isHeITSrefit;\r
716           vett4[6] =(Float_t)percentile;\r
717           vett4[7] =(Float_t)esdtrack->GetSign();\r
718           vett4[8] =(Float_t)pinTPC;\r
719           vett4[9] =(Float_t)esdtrack->GetTPCsignal();\r
720           vett4[10]=(Float_t)esdtrack->Px();\r
721           vett4[11]=(Float_t)esdtrack->Py();\r
722           vett4[12]=(Float_t)esdtrack->Pz();\r
723           vett4[13]=(Float_t)esdtrack->Eta();\r
724           vett4[14]=(Float_t)isTOF;\r
725           vett4[15]=(Float_t)poutTPC;\r
726           vett4[16]=(Float_t)timeTOF;\r
727           vett4[17]=(Float_t)trackLenghtTOF;\r
728           vett4[18]=(Float_t)impactXY;\r
729           vett4[19]=(Float_t)impactZ;\r
730           vett4[20]=(Float_t)mapITS;\r
731           vett4[21]=(Float_t)esdtrack->GetTPCNcls();\r
732           vett4[22]=(Float_t)esdtrack->GetTRDsignal();\r
733           vett4[23]=(Float_t)xPrimaryVertex;\r
734           vett4[24]=(Float_t)yPrimaryVertex;\r
735           vett4[25]=(Float_t)zPrimaryVertex;\r
736           vett4[26]=(Float_t)chi2PerClusterTPC;\r
737           \r
738           fNtuple4->Fill(vett4);\r
739         }\r
740       }\r
741     }       //! track\r
742     \r
743     \r
744     Double_t        DcaHeToPrimVertex=0;\r
745     Double_t        DcaPionToPrimVertex=0;\r
746     \r
747     impactXY=-999, impactZ=-999;\r
748     impactXYpi=-999, impactZpi=-999;\r
749     \r
750     // Track \r
751     \r
752     AliESDtrack  *PionTrack = 0x0;\r
753     AliESDtrack  *HeTrack = 0x0;\r
754     \r
755     // Vettors for il PxPyPz\r
756     \r
757     Double_t momPionVett[3];\r
758     for(Int_t i=0;i<3;i++)momPionVett[i]=0;\r
759     \r
760     Double_t momHeVett[3];\r
761     for(Int_t i=0;i<3;i++)momHeVett[i]=0;\r
762     \r
763     //At SV\r
764     \r
765     Double_t momPionVettAt[3];\r
766     for(Int_t i=0;i<3;i++)momPionVettAt[i]=0;\r
767     \r
768     Double_t momHeVettAt[3];\r
769     for(Int_t i=0;i<3;i++)momHeVettAt[i]=0;\r
770     \r
771     //---------------   LOOP PAIRS   ----------------\r
772     \r
773     for (Int_t k=0; k < nPionsTPC; k++) {                           //! Pions Loop\r
774       \r
775       DcaPionToPrimVertex=0.;\r
776       DcaHeToPrimVertex=0;\r
777       \r
778       Int_t PionIdx=PionsTPC[k];\r
779       \r
780       PionTrack=lESDevent->GetTrack(PionIdx);\r
781       \r
782       statusPi = (ULong_t)PionTrack->GetStatus();\r
783       isTOFPi  = ((((statusPi) & (AliESDtrack::kTOFout)) != 0) && (((statusPi) & (AliESDtrack::kTIME)) != 0));\r
784       IsPiITSRefit = ((statusPi) & (AliESDtrack::kITSrefit)); \r
785       \r
786       if (PionTrack) \r
787         DcaPionToPrimVertex = TMath::Abs(PionTrack->GetD(xPrimaryVertex, yPrimaryVertex,lMagneticField)); //OK\r
788       \r
789       if(DcaPionToPrimVertex<0.2)continue; //qui\r
790       \r
791       AliExternalTrackParam trackInPion(*PionTrack);  \r
792       \r
793       for (Int_t i=0; i<nHeTPC; i++){                               //! Helium Loop\r
794         \r
795         Int_t HeIdx=HeTPC[i];\r
796         \r
797         HeTrack=lESDevent->GetTrack(HeIdx);\r
798         \r
799         statusT= (ULong_t)HeTrack->GetStatus();\r
800         isTOFHe   = (((statusT & AliESDtrack::kTOFout) != 0) && ((statusT & AliESDtrack::kTIME) != 0));\r
801         IsHeITSRefit = (status & AliESDtrack::kITSrefit); \r
802         \r
803         if (HeTrack) \r
804           DcaHeToPrimVertex = TMath::Abs(HeTrack->GetD(xPrimaryVertex, yPrimaryVertex,lMagneticField)); //OK\r
805         \r
806         AliExternalTrackParam trackInHe(*HeTrack); \r
807     \r
808         if ( DcaPionToPrimVertex < fgDNmin)                //OK\r
809           if ( DcaHeToPrimVertex < fgDNmin) continue;    //OK\r
810         \r
811         Double_t xn, xp;\r
812         Double_t dca=0.;\r
813         \r
814         dca= PionTrack->GetDCA(HeTrack,lMagneticField,xn,xp); //!dca (Neg to Pos)\r
815         \r
816         if (dca > fgDCAmax) continue;\r
817         if ((xn+xp) > 2*fgRmax) continue;\r
818         if ((xn+xp) < 2*fgRmin) continue;\r
819         \r
820         //CORRECTION from AliV0Vertex\r
821         \r
822         Bool_t corrected=kFALSE;\r
823         if ((trackInPion.GetX() > 3.) && (xn < 3.)) {\r
824           //correct for the beam pipe material\r
825           corrected=kTRUE;\r
826         }\r
827         if ((trackInHe.GetX() > 3.) && (xp < 3.)) {\r
828           //correct for the beam pipe material\r
829           corrected=kTRUE;\r
830         }\r
831         if (corrected) {\r
832           dca=trackInPion.GetDCA(&trackInHe,lMagneticField,xn,xp);\r
833           if (dca > fgDCAmax) continue;\r
834           if ((xn+xp) > 2*fgRmax) continue;\r
835           if ((xn+xp) < 2*fgRmin) continue;\r
836         }\r
837         \r
838         //=============================================//\r
839         // Make "V0" with found tracks                 //\r
840         //=============================================//\r
841         \r
842         trackInPion.PropagateTo(xn,lMagneticField); \r
843         trackInHe.PropagateTo(xp,lMagneticField);\r
844         \r
845         AliESDv0 vertex(trackInPion,PionIdx,trackInHe,HeIdx);\r
846         if (vertex.GetChi2V0() > fgChi2max) continue;\r
847         \r
848         Float_t CosPointingAngle=vertex.GetV0CosineOfPointingAngle(xPrimaryVertex,yPrimaryVertex,zPrimaryVertex); //PointingAngle\r
849         if (CosPointingAngle < fgCPAmin) continue;\r
850         \r
851         vertex.SetDcaV0Daughters(dca);\r
852         vertex.SetV0CosineOfPointingAngle(CosPointingAngle);\r
853         \r
854         fPos[0]=vertex.Xv();\r
855         fPos[1]=vertex.Yv(); \r
856         fPos[2]=vertex.Zv(); \r
857         \r
858         HeTrack->PxPyPz(momHeVett);\r
859         PionTrack->PxPyPz(momPionVett); \r
860         \r
861         Double_t raggio=TMath::Sqrt(fPos[0]*fPos[0]+fPos[1]*fPos[1]+fPos[2]*fPos[2]);\r
862         HeTrack->GetPxPyPzAt(raggio,lMagneticField,momHeVettAt);\r
863         PionTrack->GetPxPyPzAt(raggio,lMagneticField,momPionVettAt); \r
864         \r
865         //------------------------------------------------------------------------//\r
866         \r
867         HeTrack->GetImpactParameters(impactXY, impactZ);\r
868         \r
869         PionTrack->GetImpactParameters(impactXYpi, impactZpi);\r
870         \r
871         if(vertex.GetD(xPrimaryVertex,yPrimaryVertex,zPrimaryVertex)>3) continue;\r
872         \r
873         //salvo solo fino a 3.1 GeV/c2\r
874         \r
875         vHelium.SetXYZM(2*momHeVettAt[0],2*momHeVettAt[1],2*momHeVettAt[2],Helium3Mass); \r
876         vPion.SetXYZM(momPionVettAt[0],momPionVettAt[1],momPionVettAt[2],PionMass);       \r
877         vSum=vHelium+vPion;\r
878         \r
879         if(vSum.M()>3.1)continue;\r
880         \r
881         //----------------------------------------------------------------------//\r
882         \r
883         vett1[0]=(Float_t)runNumber;\r
884         vett1[1]=(Float_t)BCNumber;\r
885         vett1[2]=(Float_t)OrbitNumber;\r
886         vett1[3]=(Float_t)PeriodNumber;\r
887         vett1[4]=(Float_t)eventtype;\r
888         vett1[5]=(Float_t)TrackNumber;\r
889         vett1[6]=(Float_t)percentile;\r
890         vett1[7]=(Float_t)xPrimaryVertex; //PRIMARY\r
891         vett1[8]=(Float_t)yPrimaryVertex;\r
892         vett1[9]=(Float_t)zPrimaryVertex;\r
893         vett1[10]=(Float_t)fPos[0]; //SECONDARY\r
894         vett1[11]=(Float_t)fPos[1];\r
895         vett1[12]=(Float_t)fPos[2];\r
896         vett1[13]=(Float_t)dca;           //between 2 tracks\r
897         vett1[14]=(Float_t)CosPointingAngle;          //cosPointingAngle da V0\r
898         vett1[15]=(Float_t)vertex.GetD(xPrimaryVertex,yPrimaryVertex,zPrimaryVertex);\r
899         vett1[16]=(Float_t)HeTrack->GetSign(); //He\r
900         vett1[17]=(Float_t)trackInHe.GetP();\r
901         vett1[18]=(Float_t)HeTrack->GetTPCsignal();\r
902         vett1[19]=(Float_t)DcaHeToPrimVertex;\r
903         vett1[20]=(Float_t)HeTrack->Eta();\r
904         vett1[21]=(Float_t)momHeVett[0];\r
905         vett1[22]=(Float_t)momHeVett[1];\r
906         vett1[23]=(Float_t)momHeVett[2];\r
907         vett1[24]=(Float_t)momHeVettAt[0];\r
908         vett1[25]=(Float_t)momHeVettAt[1];\r
909         vett1[26]=(Float_t)momHeVettAt[2];\r
910         vett1[27]=(Float_t)HeTrack->GetTPCNcls();\r
911         vett1[28]=(Float_t)impactXY;\r
912         vett1[29]=(Float_t)impactZ;\r
913         vett1[30]=(Float_t)HeTrack->GetITSClusterMap();\r
914         vett1[31]=(Float_t)IsHeITSRefit;\r
915         vett1[32]=(Float_t)PionTrack->GetSign(); //Pion\r
916         vett1[33]=(Float_t)trackInPion.GetP();\r
917         vett1[34]=(Float_t)PionTrack->GetTPCsignal();\r
918         vett1[35]=(Float_t)DcaPionToPrimVertex;\r
919         vett1[36]=(Float_t)PionTrack->Eta();\r
920         vett1[37]=(Float_t)momPionVett[0];\r
921         vett1[38]=(Float_t)momPionVett[1];\r
922         vett1[39]=(Float_t)momPionVett[2];\r
923         vett1[40]=(Float_t)momPionVettAt[0];\r
924         vett1[41]=(Float_t)momPionVettAt[1];\r
925         vett1[42]=(Float_t)momPionVettAt[2];\r
926         vett1[43]=(Float_t)PionTrack->GetTPCNcls();\r
927         vett1[44]=(Float_t)impactXYpi;\r
928         vett1[45]=(Float_t)impactZpi;\r
929         vett1[46]=(Float_t)PionTrack->GetITSClusterMap();\r
930         vett1[47]=(Float_t)IsPiITSRefit;\r
931         vett1[48]=(Float_t)xn;\r
932         vett1[49]=(Float_t)xp;\r
933         \r
934         fNtuple1->Fill(vett1);  \r
935         \r
936       }// positive TPC\r
937       \r
938     } //negative tpc\r
939     \r
940   }\r
941   \r
942   PostData(1,fListHistCascade);\r
943   \r
944 } //end userexec\r
945 \r
946 \r
947 //________________________________________________________________________\r
948 \r
949 void AliAnalysisTaskHelium3Pi::Terminate(Option_t *) \r
950 {\r
951   // Draw result to the screen\r
952   // Called once at the end of the query\r
953 }\r
954 \r
955 \r