]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FEMTOSCOPY/AliFemto/AliFemtoCutMonitorParticleYPt.cxx
Fix Gamma calculation^CliFemtoModelGausLCMSFreezeOutGenerator.cxx
[u/mrichter/AliRoot.git] / PWG2 / FEMTOSCOPY / AliFemto / AliFemtoCutMonitorParticleYPt.cxx
1 ////////////////////////////////////////////////////////////////////////////////
2 //                                                                            //
3 // AliFemtoCutMonitorParticleYPt - the cut monitor for particles to study    //
4 // the difference between reconstructed and true momentum                     //
5 //                                                                            //
6 ////////////////////////////////////////////////////////////////////////////////
7 #include "AliFemtoCutMonitorParticleYPt.h"
8 #include "AliFemtoModelHiddenInfo.h"
9 #include <TH1D.h>
10 #include <TH2D.h>
11 #include <TList.h>
12 #include <TMath.h>
13
14 AliFemtoCutMonitorParticleYPt::AliFemtoCutMonitorParticleYPt():
15   fYPt(0),
16   fYPhi(0),
17   fPtPhi(0),
18   fEtaPhi(0),
19   fEtaPt(0),
20   fEtaPhiW(0),
21   fEtaPtW(0),
22   fDCARPt(0),
23   fDCAZPt(0),
24   fMass(0.13957)
25 {
26   // Default constructor
27   fYPt = new TH2D("YPt", "Rapidity vs Pt",              100, -1.0, 1.0, 100, 0.0, 2.0);
28   fYPhi = new TH2D("YPhi", "Rapidity vs Phi",           100, -1.0, 1.0, 100, -TMath::Pi(), TMath::Pi());
29   fPtPhi = new TH2D("PtPhi", "Pt vs Phi",               100,  0.0, 2.0, 100, -TMath::Pi(), TMath::Pi());
30   fEtaPhi = new TH2D("EtaPhi", "Pseudorapidity vs Phi", 100, -1.0, 1.0, 100, -TMath::Pi(), TMath::Pi());
31   fEtaPt = new TH2D("EtaPt", "Pseudorapidity vs Pt",    100, -1.0, 1.0, 100, 0.0, 2.0);
32   fEtaPhiW = new TH2D("EtaPhiW", "Pseudorapidity vs Phi chi2/N weighted", 100, -1.0, 1.0, 100, -TMath::Pi(), TMath::Pi());
33   fEtaPtW = new TH2D("EtaPtW", "Pseudorapidity vs Pt chi2/N weighted",    100, -1.0, 1.0, 100, 0.0, 2.0);
34   fDCARPt = new TH2D("DCARPt", "DCA in XY vs. Pt", 400, -2.0, 2.0, 100,0.0,2.0);
35   fDCAZPt = new TH2D("DCAZPt", "DCA in Z vs. Pt", 400, -2.0, 2.0, 100,0.0,2.0);
36 }
37
38 AliFemtoCutMonitorParticleYPt::AliFemtoCutMonitorParticleYPt(const char *aName, float aMass):
39   AliFemtoCutMonitor(),
40   fYPt(0),
41   fYPhi(0),
42   fPtPhi(0),
43   fEtaPhi(0),
44   fEtaPt(0),
45   fEtaPhiW(0),
46   fEtaPtW(0),
47   fDCARPt(0),
48   fDCAZPt(0),
49   fMass(aMass)
50 {
51   // Normal constructor
52   char name[200];
53   snprintf(name, 200, "YPt%s", aName);
54   fYPt = new TH2D(name, "Rapdity vs Pt", 100, -1.0, 1.0, 100, 0.0, 2.0);
55   snprintf(name, 200, "YPhi%s", aName);
56   fYPhi = new TH2D(name, "Rapidity vs Phi",           100, -1.0, 1.0, 100, -TMath::Pi(), TMath::Pi());
57   snprintf(name, 200, "PtPhi%s", aName);
58   fPtPhi = new TH2D(name, "Pt vs Phi",               100,  0.0, 2.0, 100, -TMath::Pi(), TMath::Pi());
59   snprintf(name, 200, "EtaPhi%s", aName);
60   fEtaPhi = new TH2D(name, "Pseudorapidity vs Phi", 100, -1.0, 1.0, 100, -TMath::Pi(), TMath::Pi());
61   snprintf(name, 200, "EtaPt%s", aName);
62   fEtaPt = new TH2D(name, "Pseudorapidity vs Pt",    100, -1.0, 1.0, 100, 0.0, 2.0);
63   snprintf(name, 200, "EtaPhiW%s", aName);
64   fEtaPhiW = new TH2D(name, "Pseudorapidity vs Phi chi2/N weighted", 100, -1.0, 1.0, 100, -TMath::Pi(), TMath::Pi());
65   snprintf(name, 200, "EtaPtW%s", aName);
66   fEtaPtW = new TH2D(name, "Pseudorapidity vs Pt chi2/N weighted",    100, -1.0, 1.0, 100, 0.0, 2.0);
67   snprintf(name, 200, "DCARPt%s", aName);
68   fDCARPt = new TH2D(name, "DCA in XY vs. Pt", 400, -2.0, 2.0, 100,0.0,2.0);
69   snprintf(name, 200, "DCAZPt%s", aName);
70   fDCAZPt = new TH2D(name, "DCA in Z vs. Pt", 400, -2.0, 2.0, 100,0.0,2.0);
71 }
72
73 AliFemtoCutMonitorParticleYPt::AliFemtoCutMonitorParticleYPt(const AliFemtoCutMonitorParticleYPt &aCut):
74   AliFemtoCutMonitor(),
75   fYPt(0),
76   fYPhi(0),
77   fPtPhi(0),
78   fEtaPhi(0),
79   fEtaPt(0),
80   fEtaPhiW(0),
81   fEtaPtW(0),
82   fDCARPt(0),
83   fDCAZPt(0),
84   fMass(0.13957)
85 {
86   // copy constructor
87   if (fYPt) delete fYPt;
88   fYPt = new TH2D(*aCut.fYPt);
89   fYPhi = new TH2D(*aCut.fYPhi);
90   fPtPhi = new TH2D(*aCut.fPtPhi);
91   fEtaPhi = new TH2D(*aCut.fEtaPhi);
92   fEtaPt = new TH2D(*aCut.fEtaPt);
93   fEtaPhiW = new TH2D(*aCut.fEtaPhiW);
94   fEtaPtW = new TH2D(*aCut.fEtaPtW);
95   fDCARPt = new TH2D(*aCut.fDCARPt);
96   fDCAZPt = new TH2D(*aCut.fDCAZPt);
97   fMass = aCut.fMass; 
98 }
99
100 AliFemtoCutMonitorParticleYPt::~AliFemtoCutMonitorParticleYPt()
101 {
102   // Destructor
103   delete fYPt;
104   delete fYPhi;
105   delete fPtPhi;
106   delete fEtaPhi;
107   delete fEtaPt;
108   delete fEtaPhiW;
109   delete fEtaPtW;
110   delete fDCARPt;
111   delete fDCAZPt;
112 }
113
114 AliFemtoCutMonitorParticleYPt& AliFemtoCutMonitorParticleYPt::operator=(const AliFemtoCutMonitorParticleYPt& aCut)
115 {
116   // assignment operator
117   if (this == &aCut) 
118     return *this;
119
120   if (fYPt) delete fYPt;
121   fYPt = new TH2D(*aCut.fYPt);
122   if (fYPhi) delete fYPhi;
123   fYPhi = new TH2D(*aCut.fYPhi);
124   if (fPtPhi) delete fPtPhi;
125   fPtPhi = new TH2D(*aCut.fPtPhi);
126   if (fEtaPhi) delete fEtaPhi;
127   fEtaPhi = new TH2D(*aCut.fEtaPhi);
128   if (fEtaPt) delete fEtaPt;
129   fEtaPt = new TH2D(*aCut.fEtaPt);
130   if (fEtaPhiW) delete fEtaPhiW;
131   fEtaPhiW = new TH2D(*aCut.fEtaPhiW);
132   if (fEtaPtW) delete fEtaPtW;
133   fEtaPtW = new TH2D(*aCut.fEtaPtW);
134   if (fDCARPt) delete fDCARPt;
135   fDCARPt = new TH2D(*aCut.fDCARPt);
136   if (fDCAZPt) delete fDCAZPt;
137   fDCAZPt = new TH2D(*aCut.fDCAZPt);
138   
139   return *this;
140 }
141
142 AliFemtoString AliFemtoCutMonitorParticleYPt::Report(){ 
143   // Prepare report from the execution
144   string stemp = "*** AliFemtoCutMonitorParticleYPt report"; 
145   AliFemtoString returnThis = stemp;
146   return returnThis; 
147 }
148
149 void AliFemtoCutMonitorParticleYPt::Fill(const AliFemtoTrack* aTrack)
150 {
151   // Fill in the monitor histograms with the values from the current track
152   float tEnergy = ::sqrt(aTrack->P().mag2()+fMass*fMass);
153   float tRapidity = 0.5*::log((tEnergy+aTrack->P().z())/(tEnergy-aTrack->P().z()));
154   float tPt = ::sqrt((aTrack->P().x())*(aTrack->P().x())+(aTrack->P().y())*(aTrack->P().y()));
155   float tEta = -TMath::Log(TMath::Tan(aTrack->P().theta()/2.0));
156   float tPhi = aTrack->P().phi();
157   float chi2w;
158   float dcar = aTrack->ImpactD();
159   float dcaz = aTrack->ImpactZ();
160   if (aTrack->TPCncls() > 0)
161     chi2w = aTrack->TPCchi2()/aTrack->TPCncls();
162   else
163     chi2w = 6.0;
164
165   fYPt->Fill(tRapidity, tPt);
166   fYPhi->Fill(tRapidity, tPhi);
167   fPtPhi->Fill(tPt, tPhi);
168   fEtaPhi->Fill(tEta, tPhi);
169   fEtaPt->Fill(tEta, tPt);
170   fEtaPhiW->Fill(tEta, tPhi, chi2w);
171   fEtaPtW->Fill(tEta, tPt, chi2w);
172   fDCARPt->Fill(dcar, tPt);
173   fDCAZPt->Fill(dcaz, tPt);
174 }
175
176 void AliFemtoCutMonitorParticleYPt::Write()
177 {
178   // Write out the relevant histograms
179   fYPt->Write();
180   fYPhi->Write();
181   fPtPhi->Write();
182   fEtaPhi->Write();
183   fEtaPt->Write();
184   fEtaPhiW->Write();
185   fEtaPtW->Write();
186   fDCARPt->Write();
187   fDCAZPt->Write();
188 }
189
190 TList *AliFemtoCutMonitorParticleYPt::GetOutputList()
191 {
192   TList *tOutputList = new TList();
193   tOutputList->Add(fYPt);
194   tOutputList->Add(fYPhi);
195   tOutputList->Add(fPtPhi);
196   tOutputList->Add(fEtaPhi);
197   tOutputList->Add(fEtaPt);
198   tOutputList->Add(fEtaPhiW);
199   tOutputList->Add(fEtaPtW);
200   tOutputList->Add(fDCARPt);
201   tOutputList->Add(fDCAZPt);
202
203   return tOutputList;
204 }