]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FEMTOSCOPY/AliFemto/AliFemtoV0PairCut.cxx
updates in macros for Femto QA in train
[u/mrichter/AliRoot.git] / PWGCF / FEMTOSCOPY / AliFemto / AliFemtoV0PairCut.cxx
1 /////////////////////////////////////////////////////////////////////////////
2 //                                                                         //
3 // AliFemtoShareQualityPairCut - a pair cut which checks for some pair     //
4 // qualities that attempt to identify slit/doubly reconstructed tracks     //
5 //                                                                         //
6 /////////////////////////////////////////////////////////////////////////////
7 /***************************************************************************
8  *
9  * $Id: AliFemtoShareQualityPairCut.cxx 50722 2011-07-21 15:18:38Z akisiel $
10  *
11  * Author: Adam Kisiel, Ohio State, kisiel@mps.ohio-state.edu
12  ***************************************************************************
13  *
14  * Description: part of STAR HBT Framework: AliFemtoMaker package
15  *   a cut to remove "shared" and "split" pairs
16  *
17  ***************************************************************************
18  *
19  *
20  **************************************************************************/
21
22 #include "AliFemtoV0PairCut.h"
23 #include <string>
24 #include <cstdio>
25
26 #ifdef __ROOT__
27 ClassImp(AliFemtoV0PairCut)
28 #endif
29
30 //__________________
31 AliFemtoV0PairCut::AliFemtoV0PairCut():
32   fNPairsPassed(0),
33   fNPairsFailed(0),
34   fV0Max(1.0),
35   fShareFractionMax(1.0),
36   fRemoveSameLabel(0),
37   fDataType(kAOD),
38   fDTPCMin(0),
39   fDTPCExitMin(0),
40   fMinAvgSepPosPos(0),
41   fMinAvgSepPosNeg(0),
42   fMinAvgSepNegPos(0),
43   fMinAvgSepNegNeg(0)
44 {
45   // Default constructor
46   // Nothing to do
47 }
48 //__________________
49 AliFemtoV0PairCut::~AliFemtoV0PairCut(){
50   /* no-op */
51 }
52
53 AliFemtoV0PairCut& AliFemtoV0PairCut::operator=(const AliFemtoV0PairCut& cut) 
54 {
55   if (this != &cut) {
56    
57     AliFemtoPairCut::operator=(cut);
58     fNPairsPassed = 0;
59     fNPairsFailed =0;
60     fV0Max = 1.0;
61     fShareFractionMax = 1.0;
62     fRemoveSameLabel = 0;
63     fDataType = kAOD;
64     fDTPCMin = 0;
65     fDTPCExitMin = 0;
66     fMinAvgSepPosPos = 0;
67     fMinAvgSepPosNeg = 0;
68     fMinAvgSepNegPos = 0;
69     fMinAvgSepNegNeg = 0;
70   }
71
72   return *this;
73 }
74
75
76 //__________________
77 bool AliFemtoV0PairCut::Pass(const AliFemtoPair* pair){
78   // Check for pairs that are possibly shared/double reconstruction
79
80   bool temp = true;
81
82   /*cout<<"pair->Track1(): "<<pair->Track1()<<endl;
83   cout<<"pair->Track2(): "<<pair->Track2()<<endl;
84   cout<<"pair->Track1()->V0(): "<<pair->Track1()->V0()<<endl;
85   cout<<"pair->Track2()->V0(): "<<pair->Track2()->V0()<<endl;
86   cout<<"pair->Track1()->V0()->IdNeg(): "<<pair->Track1()->V0()->IdNeg()<<endl;
87   cout<<"pair->Track2()->V0()->IdNeg(): "<<pair->Track2()->V0()->IdNeg()<<endl;
88   cout<<"pair->Track1()->V0()->IdPos(): "<<pair->Track1()->V0()->IdPos()<<endl;
89   cout<<"pair->Track2()->V0()->IdPos(): "<<pair->Track2()->V0()->IdPos()<<endl;*/
90
91   bool tempTPCEntrancePos = true;
92   bool tempTPCEntranceNeg = true;
93   bool tempTPCExitPos = true;
94   bool tempTPCExitNeg = true;
95   if(fDataType==kESD || fDataType==kAOD)
96     {
97       double distx = pair->Track1()->V0()->NominalTpcEntrancePointPos().x() - pair->Track2()->V0()->NominalTpcEntrancePointPos().x();
98       double disty = pair->Track1()->V0()->NominalTpcEntrancePointPos().y() - pair->Track2()->V0()->NominalTpcEntrancePointPos().y();
99       double distz = pair->Track1()->V0()->NominalTpcEntrancePointPos().z() - pair->Track2()->V0()->NominalTpcEntrancePointPos().z();
100       double distPos = sqrt(distx*distx + disty*disty + distz*distz);
101
102       distx = pair->Track1()->V0()->NominalTpcEntrancePointNeg().x() - pair->Track2()->V0()->NominalTpcEntrancePointNeg().x();
103       disty = pair->Track1()->V0()->NominalTpcEntrancePointNeg().y() - pair->Track2()->V0()->NominalTpcEntrancePointNeg().y();
104       distz = pair->Track1()->V0()->NominalTpcEntrancePointNeg().z() - pair->Track2()->V0()->NominalTpcEntrancePointNeg().z();
105       double distNeg = sqrt(distx*distx + disty*disty + distz*distz);
106
107       double distExitX = pair->Track1()->V0()->NominalTpcExitPointPos().x() - pair->Track2()->V0()->NominalTpcExitPointPos().x();
108       double distExitY = pair->Track1()->V0()->NominalTpcExitPointPos().y() - pair->Track2()->V0()->NominalTpcExitPointPos().y();
109       double distExitZ = pair->Track1()->V0()->NominalTpcExitPointPos().z() - pair->Track2()->V0()->NominalTpcExitPointPos().z();
110       double distExitPos = sqrt(distExitX*distExitX + distExitY*distExitY + distExitZ*distExitZ);
111
112       distExitX = pair->Track1()->V0()->NominalTpcExitPointNeg().x() - pair->Track2()->V0()->NominalTpcExitPointNeg().x();
113       distExitY = pair->Track1()->V0()->NominalTpcExitPointNeg().y() - pair->Track2()->V0()->NominalTpcExitPointNeg().y();
114       distExitZ = pair->Track1()->V0()->NominalTpcExitPointNeg().z() - pair->Track2()->V0()->NominalTpcExitPointNeg().z();
115       double distExitNeg = sqrt(distExitX*distExitX + distExitY*distExitY + distExitZ*distExitZ);
116
117       tempTPCEntrancePos = distPos > fDTPCMin;
118       tempTPCEntranceNeg = distNeg > fDTPCMin;
119
120       tempTPCExitPos = distExitPos > fDTPCExitMin;
121       tempTPCExitNeg = distExitNeg > fDTPCExitMin;
122     }
123  
124
125   if(!(pair->Track1()->V0() && pair->Track2()->V0()))
126     {
127       return false;
128     }
129   if(pair->Track1()->V0()->IdNeg()==pair->Track2()->V0()->IdNeg() || pair->Track1()->V0()->IdPos()==pair->Track2()->V0()->IdPos())
130     {
131
132       return false;
133     }
134
135   if(!tempTPCEntrancePos || !tempTPCEntranceNeg || !tempTPCExitPos || !tempTPCExitNeg) return false;  
136
137   double avgSep=0;
138   AliFemtoThreeVector first, second, tmp;
139   for(int i=0; i<8 ;i++)
140     {
141       tmp = pair->Track1()->V0()->NominalTpcPointPos(i);
142       //cout<<"X pos: "<<tmp.x()<<endl;
143       first.SetX((double)(tmp.x()));
144       first.SetY((double)tmp.y());
145       first.SetZ((double)tmp.z());
146       
147       tmp = pair->Track2()->V0()->NominalTpcPointPos(i);
148       second.SetX((double)tmp.x());
149       second.SetY((double)tmp.y());
150       second.SetZ((double)tmp.z()); 
151
152       avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
153     }
154   avgSep /= 8;
155
156   if(avgSep<fMinAvgSepPosPos) return false;
157
158   avgSep = 0;
159       
160
161   for(int i=0; i<8 ;i++)
162     {
163       tmp = pair->Track1()->V0()->NominalTpcPointPos(i);
164       first.SetX((double)(tmp.x()));
165       first.SetY((double)tmp.y());
166       first.SetZ((double)tmp.z());
167       
168       tmp = pair->Track2()->V0()->NominalTpcPointNeg(i);
169       //cout<<"X neg: "<<tmp.x()<<endl;
170       second.SetX((double)tmp.x());
171       second.SetY((double)tmp.y());
172       second.SetZ((double)tmp.z()); 
173
174       avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
175     }
176   avgSep /= 8;
177   if(avgSep<fMinAvgSepPosNeg) return false;
178
179   avgSep = 0;
180
181   for(int i=0; i<8 ;i++)
182     {
183       tmp = pair->Track1()->V0()->NominalTpcPointNeg(i);
184       first.SetX((double)(tmp.x()));
185       first.SetY((double)tmp.y());
186       first.SetZ((double)tmp.z());
187       
188       tmp = pair->Track2()->V0()->NominalTpcPointPos(i);
189       second.SetX((double)tmp.x());
190       second.SetY((double)tmp.y());
191       second.SetZ((double)tmp.z()); 
192
193       avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
194     }
195   avgSep /= 8;
196   if(avgSep<fMinAvgSepNegPos) return false;
197       
198   avgSep = 0;
199  
200   for(int i=0; i<8 ;i++)
201     {
202       tmp = pair->Track1()->V0()->NominalTpcPointNeg(i);
203       first.SetX((double)(tmp.x()));
204       first.SetY((double)tmp.y());
205       first.SetZ((double)tmp.z());
206       
207       tmp = pair->Track2()->V0()->NominalTpcPointNeg(i);
208       second.SetX((double)tmp.x());
209       second.SetY((double)tmp.y());
210       second.SetZ((double)tmp.z()); 
211
212       avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
213     }
214   avgSep /= 8;
215   if(avgSep<fMinAvgSepNegNeg) return false;
216
217   avgSep = 0;
218
219
220   return temp;
221 }
222 //__________________
223 AliFemtoString AliFemtoV0PairCut::Report(){
224   // Prepare the report from the execution
225   string stemp = "AliFemtoV0 Pair Cut - remove shared and split pairs\n";  char ctemp[100];
226   snprintf(ctemp , 100, "Number of pairs which passed:\t%ld  Number which failed:\t%ld\n",fNPairsPassed,fNPairsFailed);
227   stemp += ctemp;
228   AliFemtoString returnThis = stemp;
229   return returnThis;}
230 //__________________
231
232 void AliFemtoV0PairCut::SetV0Max(Double_t aV0Max) {
233   fV0Max = aV0Max;
234 }
235
236 Double_t AliFemtoV0PairCut::GetAliFemtoV0Max() const {
237   return fV0Max;
238 }
239
240
241 TList *AliFemtoV0PairCut::ListSettings()
242 {
243   // return a list of settings in a writable form
244   TList *tListSetttings = new TList();
245   char buf[200];
246   snprintf(buf, 200, "AliFemtoV0PairCut.sharequalitymax=%f", fV0Max);
247   snprintf(buf, 200, "AliFemtoV0PairCut.sharefractionmax=%f", fShareFractionMax);
248   tListSetttings->AddLast(new TObjString(buf));
249
250   return tListSetttings;
251 }
252
253 void     AliFemtoV0PairCut::SetRemoveSameLabel(Bool_t aRemove)
254 {
255   fRemoveSameLabel = aRemove;
256 }
257
258 void AliFemtoV0PairCut::SetDataType(AliFemtoDataType type)
259 {
260   fDataType = type;
261 }
262
263 void AliFemtoV0PairCut::SetTPCEntranceSepMinimum(double dtpc)
264 {
265   fDTPCMin = dtpc;
266 }
267
268 void AliFemtoV0PairCut::SetTPCExitSepMinimum(double dtpc)
269 {
270   fDTPCExitMin = dtpc;
271 }
272
273 void AliFemtoV0PairCut::SetMinAvgSeparation(int type, double minSep)
274 {
275   if(type == 0) //Pos-Pos
276     fMinAvgSepPosPos = minSep;
277   else if(type == 1) //Pos-Neg
278     fMinAvgSepPosNeg = minSep;
279   else if(type == 2) //Neg-Pos
280     fMinAvgSepNegPos = minSep;
281   else if(type == 3) //Neg-Neg
282     fMinAvgSepNegNeg = minSep;
283 }