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