]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/PartCorrBase/AliFidutialCut.cxx
In case of merging AOD inputs, recover vertex from both inputs, not from first.
[u/mrichter/AliRoot.git] / PWG4 / PartCorrBase / AliFidutialCut.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15 /* $Id: AliFidutialCut.cxx 21839 2007-10-29 13:49:42Z gustavo $ */
16
17 //_________________________________________________________________________
18 // Class for track/cluster acceptance selection
19 // Selection in Central barrel, EMCAL and PHOS
20 //                
21 //*-- Author: Gustavo Conesa (LNF-INFN) 
22 //////////////////////////////////////////////////////////////////////////////
23
24
25 // --- ROOT system ---
26 #include <TMath.h>
27 #include <TLorentzVector.h>
28 #include <TString.h>
29 //#include <TArrayF.h>
30   
31 //---- ANALYSIS system ----
32 #include "AliFidutialCut.h"
33
34 ClassImp(AliFidutialCut)
35
36
37 //____________________________________________________________________________
38 AliFidutialCut::AliFidutialCut() : 
39   TObject(),
40   fEMCALFidutialCut(0),  fPHOSFidutialCut(0),  fCTSFidutialCut(0),
41   fCTSFidCutMinEta(0x0),fCTSFidCutMinPhi(0x0),fCTSFidCutMaxEta(0x0), fCTSFidCutMaxPhi(0x0),
42   fEMCALFidCutMinEta(0x0),fEMCALFidCutMinPhi(0x0),fEMCALFidCutMaxEta(0x0), fEMCALFidCutMaxPhi(0x0),
43   fPHOSFidCutMinEta(0x0),fPHOSFidCutMinPhi(0x0),fPHOSFidCutMaxEta(0x0), fPHOSFidCutMaxPhi(0x0)
44
45 {
46   //Ctor
47
48   //Initialize parameters
49   InitParameters();
50 }
51
52 //____________________________________________________________________________
53 AliFidutialCut::AliFidutialCut(const AliFidutialCut & g) :   
54   TObject(g), 
55   fEMCALFidutialCut(g.fEMCALFidutialCut),  fPHOSFidutialCut(g.fPHOSFidutialCut),  fCTSFidutialCut(g. fCTSFidutialCut),
56   fCTSFidCutMinEta(g.fCTSFidCutMinEta?new TArrayF(*g.fCTSFidCutMinEta):0x0),
57   fCTSFidCutMinPhi(g.fCTSFidCutMinPhi?new TArrayF(*g.fCTSFidCutMinPhi):0x0),
58   fCTSFidCutMaxEta(g.fCTSFidCutMaxEta?new TArrayF(*g.fCTSFidCutMaxEta):0x0),
59   fCTSFidCutMaxPhi(g.fCTSFidCutMaxPhi?new TArrayF(*g.fCTSFidCutMaxPhi):0x0),
60   fEMCALFidCutMinEta(g.fEMCALFidCutMinEta?new TArrayF(*g.fEMCALFidCutMinEta):0x0),
61   fEMCALFidCutMinPhi(g.fEMCALFidCutMinPhi?new TArrayF(*g.fEMCALFidCutMinPhi):0x0),
62   fEMCALFidCutMaxEta(g.fEMCALFidCutMaxEta?new TArrayF(*g.fEMCALFidCutMaxEta):0x0),
63   fEMCALFidCutMaxPhi(g.fEMCALFidCutMaxPhi?new TArrayF(*g.fEMCALFidCutMaxPhi):0x0),
64   fPHOSFidCutMinEta(g.fPHOSFidCutMinEta?new TArrayF(*g.fPHOSFidCutMinEta):0x0),
65   fPHOSFidCutMinPhi(g.fPHOSFidCutMinPhi?new TArrayF(*g.fPHOSFidCutMinPhi):0x0),
66   fPHOSFidCutMaxEta(g.fPHOSFidCutMaxEta?new TArrayF(*g.fPHOSFidCutMaxEta):0x0),
67   fPHOSFidCutMaxPhi(g.fPHOSFidCutMaxPhi?new TArrayF(*g.fPHOSFidCutMaxPhi):0x0)
68
69 {
70   // cpy ctor
71
72 }
73
74 //_________________________________________________________________________
75 AliFidutialCut & AliFidutialCut::operator = (const AliFidutialCut & source)
76 {
77   // assignment operator
78   
79   if(&source == this) return *this;
80   
81   fEMCALFidutialCut = source.fEMCALFidutialCut;  
82   fPHOSFidutialCut = source.fPHOSFidutialCut;
83   fCTSFidutialCut = source.fCTSFidutialCut;
84   
85   fCTSFidCutMinEta = source.fCTSFidCutMinEta?new TArrayF(*source.fCTSFidCutMinEta):0x0;
86   fCTSFidCutMinPhi = source.fCTSFidCutMinPhi?new TArrayF(*source.fCTSFidCutMinPhi):0x0;
87   fCTSFidCutMaxEta = source.fCTSFidCutMaxEta?new TArrayF(*source.fCTSFidCutMaxEta):0x0;
88   fCTSFidCutMaxPhi = source.fCTSFidCutMaxPhi?new TArrayF(*source.fCTSFidCutMaxPhi):0x0;
89   fEMCALFidCutMinEta = source.fEMCALFidCutMinEta?new TArrayF(*source.fEMCALFidCutMinEta):0x0;
90   fEMCALFidCutMinPhi = source.fEMCALFidCutMinPhi?new TArrayF(*source.fEMCALFidCutMinPhi):0x0;
91   fEMCALFidCutMaxEta = source.fEMCALFidCutMaxEta?new TArrayF(*source.fEMCALFidCutMaxEta):0x0;
92   fEMCALFidCutMaxPhi = source.fEMCALFidCutMaxPhi?new TArrayF(*source.fEMCALFidCutMaxPhi):0x0;
93   fPHOSFidCutMinEta = source.fPHOSFidCutMinEta?new TArrayF(*source.fPHOSFidCutMinEta):0x0;
94   fPHOSFidCutMinPhi = source.fPHOSFidCutMinPhi?new TArrayF(*source.fPHOSFidCutMinPhi):0x0;
95   fPHOSFidCutMaxEta = source.fPHOSFidCutMaxEta?new TArrayF(*source.fPHOSFidCutMaxEta):0x0;
96   fPHOSFidCutMaxPhi = source.fPHOSFidCutMaxPhi?new TArrayF(*source.fPHOSFidCutMaxPhi):0x0;
97
98   return *this;
99
100 }
101
102 //_________________________________
103 AliFidutialCut::~AliFidutialCut() {
104   //Dtor
105
106   delete fCTSFidCutMinEta ;
107   delete fCTSFidCutMinPhi ;
108   delete fCTSFidCutMaxEta ;
109   delete fCTSFidCutMaxPhi ;
110   
111   delete fEMCALFidCutMinEta ;
112   delete fEMCALFidCutMinPhi ;
113   delete fEMCALFidCutMaxEta ; 
114   delete fEMCALFidCutMaxPhi ;
115   
116   delete fPHOSFidCutMinEta ; 
117   delete fPHOSFidCutMinPhi ; 
118   delete fPHOSFidCutMaxEta ;
119   delete fPHOSFidCutMaxPhi ;
120
121 }
122
123
124 //_______________________________________________________________
125 Bool_t AliFidutialCut::IsInFidutialCut(TLorentzVector momentum, TString det) const
126 {
127   //Selects EMCAL or PHOS cluster or CTS track if it is inside eta-phi defined regions
128
129   Double_t phi = momentum.Phi();
130   if(phi<0) phi+=TMath::TwoPi() ;
131   Double_t eta =  momentum.Eta();
132   //printf("IsInFidutialCut::Det: %s, phi = %f, eta = %f\n", det.Data(),phi*TMath::RadToDeg(), eta);
133
134   Float_t * maxeta = new Float_t;
135   Float_t * maxphi = new Float_t;
136   Float_t * mineta = new Float_t;
137   Float_t * minphi = new Float_t;
138   Int_t nphiregions = 0;
139   Int_t netaregions = 0;
140   Bool_t selection = kFALSE;
141
142   if(det == "CTS"){
143     netaregions =  fCTSFidCutMaxEta->GetSize();
144     nphiregions =  fCTSFidCutMaxPhi->GetSize();
145     if(netaregions !=  fCTSFidCutMinEta->GetSize() || nphiregions !=  fCTSFidCutMinPhi->GetSize())
146       printf("AliFidutialCut::IsInFidutialCut() - Wrong number of CTS fidutial cut regions: nmaxeta %d != nmineta %d; nmaxphi %d != nminphi %d\n",
147                     netaregions, fCTSFidCutMinEta->GetSize(),  nphiregions, fCTSFidCutMinPhi->GetSize());
148     
149     maxeta = fCTSFidCutMaxEta->GetArray();
150     maxphi = fCTSFidCutMaxPhi->GetArray();
151     mineta = fCTSFidCutMinEta->GetArray();
152     minphi = fCTSFidCutMinPhi->GetArray();
153     selection =  fCTSFidutialCut ; 
154   }
155   else   if(det == "EMCAL"){
156     netaregions =  fEMCALFidCutMaxEta->GetSize();
157     nphiregions =  fEMCALFidCutMaxPhi->GetSize();
158     if(netaregions !=  fEMCALFidCutMinEta->GetSize() || nphiregions !=  fEMCALFidCutMinPhi->GetSize())
159       printf("AliFidutialCut::IsInFidutialCut() - Wrong number of EMCAL fidutial cut regions: nmaxeta %d != nmineta %d; nmaxphi %d != nminphi %d\n",
160                     netaregions, fEMCALFidCutMinEta->GetSize(),  nphiregions, fEMCALFidCutMinPhi->GetSize());
161     
162     maxeta = fEMCALFidCutMaxEta->GetArray();
163     maxphi = fEMCALFidCutMaxPhi->GetArray();
164     mineta = fEMCALFidCutMinEta->GetArray();
165     minphi = fEMCALFidCutMinPhi->GetArray();
166     selection =  fEMCALFidutialCut ; 
167   }
168   else   if(det == "PHOS"){
169     netaregions =  fPHOSFidCutMaxEta->GetSize();
170     nphiregions =  fPHOSFidCutMaxPhi->GetSize();
171     if(netaregions !=  fPHOSFidCutMinEta->GetSize() || nphiregions !=  fPHOSFidCutMinPhi->GetSize())
172       printf("AliFidutialCut::IsInFidutialCut() - Wrong number of PHOS fidutial cut regions: nmaxeta %d != nmineta %d; nmaxphi %d != nminphi %d\n",
173                     netaregions, fPHOSFidCutMinEta->GetSize(),  nphiregions, fPHOSFidCutMinPhi->GetSize());
174     
175     maxeta = fPHOSFidCutMaxEta->GetArray();
176     maxphi = fPHOSFidCutMaxPhi->GetArray();
177     mineta = fPHOSFidCutMinEta->GetArray();
178     minphi = fPHOSFidCutMinPhi->GetArray();
179     selection =  fPHOSFidutialCut ; 
180   }
181   else
182     printf("AliFidutialCut::IsInFidutialCut() - Wrong detector name = %s\n", det.Data());
183
184   //printf("IsInFidutialCut::nphiregions = %d, netaregions = %d\n",nphiregions, netaregions);
185
186   if(!selection) return kTRUE; //No cuts applied, all tracks/clusters used
187
188   //Eta fidutial cut
189   Bool_t bInEtaFidCut = kFALSE;
190   for(Int_t ieta = 0; ieta < netaregions; ieta++)
191     if(eta > mineta[ieta] && eta < maxeta[ieta]) bInEtaFidCut = kTRUE;
192
193   if(bInEtaFidCut){
194     //printf("Eta cut passed\n");
195     //Phi fidutial cut
196     Bool_t bInPhiFidCut = kFALSE;
197     for(Int_t iphi = 0; iphi < nphiregions; iphi++)
198       if(phi > minphi[iphi] *TMath::DegToRad()&& phi < maxphi[iphi]*TMath::DegToRad()) bInPhiFidCut = kTRUE ;
199     
200     if(bInPhiFidCut) {
201       //printf("IsInFidutialCut:: %s cluster/track accepted\n",det.Data());
202       return kTRUE;
203     }
204     else return kFALSE;
205     
206   }//In eta fid cut
207   else
208     return kFALSE;
209   
210 }
211
212
213 //_______________________________________________________________
214 void AliFidutialCut::InitParameters()
215 {
216  
217   //Initialize the parameters of the analysis.
218
219   fEMCALFidutialCut = kTRUE ;  
220   fPHOSFidutialCut = kTRUE ;
221   fCTSFidutialCut = kTRUE ;
222
223   fCTSFidCutMinEta = new TArrayF(1);
224   fCTSFidCutMinEta->SetAt(-0.9,0); 
225   fCTSFidCutMaxEta = new TArrayF(1);
226   fCTSFidCutMaxEta->SetAt(0.9,0); 
227
228   fCTSFidCutMinPhi = new TArrayF(1);
229   fCTSFidCutMinPhi->SetAt(0,0); 
230   fCTSFidCutMaxPhi = new TArrayF(1);
231   fCTSFidCutMaxPhi->SetAt(360.,0); 
232
233   fEMCALFidCutMinEta = new TArrayF(1);
234   fEMCALFidCutMinEta->SetAt(-0.7,0); 
235   fEMCALFidCutMaxEta = new TArrayF(1);
236   fEMCALFidCutMaxEta->SetAt(0.7,0); 
237
238   fEMCALFidCutMinPhi = new TArrayF(1);
239   fEMCALFidCutMinPhi->SetAt(80.,0); 
240   fEMCALFidCutMaxPhi = new TArrayF(1);
241   fEMCALFidCutMaxPhi->SetAt(190.,0); 
242
243   fPHOSFidCutMinEta = new TArrayF(1);
244   fPHOSFidCutMinEta->SetAt(-0.13,0); 
245   fPHOSFidCutMaxEta = new TArrayF(1);
246   fPHOSFidCutMaxEta->SetAt(0.13,0); 
247
248   fPHOSFidCutMinPhi = new TArrayF(1);
249   fPHOSFidCutMinPhi->SetAt(220.,0); 
250   fPHOSFidCutMaxPhi = new TArrayF(1);
251   fPHOSFidCutMaxPhi->SetAt(320.,0); 
252
253 }
254
255
256 //________________________________________________________________
257 void AliFidutialCut::Print(const Option_t * opt) const
258 {
259
260   //Print some relevant parameters set for the analysis
261   if(! opt)
262     return;
263
264   printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
265
266   if(fCTSFidutialCut){
267     Int_t netaregions =  fCTSFidCutMaxEta->GetSize();
268     Int_t nphiregions =  fCTSFidCutMaxPhi->GetSize();
269     printf(">> CTS Fidutial regions : phi %d eta %d\n", netaregions, nphiregions) ;
270     for(Int_t ieta = 0; ieta < netaregions; ieta++)
271       printf(" region %d : %3.2f < eta < %3.2f\n", ieta, fCTSFidCutMinEta->GetAt(ieta), fCTSFidCutMaxEta->GetAt(ieta)) ;
272     for(Int_t iphi = 0; iphi < nphiregions; iphi++)
273       printf(" region %d : %3.1f < phi < %3.1f\n", iphi, fCTSFidCutMinPhi->GetAt(iphi), fCTSFidCutMaxPhi->GetAt(iphi)) ; 
274   }
275   else printf(">>No fidutial cuts in CTS\n");
276
277   if(fEMCALFidutialCut){
278     Int_t netaregions =  fEMCALFidCutMaxEta->GetSize();
279     Int_t nphiregions =  fEMCALFidCutMaxPhi->GetSize();
280     printf(">>EMCAL Fidutial regions : phi %d eta %d\n", netaregions, nphiregions) ;
281     for(Int_t ieta = 0; ieta < netaregions; ieta++)
282       printf(" region %d : %3.2f < eta < %3.2f\n", ieta, fEMCALFidCutMinEta->GetAt(ieta), fEMCALFidCutMaxEta->GetAt(ieta)) ;
283     for(Int_t iphi = 0; iphi < nphiregions; iphi++)
284       printf(" region %d : %3.1f < phi < %3.1f\n", iphi, fEMCALFidCutMinPhi->GetAt(iphi), fEMCALFidCutMaxPhi->GetAt(iphi)) ; 
285   }
286   else printf(">>No fidutial cuts in EMCAL\n");
287
288   if(fPHOSFidutialCut){
289     Int_t netaregions =  fPHOSFidCutMaxEta->GetSize();
290     Int_t nphiregions =  fPHOSFidCutMaxPhi->GetSize();
291     printf(">>PHOS Fidutial regions : phi %d eta %d\n", netaregions, nphiregions) ;
292     for(Int_t ieta = 0; ieta < netaregions; ieta++)
293       printf(" region %d : %3.2f < eta < %3.2f\n", ieta, fPHOSFidCutMinEta->GetAt(ieta), fPHOSFidCutMaxEta->GetAt(ieta)) ;
294     for(Int_t iphi = 0; iphi < nphiregions; iphi++)
295       printf(" region %d : %3.1f < phi < %3.1f\n", iphi, fPHOSFidCutMinPhi->GetAt(iphi), fPHOSFidCutMaxPhi->GetAt(iphi)) ; 
296   }
297   else printf(">>No fidutial cuts in PHOS\n");
298   printf("    \n") ;
299
300
301
302 //_______________________________________________________________
303 void AliFidutialCut::SetSimpleCTSFidutialCut(const Float_t eta, const Float_t minphi, const Float_t maxphi){
304
305   //Method to set simple acceptance cut to CTS
306   fCTSFidCutMinEta->Set(1);
307   fCTSFidCutMaxEta->Set(1);
308   fCTSFidCutMinPhi->Set(1);
309   fCTSFidCutMaxPhi->Set(1);
310
311   fCTSFidCutMinEta->SetAt(-eta,0);
312   fCTSFidCutMaxEta->SetAt(eta,0);
313   fCTSFidCutMinPhi->SetAt(minphi,0);
314   fCTSFidCutMaxPhi->SetAt(maxphi,0);
315
316
317 }
318
319
320 //_______________________________________________________________
321 void AliFidutialCut::SetSimpleEMCALFidutialCut(const Float_t eta, const Float_t minphi, const Float_t maxphi){
322   //Method to set simple acceptance cut to EMCAL
323
324   fEMCALFidCutMinEta->Set(1);
325   fEMCALFidCutMaxEta->Set(1);
326   fEMCALFidCutMinPhi->Set(1);
327   fEMCALFidCutMaxPhi->Set(1);
328
329   fEMCALFidCutMinEta->SetAt(-eta,0);
330   fEMCALFidCutMaxEta->SetAt(eta,0);
331   fEMCALFidCutMinPhi->SetAt(minphi,0);
332   fEMCALFidCutMaxPhi->SetAt(maxphi,0);
333
334
335 }
336
337 //_______________________________________________________________
338 void AliFidutialCut::SetSimplePHOSFidutialCut(const Float_t eta, const Float_t minphi, const Float_t maxphi){
339
340   //Method to set simple acceptance cut to PHOS
341   fPHOSFidCutMinEta->Set(1);
342   fPHOSFidCutMaxEta->Set(1);
343   fPHOSFidCutMinPhi->Set(1);
344   fPHOSFidCutMaxPhi->Set(1);
345
346   fPHOSFidCutMinEta->SetAt(-eta,0);
347   fPHOSFidCutMaxEta->SetAt(eta,0);
348   fPHOSFidCutMinPhi->SetAt(minphi,0);
349   fPHOSFidCutMaxPhi->SetAt(maxphi,0);
350
351
352 }