]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FEMTOSCOPY/AliFemto/AliFemtoCutMonitorParticleYPt.cxx
3f200978be0c306d5a38d0d5f4f6f494caa7a913
[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   fMass(0.13957)
23 {
24   // Default constructor
25   fYPt = new TH2D("YPt", "Rapidity vs Pt",              100, -1.0, 1.0, 100, 0.1, 2.0);
26   fYPhi = new TH2D("YPhi", "Rapidity vs Phi",           100, -1.0, 1.0, 100, -TMath::Pi(), TMath::Pi());
27   fPtPhi = new TH2D("PtPhi", "Pt vs Phi",               100,  0.1, 2.0, 100, -TMath::Pi(), TMath::Pi());
28   fEtaPhi = new TH2D("EtaPhi", "Pseudorapidity vs Phi", 100, -1.0, 1.0, 100, -TMath::Pi(), TMath::Pi());
29   fEtaPt = new TH2D("EtaPt", "Pseudorapidity vs Pt",    100, -1.0, 1.0, 100, 0.1, 2.0);
30   fEtaPhiW = new TH2D("EtaPhiW", "Pseudorapidity vs Phi chi2/N weighted", 100, -1.0, 1.0, 100, -TMath::Pi(), TMath::Pi());
31   fEtaPtW = new TH2D("EtaPtW", "Pseudorapidity vs Pt chi2/N weighted",    100, -1.0, 1.0, 100, 0.1, 2.0);
32 }
33
34 AliFemtoCutMonitorParticleYPt::AliFemtoCutMonitorParticleYPt(const char *aName, float aMass):
35   AliFemtoCutMonitor(),
36   fYPt(0),
37   fYPhi(0),
38   fPtPhi(0),
39   fEtaPhi(0),
40   fEtaPt(0),
41   fEtaPhiW(0),
42   fEtaPtW(0),
43   fMass(aMass)
44 {
45   // Normal constructor
46   char name[200];
47   snprintf(name, 200, "YPt%s", aName);
48   fYPt = new TH2D(name, "Rapdity vs Pt", 100, -1.0, 1.0, 100, 0.1, 2.0);
49   snprintf(name, 200, "YPhi%s", aName);
50   fYPhi = new TH2D(name, "Rapidity vs Phi",           100, -1.0, 1.0, 100, -TMath::Pi(), TMath::Pi());
51   snprintf(name, 200, "PtPhi%s", aName);
52   fPtPhi = new TH2D(name, "Pt vs Phi",               100,  0.1, 2.0, 100, -TMath::Pi(), TMath::Pi());
53   snprintf(name, 200, "EtaPhi%s", aName);
54   fEtaPhi = new TH2D(name, "Pseudorapidity vs Phi", 100, -1.0, 1.0, 100, -TMath::Pi(), TMath::Pi());
55   snprintf(name, 200, "EtaPt%s", aName);
56   fEtaPt = new TH2D(name, "Pseudorapidity vs Pt",    100, -1.0, 1.0, 100, 0.1, 2.0);
57   snprintf(name, 200, "EtaPhiW%s", aName);
58   fEtaPhiW = new TH2D(name, "Pseudorapidity vs Phi chi2/N weighted", 100, -1.0, 1.0, 100, -TMath::Pi(), TMath::Pi());
59   snprintf(name, 200, "EtaPtW%s", aName);
60   fEtaPtW = new TH2D(name, "Pseudorapidity vs Pt chi2/N weighted",    100, -1.0, 1.0, 100, 0.1, 2.0);
61 }
62
63 AliFemtoCutMonitorParticleYPt::AliFemtoCutMonitorParticleYPt(const AliFemtoCutMonitorParticleYPt &aCut):
64   AliFemtoCutMonitor(),
65   fYPt(0),
66   fYPhi(0),
67   fPtPhi(0),
68   fEtaPhi(0),
69   fEtaPt(0),
70   fEtaPhiW(0),
71   fEtaPtW(0),
72   fMass(0.13957)
73 {
74   // copy constructor
75   if (fYPt) delete fYPt;
76   fYPt = new TH2D(*aCut.fYPt);
77   fYPhi = new TH2D(*aCut.fYPhi);
78   fPtPhi = new TH2D(*aCut.fPtPhi);
79   fEtaPhi = new TH2D(*aCut.fEtaPhi);
80   fEtaPt = new TH2D(*aCut.fEtaPt);
81   fEtaPhiW = new TH2D(*aCut.fEtaPhiW);
82   fEtaPtW = new TH2D(*aCut.fEtaPtW);
83   fMass = aCut.fMass; 
84 }
85
86 AliFemtoCutMonitorParticleYPt::~AliFemtoCutMonitorParticleYPt()
87 {
88   // Destructor
89   delete fYPt;
90   delete fYPhi;
91   delete fPtPhi;
92   delete fEtaPhi;
93   delete fEtaPt;
94   delete fEtaPhiW;
95   delete fEtaPtW;
96 }
97
98 AliFemtoCutMonitorParticleYPt& AliFemtoCutMonitorParticleYPt::operator=(const AliFemtoCutMonitorParticleYPt& aCut)
99 {
100   // assignment operator
101   if (this == &aCut) 
102     return *this;
103
104   if (fYPt) delete fYPt;
105   fYPt = new TH2D(*aCut.fYPt);
106   if (fYPhi) delete fYPhi;
107   fYPhi = new TH2D(*aCut.fYPhi);
108   if (fPtPhi) delete fPtPhi;
109   fPtPhi = new TH2D(*aCut.fPtPhi);
110   if (fEtaPhi) delete fEtaPhi;
111   fEtaPhi = new TH2D(*aCut.fEtaPhi);
112   if (fEtaPt) delete fEtaPt;
113   fEtaPt = new TH2D(*aCut.fEtaPt);
114   if (fEtaPhiW) delete fEtaPhiW;
115   fEtaPhiW = new TH2D(*aCut.fEtaPhiW);
116   if (fEtaPtW) delete fEtaPtW;
117   fEtaPtW = new TH2D(*aCut.fEtaPtW);
118   
119   return *this;
120 }
121
122 AliFemtoString AliFemtoCutMonitorParticleYPt::Report(){ 
123   // Prepare report from the execution
124   string stemp = "*** AliFemtoCutMonitorParticleYPt report"; 
125   AliFemtoString returnThis = stemp;
126   return returnThis; 
127 }
128
129 void AliFemtoCutMonitorParticleYPt::Fill(const AliFemtoTrack* aTrack)
130 {
131   // Fill in the monitor histograms with the values from the current track
132   float tEnergy = ::sqrt(aTrack->P().mag2()+fMass*fMass);
133   float tRapidity = 0.5*::log((tEnergy+aTrack->P().z())/(tEnergy-aTrack->P().z()));
134   float tPt = ::sqrt((aTrack->P().x())*(aTrack->P().x())+(aTrack->P().y())*(aTrack->P().y()));
135   float tEta = -TMath::Log(TMath::Tan(aTrack->P().theta()/2.0));
136   float tPhi = aTrack->P().phi();
137   float chi2w;
138   if (aTrack->TPCncls() > 0)
139     chi2w = aTrack->TPCchi2()/aTrack->TPCncls();
140   else
141     chi2w = 6.0;
142
143   fYPt->Fill(tRapidity, tPt);
144   fYPhi->Fill(tRapidity, tPhi);
145   fPtPhi->Fill(tPt, tPhi);
146   fEtaPhi->Fill(tEta, tPhi);
147   fEtaPt->Fill(tEta, tPt);
148   fEtaPhiW->Fill(tEta, tPhi, chi2w);
149   fEtaPtW->Fill(tEta, tPt, chi2w);
150 }
151
152 void AliFemtoCutMonitorParticleYPt::Write()
153 {
154   // Write out the relevant histograms
155   fYPt->Write();
156   fYPhi->Write();
157   fPtPhi->Write();
158   fEtaPhi->Write();
159   fEtaPt->Write();
160   fEtaPhiW->Write();
161   fEtaPtW->Write();
162 }
163
164 TList *AliFemtoCutMonitorParticleYPt::GetOutputList()
165 {
166   TList *tOutputList = new TList();
167   tOutputList->Add(fYPt);
168   tOutputList->Add(fYPhi);
169   tOutputList->Add(fPtPhi);
170   tOutputList->Add(fEtaPhi);
171   tOutputList->Add(fEtaPt);
172   tOutputList->Add(fEtaPhiW);
173   tOutputList->Add(fEtaPtW);
174
175   return tOutputList;
176 }