]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FEMTOSCOPY/AliFemtoUser/AliFemtoCorrFctnMinvMonitor.cxx
Split: fixed incpaths for ANALYSISalice -> OADB
[u/mrichter/AliRoot.git] / PWGCF / FEMTOSCOPY / AliFemtoUser / AliFemtoCorrFctnMinvMonitor.cxx
1 ////////////////////////////////////////////////////////////////////////////////
2 //                                                                            //
3 // AliFemtoCorrFctnGammaMonitor - A correlation function that analyzes        //
4 // two particle mass minvariant with various mass assumptions                 //
5 //                                                                            //
6 // Authors: MaƂgorzata Janik majanik@cern.ch
7 //          Anna Zaborowska azaborow@cern.ch                                   //
8 //                                                                            //
9 ////////////////////////////////////////////////////////////////////////////////
10
11 #include "AliFemtoCorrFctnMinvMonitor.h"
12 #include <cstdio>
13 #include <TMath.h>
14
15 #ifdef __ROOT__
16 ClassImp(AliFemtoCorrFctnMinvMonitor)
17 #endif
18
19 //____________________________
20 AliFemtoCorrFctnMinvMonitor::AliFemtoCorrFctnMinvMonitor(char* title):
21    AliFemtoCorrFctn(),
22    fMinveeFail(0),
23    fMinvee(0),
24    fMinv2piFail(0),
25    fMinv2pi(0),
26    fMinvppiFail(0),
27    fMinvppi(0)
28 {
29   fMinveeFail = new TH1D(Form("MinveeGamma%s",title), "ee mass assumption GAMMA, minv",1000, 0.0, 10.0);
30    fMinvee = new TH1D(Form("Minvee%s",title), "ee mass assumption, minv",1000, 0.0, 10.0);
31    fMinv2piFail = new TH1D(Form("Minv2piResonances%s",title), "pipi mass assumption RESONANCES, minv",1000, 0.0, 10.0);
32    fMinv2pi = new TH1D(Form("Minv2pi%s",title), "pipi mass assumption, minv",1000, 0.0, 10.0);
33    fMinvppiFail = new TH1D(Form("MinvppiResonances%s",title), "ppi mass assumption RESONANCES, minv",1000, 0.0, 10.0);
34    fMinvppi = new TH1D(Form("Minvppi%s",title), "ppi mass assumption, minv",1000, 0.0, 10.0);
35 }
36
37 //____________________________
38 AliFemtoCorrFctnMinvMonitor::AliFemtoCorrFctnMinvMonitor(const AliFemtoCorrFctnMinvMonitor& aCorrFctn) :
39   AliFemtoCorrFctn(),
40    fMinveeFail(0),
41    fMinvee(0),
42    fMinv2piFail(0),
43    fMinv2pi(0),
44    fMinvppiFail(0),
45    fMinvppi(0)
46 {
47   // copy constructor
48   if (fMinveeFail) delete fMinveeFail;
49   fMinveeFail = new TH1D(*aCorrFctn.fMinveeFail);
50   if (fMinvee) delete fMinvee;
51   fMinvee = new TH1D(*aCorrFctn.fMinvee);
52   if (fMinv2piFail) delete fMinv2piFail;
53   fMinv2piFail = new TH1D(*aCorrFctn.fMinv2piFail);
54   if (fMinv2pi) delete fMinv2pi;
55   fMinv2pi = new TH1D(*aCorrFctn.fMinv2pi);
56   if (fMinvppiFail) delete fMinvppiFail;
57   fMinvppiFail = new TH1D(*aCorrFctn.fMinvppiFail);
58   if (fMinvppi) delete fMinvppi;
59   fMinvppi = new TH1D(*aCorrFctn.fMinvppi);
60 }
61 //____________________________
62 AliFemtoCorrFctnMinvMonitor::~AliFemtoCorrFctnMinvMonitor(){
63   // destructor
64     delete fMinveeFail;
65   delete fMinvee;
66   delete fMinv2piFail;
67   delete fMinv2pi;
68   delete fMinvppiFail;
69   delete fMinvppi;
70 }
71 //_________________________
72 AliFemtoCorrFctnMinvMonitor& AliFemtoCorrFctnMinvMonitor::operator=(const AliFemtoCorrFctnMinvMonitor& aCorrFctn)
73 {
74   // assignment operator
75   if (this == &aCorrFctn)
76     return *this;
77
78   if (fMinveeFail) delete fMinveeFail;
79   fMinveeFail = new TH1D(*aCorrFctn.fMinveeFail);
80   if (fMinvee) delete fMinvee;
81   fMinvee = new TH1D(*aCorrFctn.fMinvee);
82   if (fMinv2piFail) delete fMinv2piFail;
83   fMinv2piFail = new TH1D(*aCorrFctn.fMinv2piFail);
84   if (fMinv2pi) delete fMinv2pi;
85   fMinv2pi = new TH1D(*aCorrFctn.fMinv2pi);
86   if (fMinvppiFail) delete fMinvppiFail;
87   fMinvppiFail = new TH1D(*aCorrFctn.fMinvppiFail);
88   if (fMinvppi) delete fMinvppi;
89   fMinvppi = new TH1D(*aCorrFctn.fMinvppi);
90
91   return *this;
92 }
93 //_________________________
94 void AliFemtoCorrFctnMinvMonitor::Finish(){
95   // here is where we should normalize, fit, etc...
96   // we should NOT Draw() the histos (as I had done it below),
97   // since we want to insulate ourselves from root at this level
98   // of the code.  Do it instead at root command line with browser.
99   //  mShareNumerator->Draw();
100   //mShareDenominator->Draw();
101   //mRatio->Draw();
102
103 }
104
105 //____________________________
106 AliFemtoString AliFemtoCorrFctnMinvMonitor::Report(){
107   // create report
108   string stemp = "Mass invariant Monitor Function Report\n";
109   AliFemtoString returnThis = stemp;
110   return returnThis;
111 }
112 //____________________________
113 void AliFemtoCorrFctnMinvMonitor::AddRealPair( AliFemtoPair* pair){
114    double me = 0.000511;
115   double mPi = 0.13957018;
116   double mp = 0.938272046;
117
118   double mgammamax = 0.04;
119   double mK0min = 0.00049;
120   double mK0max = 0.00051;
121   //double mK0 = 0.000497614;
122   double mRhomin = 0.000765;
123   double mRhomax = 0.000785;
124   //double mRho = 0.00077526;
125   double mLmin = 1.095;
126   double mLmax = 1.135;
127   //double mL = 1.115683;
128
129   if ((pair->Track1()->Track()->Charge() * pair->Track2()->Track()->Charge()) < 0.0) {
130
131     // check on ee pairs (gamma)
132     double e1 = TMath::Sqrt(me*me + pair->Track1()->Track()->P().Mag2());
133     double e2 = TMath::Sqrt(me*me + pair->Track2()->Track()->P().Mag2());
134     double minvGamma = 2*me*me + 2*(e1*e2 -
135                                pair->Track1()->Track()->P().x()*pair->Track2()->Track()->P().x() -
136                                pair->Track1()->Track()->P().y()*pair->Track2()->Track()->P().y() -
137                                pair->Track1()->Track()->P().z()*pair->Track2()->Track()->P().z());
138     if ( minvGamma < mgammamax )
139     {
140        fMinveeFail->Fill(minvGamma);
141     }
142     else fMinvee->Fill(minvGamma);
143     //check on resonances
144     double pi1 =  TMath::Sqrt(mPi*mPi + pair->Track1()->Track()->P().Mag2());
145     double pi2 =  TMath::Sqrt(mPi*mPi + pair->Track2()->Track()->P().Mag2());
146     double p1 =  TMath::Sqrt(mp*mp + pair->Track1()->Track()->P().Mag2());
147     double p2 =  TMath::Sqrt(mp*mp + pair->Track2()->Track()->P().Mag2());
148     //check on K0 and Rho
149     double minv2pi = 2*mPi*mPi + 2*(pi1*pi2 -
150                                pair->Track1()->Track()->P().x()*pair->Track2()->Track()->P().x() -
151                                pair->Track1()->Track()->P().y()*pair->Track2()->Track()->P().y() -
152                                pair->Track1()->Track()->P().z()*pair->Track2()->Track()->P().z());
153     if ( ((minv2pi>mK0min && minv2pi<mK0max) || (minv2pi>mRhomin && minv2pi<mRhomax)) )
154     {
155        fMinv2piFail->Fill(minv2pi);
156     }
157     else fMinv2pi->Fill(minv2pi);
158     //check on L0
159     double minvpPi = 2*mp*mPi + 2*(p1*pi2 -
160                                pair->Track1()->Track()->P().x()*pair->Track2()->Track()->P().x() -
161                                pair->Track1()->Track()->P().y()*pair->Track2()->Track()->P().y() -
162                                pair->Track1()->Track()->P().z()*pair->Track2()->Track()->P().z());
163     double minvPip = 2*mPi*mp + 2*(pi1*p2 -
164                                pair->Track1()->Track()->P().x()*pair->Track2()->Track()->P().x() -
165                                pair->Track1()->Track()->P().y()*pair->Track2()->Track()->P().y() -
166                                pair->Track1()->Track()->P().z()*pair->Track2()->Track()->P().z());
167     if( (minvpPi>mLmin) && (minvpPi<mLmax) )
168     {
169        fMinvppiFail->Fill(minvpPi);
170     }
171     else fMinvppi->Fill(minvpPi);
172     if( (minvPip>mLmin) && (minvPip<mLmax) )
173     {
174        fMinvppiFail->Fill(minvPip);
175     }
176     else fMinvppi->Fill(minvPip);
177   }
178 }
179 //____________________________
180
181 void AliFemtoCorrFctnMinvMonitor::WriteHistos()
182 {
183   // Write out result histograms
184   fMinveeFail->Write();
185   fMinvee->Write();
186   fMinv2piFail->Write();
187   fMinv2pi->Write();
188   fMinvppiFail->Write();
189   fMinvppi->Write();
190 }
191
192 TList* AliFemtoCorrFctnMinvMonitor::GetOutputList()
193 {
194   // Prepare the list of objects to be written to the output
195   TList *tOutputList = new TList();
196
197   tOutputList->Add(fMinveeFail);
198   tOutputList->Add(fMinvee);
199   tOutputList->Add(fMinv2piFail);
200   tOutputList->Add(fMinv2pi);
201   tOutputList->Add(fMinvppiFail);
202   tOutputList->Add(fMinvppi);
203
204   return tOutputList;
205 }