]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FEMTOSCOPY/AliFemto/AliFemtoCutMonitorV0.cxx
Adding V0 femtoscopy analysis
[u/mrichter/AliRoot.git] / PWGCF / FEMTOSCOPY / AliFemto / AliFemtoCutMonitorV0.cxx
1 ////////////////////////////////////////////////////////////////////////////////
2 //                                                                            //
3 // AliFemtoCutMonitorV0                                                       //
4 //                                                                            //
5 ////////////////////////////////////////////////////////////////////////////////
6 #include "AliFemtoCutMonitorV0.h"
7 #include <TH1D.h>
8 #include <TH2D.h>
9 #include <TH3D.h>
10 #include <TList.h>
11 #include "AliFemtoModelHiddenInfo.h"
12
13 AliFemtoCutMonitorV0::AliFemtoCutMonitorV0():
14   fLambdaMass(0),
15   fAntiLambdaMass(0),
16   fK0ShortMass(0),
17   fDcaDaughters(0), 
18   fDcaV0ToPrimVertex(0),
19   fCosPointingAngle(0),
20   fEtaV0(0),
21   fPtV0(0),
22   fPtPosDaughter(0),
23   fPtNegDaughter(0)
24 {
25   // Default constructor
26   fLambdaMass = new TH1F("LambdaMass", "Mass Assuming Lambda Hypothesis", 10000, 0, 5);
27   fAntiLambdaMass = new TH1F("AntiLambdaMass", "Mass Assuming AntiLambda Hypothesis", 10000, 0, 5);
28   fK0ShortMass= new TH1F("K0ShortMass", "Mass Assuming K0 short Hypothesis", 500, 0, 5);
29   fDcaDaughters = new TH1F("DcaDaughters", "DCA Daughters", 500, 0, 2);
30   fDcaV0ToPrimVertex = new TH1F("DcaV0ToPrimVertex", "DCA V0 to primary vertex", 500, 0, 3);
31   fCosPointingAngle = new TH1F("CosPointingAngle","Cosinus Pointing Angle",500,0,1);
32   fEtaV0 = new TH1F("EtaV0", "|Eta| distribution of V0s", 500, 0.0, 8.);
33   fPtV0 = new TH1F("PtV0", "Pt distribution of V0s", 500, 0.0, 8.);
34   fPtPosDaughter = new TH1F("PtPosDaughter", "Pt distribution of positive daughters", 500, 0.0, 5.);
35   fPtNegDaughter = new TH1F("PtNegDaughter", "Pt distribution of negative daughters", 500, 0.0, 5.);
36
37
38   fLambdaMass->Sumw2();
39   fAntiLambdaMass->Sumw2();
40   fK0ShortMass->Sumw2();
41   fDcaDaughters->Sumw2();
42   fDcaV0ToPrimVertex->Sumw2();
43   fCosPointingAngle->Sumw2();
44   fEtaV0->Sumw2();
45   fPtPosDaughter->Sumw2();
46   fPtNegDaughter->Sumw2();
47
48 }
49
50 AliFemtoCutMonitorV0::AliFemtoCutMonitorV0(const char *aName):
51   AliFemtoCutMonitor(),
52   fLambdaMass(0),
53   fAntiLambdaMass(0),
54   fK0ShortMass(0),
55   fDcaDaughters(0), 
56   fDcaV0ToPrimVertex(0),
57   fCosPointingAngle(0),
58   fEtaV0(0),
59   fPtV0(0),
60   fPtPosDaughter(0),
61   fPtNegDaughter(0)
62 {
63   // Normal constructor
64   char name[200];
65   snprintf(name, 200, "LambdaMass%s", aName);
66   fLambdaMass = new TH1F(name, "Mass Assuming Lambda Hypothesis", 10000, 0, 5);
67   snprintf(name, 200, "AntiLambdaMass%s", aName);
68   fAntiLambdaMass = new TH1F(name, "Mass Assuming AntiLambda Hypothesis", 10000, 0, 5);
69   snprintf(name, 200, "K0ShortMass%s", aName);
70   fK0ShortMass = new TH1F(name, "Mass Assuming K0 short Hypothesis", 500, 0, 5);
71   snprintf(name, 200, "DcaDaughters%s", aName);
72   fDcaDaughters = new TH1F(name, "DCA Daughters", 500, 0, 2);
73   snprintf(name, 200, "DcaV0ToPrimVertex%s", aName);
74   fDcaV0ToPrimVertex = new TH1F(name, "DCA V0 to primary vertex", 500, 0, 3);
75   snprintf(name, 200, "CosPointingAngle%s", aName);
76   fCosPointingAngle = new TH1F(name,"Cosinus Pointing Angle",500,0,1);
77   snprintf(name, 200, "EtaV0%s", aName);
78   fEtaV0 = new TH1F(name, "|Eta| distribution of V0s", 500, 0.0, 1.); 
79   snprintf(name, 200, "PtV0%s", aName);
80   fPtV0 = new TH1F(name, "Pt distribution of V0s", 500, 0.0, 8.); 
81   snprintf(name, 200, "fPtPosDaughter%s", aName);
82   fPtPosDaughter = new TH1F(name, "Pt distribution of positive daughters", 500, 0.0, 5.);
83   snprintf(name, 200, "fPtNegDaughter%s", aName);
84   fPtNegDaughter = new TH1F(name, "Pt distribution of negative daughters", 500, 0.0, 5.);
85
86   fLambdaMass->Sumw2();
87   fAntiLambdaMass->Sumw2();
88   fK0ShortMass->Sumw2();
89   fDcaDaughters->Sumw2();
90   fDcaV0ToPrimVertex->Sumw2();
91   fCosPointingAngle->Sumw2();
92   fEtaV0->Sumw2();
93   fPtPosDaughter->Sumw2();
94   fPtNegDaughter->Sumw2();
95 }
96
97 AliFemtoCutMonitorV0::AliFemtoCutMonitorV0(const AliFemtoCutMonitorV0 &aCut):
98   AliFemtoCutMonitor(),
99   fLambdaMass(0),
100   fAntiLambdaMass(0),
101   fK0ShortMass(0),
102   fDcaDaughters(0), 
103   fDcaV0ToPrimVertex(0),
104   fCosPointingAngle(0),
105   fEtaV0(0),
106   fPtV0(0),
107   fPtPosDaughter(0),
108   fPtNegDaughter(0)
109 {
110   // copy constructor
111   if (fLambdaMass) delete fLambdaMass;
112   fLambdaMass = new TH1F(*aCut.fLambdaMass);
113   if (fAntiLambdaMass) delete fAntiLambdaMass;
114   fAntiLambdaMass = new TH1F(*aCut.fAntiLambdaMass);
115   if (fK0ShortMass) delete fK0ShortMass;
116   fK0ShortMass = new TH1F(*aCut.fK0ShortMass);
117   if (fDcaDaughters) delete fDcaDaughters;
118   fDcaDaughters = new TH1F(*aCut.fDcaDaughters);
119   if (fDcaV0ToPrimVertex) delete fDcaV0ToPrimVertex;
120   fDcaV0ToPrimVertex = new TH1F(*aCut.fDcaV0ToPrimVertex);
121   if(fCosPointingAngle) delete fCosPointingAngle;
122   fCosPointingAngle = new TH1F(*aCut.fCosPointingAngle);
123   if(fEtaV0) delete fEtaV0;
124   fEtaV0 = new TH1F(*aCut.fEtaV0);
125   if(fPtV0) delete fPtV0;
126   fPtV0 = new TH1F(*aCut.fPtV0);
127   if(fPtPosDaughter) delete fPtPosDaughter;
128   fPtPosDaughter = new TH1F(*aCut.fPtPosDaughter);
129   if(fPtNegDaughter) delete fPtNegDaughter;
130   fPtNegDaughter = new TH1F(*aCut.fPtNegDaughter);
131
132   fLambdaMass->Sumw2();
133   fAntiLambdaMass->Sumw2();
134   fK0ShortMass->Sumw2();
135   fDcaDaughters->Sumw2();
136   fDcaV0ToPrimVertex->Sumw2();
137   fCosPointingAngle->Sumw2();
138   fEtaV0->Sumw2();
139   fPtPosDaughter->Sumw2();
140   fPtNegDaughter->Sumw2();
141 }
142
143 AliFemtoCutMonitorV0::~AliFemtoCutMonitorV0()
144 {
145   // Destructor
146   delete fLambdaMass;
147   delete fAntiLambdaMass;
148   delete fK0ShortMass;
149   delete fDcaDaughters; 
150   delete fDcaV0ToPrimVertex;
151   delete fCosPointingAngle;
152   delete fEtaV0;
153   delete fPtV0;
154   delete fPtPosDaughter;
155   delete fPtNegDaughter;
156 }
157
158 AliFemtoCutMonitorV0& AliFemtoCutMonitorV0::operator=(const AliFemtoCutMonitorV0& aCut)
159 {
160   // assignment operator
161   if (this == &aCut) 
162     return *this;
163
164   if (fLambdaMass) delete fLambdaMass;
165   fLambdaMass = new TH1F(*aCut.fLambdaMass);
166   if (fAntiLambdaMass) delete fAntiLambdaMass;
167   fAntiLambdaMass = new TH1F(*aCut.fAntiLambdaMass);
168   if (fK0ShortMass) delete fK0ShortMass;
169   fK0ShortMass = new TH1F(*aCut.fK0ShortMass);
170   if (fDcaDaughters) delete fDcaDaughters;
171   fDcaDaughters = new TH1F(*aCut.fDcaDaughters);
172   if (fDcaV0ToPrimVertex) delete fDcaV0ToPrimVertex;
173   fDcaV0ToPrimVertex = new TH1F(*aCut.fDcaV0ToPrimVertex);
174   if(fCosPointingAngle) delete fCosPointingAngle;
175   fCosPointingAngle = new TH1F(*aCut.fCosPointingAngle);
176   if(fEtaV0) delete fEtaV0;
177   fEtaV0 = new TH1F(*aCut.fEtaV0);
178   if(fPtV0) delete fPtV0;
179   fPtV0 = new TH1F(*aCut.fPtV0);
180   if(fPtPosDaughter) delete fPtPosDaughter;
181   fPtPosDaughter = new TH1F(*aCut.fPtPosDaughter);
182   if(fPtNegDaughter) delete fPtNegDaughter;
183   fPtNegDaughter = new TH1F(*aCut.fPtNegDaughter);
184
185   fLambdaMass->Sumw2();
186   fAntiLambdaMass->Sumw2();
187   fK0ShortMass->Sumw2();
188   fDcaDaughters->Sumw2();
189   fDcaV0ToPrimVertex->Sumw2();
190   fCosPointingAngle->Sumw2();
191   fEtaV0->Sumw2();
192   fPtPosDaughter->Sumw2();
193   fPtNegDaughter->Sumw2();
194
195   return *this;
196 }
197
198 AliFemtoString AliFemtoCutMonitorV0::Report(){ 
199   // Prepare report from the execution
200   string stemp = "*** AliFemtoCutMonitorV0 report"; 
201   AliFemtoString returnThis = stemp;
202   return returnThis; 
203 }
204
205 void AliFemtoCutMonitorV0::Fill(const AliFemtoV0* aV0)
206 {
207   // Fill momentum resolution histograms for the particle
208   fLambdaMass->Fill(aV0->MassLambda());
209   fAntiLambdaMass->Fill(aV0->MassAntiLambda());
210   fK0ShortMass->Fill(aV0->MassK0Short());
211   fDcaDaughters->Fill(aV0->DcaV0Daughters()); 
212   fDcaV0ToPrimVertex->Fill(aV0->DcaV0ToPrimVertex());
213   fCosPointingAngle->Fill(aV0->CosPointingAngle());
214   fEtaV0->Fill(aV0->EtaV0());
215   fPtV0->Fill(aV0->PtV0());
216   fPtPosDaughter->Fill(aV0->PtPos());
217   fPtNegDaughter->Fill(aV0->PtNeg());
218 }
219
220 void AliFemtoCutMonitorV0::Write()
221 {
222   // Write out the relevant histograms
223   fLambdaMass->Write();
224   fAntiLambdaMass->Write();
225   fK0ShortMass->Write();
226   fDcaDaughters->Write();
227   fDcaV0ToPrimVertex->Write();
228   fCosPointingAngle->Write();
229   fEtaV0->Write();
230   fPtV0->Write();
231   fPtPosDaughter->Write();
232   fPtNegDaughter->Write();
233 }
234
235 TList *AliFemtoCutMonitorV0::GetOutputList()
236 {
237   // Get the list of histograms to write
238   TList *tOutputList = new TList();
239   tOutputList->Add(fLambdaMass);
240   tOutputList->Add(fAntiLambdaMass);
241   tOutputList->Add(fK0ShortMass);
242   tOutputList->Add(fDcaDaughters);
243   tOutputList->Add(fDcaV0ToPrimVertex);
244   tOutputList->Add(fCosPointingAngle);
245   tOutputList->Add(fEtaV0);
246   tOutputList->Add(fPtV0);
247   tOutputList->Add(fPtPosDaughter);
248   tOutputList->Add(fPtNegDaughter);
249
250   return tOutputList;
251 }