]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/CaloTrackCorrBase/AliFiducialCut.cxx
Avoid names that only differ in use of capital/small letters
[u/mrichter/AliRoot.git] / PWG / CaloTrackCorrBase / AliFiducialCut.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
16 //_________________________________________________________________________
17 // Class for track/cluster acceptance selection
18 // Selection in Central barrel, EMCAL and PHOS
19 //  
20 // Several selection regions possible for the different
21 // detectors
22 //
23 //*-- Author: Gustavo Conesa (LNF-INFN) 
24 //////////////////////////////////////////////////////////////////////////////
25
26
27 // --- ROOT system ---
28 #include <TMath.h>
29 #include <TLorentzVector.h>
30 #include <TString.h>
31
32 //---- ANALYSIS system ----
33 #include "AliFiducialCut.h"
34
35 ClassImp(AliFiducialCut)
36
37
38 //________________________________
39 AliFiducialCut::AliFiducialCut() : 
40 TObject(),
41 fEMCALFiducialCut(0),   fPHOSFiducialCut(0),    fCTSFiducialCut(0),
42 fCTSFidCutMinEta(0x0),  fCTSFidCutMinPhi(0x0),  fCTSFidCutMaxEta(0x0),   fCTSFidCutMaxPhi(0x0),
43 fEMCALFidCutMinEta(0x0),fEMCALFidCutMinPhi(0x0),fEMCALFidCutMaxEta(0x0), fEMCALFidCutMaxPhi(0x0),
44 fPHOSFidCutMinEta(0x0), fPHOSFidCutMinPhi(0x0), fPHOSFidCutMaxEta(0x0),  fPHOSFidCutMaxPhi(0x0)
45
46 {
47   //Ctor
48   
49   //Initialize parameters
50   InitParameters();
51   
52 }
53
54 //_______________________________
55 AliFiducialCut::~AliFiducialCut()
56 {
57   //Dtor
58   
59   if(fCTSFidCutMinEta)   delete fCTSFidCutMinEta ;
60   if(fCTSFidCutMinPhi)   delete fCTSFidCutMinPhi ;
61   if(fCTSFidCutMaxEta)   delete fCTSFidCutMaxEta ;
62   if(fCTSFidCutMaxPhi)   delete fCTSFidCutMaxPhi ;
63   
64   if(fEMCALFidCutMinEta) delete fEMCALFidCutMinEta ;
65   if(fEMCALFidCutMinPhi) delete fEMCALFidCutMinPhi ;
66   if(fEMCALFidCutMaxEta) delete fEMCALFidCutMaxEta ; 
67   if(fEMCALFidCutMaxPhi) delete fEMCALFidCutMaxPhi ;
68   
69   if(fPHOSFidCutMinEta)  delete fPHOSFidCutMinEta ; 
70   if(fPHOSFidCutMinPhi)  delete fPHOSFidCutMinPhi ; 
71   if(fPHOSFidCutMaxEta)  delete fPHOSFidCutMaxEta ;
72   if(fPHOSFidCutMaxPhi)  delete fPHOSFidCutMaxPhi ;
73   
74 }
75
76
77 //________________________________________________________________________________
78 Bool_t AliFiducialCut::IsInFiducialCut(TLorentzVector momentum, TString det) const
79 {
80   //Selects EMCAL or PHOS cluster or CTS track if it is inside eta-phi defined regions
81   
82   if(det == "CTS")
83   {
84           if(!fCTSFiducialCut)  
85       return kTRUE; //Fiducial cut not requested, accept all tracks  
86           else 
87       return CheckFiducialRegion(momentum, fCTSFidCutMinPhi  , fCTSFidCutMaxPhi , fCTSFidCutMinEta  , fCTSFidCutMaxEta  );
88   }
89   else   if(det == "EMCAL") 
90   {
91           if(!fEMCALFiducialCut) 
92       return kTRUE; //Fiducial cut not requested, accept all clusters  
93           else                   
94       return CheckFiducialRegion(momentum, fEMCALFidCutMinPhi, fEMCALFidCutMaxPhi, fEMCALFidCutMinEta, fEMCALFidCutMaxEta);
95   }
96   else   if(det == "PHOS")  
97   {
98           if(!fPHOSFiducialCut) 
99       return kTRUE; //Fiducial cut not requested, accept all clusters 
100           else 
101       return CheckFiducialRegion(momentum, fPHOSFidCutMinPhi , fPHOSFidCutMaxPhi , fPHOSFidCutMinEta , fPHOSFidCutMaxEta );
102   }
103   else
104   {
105           printf("AliFiducialCut::IsInFiducialCut() - Wrong detector name = %s\n", det.Data());
106           return kFALSE;
107   }
108   
109 }
110
111 //___________________________________________________________________________________________
112 Bool_t AliFiducialCut::CheckFiducialRegion(TLorentzVector momentum,
113                                            const TArrayF* minphi, const TArrayF* maxphi, 
114                                            const TArrayF* mineta, const TArrayF* maxeta) const 
115 {
116   //Given the selection regions in Eta and Phi, check if particle is in this region.
117   
118   Double_t phi = momentum.Phi();
119         if(phi < 0) phi+=TMath::TwoPi() ;
120         Double_t eta =  momentum.Eta();
121         //printf("IsInFiducialCut::Det: %s, phi = %f, eta = %f\n", det.Data(),phi*TMath::RadToDeg(), eta);
122   
123   Int_t netaregions = maxeta->GetSize();
124   Int_t nphiregions = maxphi->GetSize();
125   if(netaregions !=  mineta->GetSize() || nphiregions !=  minphi->GetSize())
126                 printf("AliFiducialCut::IsInFiducialCut() - Wrong number of fiducial cut regions: nmaxeta %d != nmineta %d; nmaxphi %d != nminphi %d\n",
127            netaregions, mineta->GetSize(),  nphiregions, minphi->GetSize());
128         
129         //Eta fiducial cut
130         Bool_t bInEtaFidCut = kFALSE;
131         for(Int_t ieta = 0; ieta < netaregions; ieta++)
132                 if(eta > mineta->GetAt(ieta) && eta < maxeta->GetAt(ieta)) bInEtaFidCut = kTRUE;
133   
134         if(bInEtaFidCut){
135                 //printf("Eta cut passed\n");
136                 //Phi fiducial cut
137                 Bool_t bInPhiFidCut = kFALSE;
138                 for(Int_t iphi = 0; iphi < nphiregions; iphi++)
139                         if(phi > minphi->GetAt(iphi) *TMath::DegToRad()&& phi < maxphi->GetAt(iphi)*TMath::DegToRad()) bInPhiFidCut = kTRUE ;
140           
141                 if(bInPhiFidCut) {
142                         //printf("IsInFiducialCut:: %s cluster/track accepted\n",det.Data());
143                         return kTRUE;
144                 }
145                 else return kFALSE;
146     
147         }//In eta fid cut
148         else
149     return kFALSE;
150 }
151
152
153 //_______________________________________________________________
154 void AliFiducialCut::InitParameters()
155 {
156   
157   //Initialize the parameters of the analysis.
158   
159   fEMCALFiducialCut = kTRUE ;  
160   fPHOSFiducialCut  = kTRUE ;
161   fCTSFiducialCut   = kTRUE ;
162   
163   fCTSFidCutMinEta = new TArrayF(1);
164   fCTSFidCutMinEta->SetAt(-0.9,0); 
165   fCTSFidCutMaxEta = new TArrayF(1);
166   fCTSFidCutMaxEta->SetAt( 0.9,0); 
167   
168   fCTSFidCutMinPhi = new TArrayF(1);
169   fCTSFidCutMinPhi->SetAt(0.  ,0); 
170   fCTSFidCutMaxPhi = new TArrayF(1);
171   fCTSFidCutMaxPhi->SetAt(360.,0); 
172   
173   fEMCALFidCutMinEta = new TArrayF(1);
174   fEMCALFidCutMinEta->SetAt(-0.7,0); 
175   fEMCALFidCutMaxEta = new TArrayF(1);
176   fEMCALFidCutMaxEta->SetAt( 0.7,0); 
177   
178   fEMCALFidCutMinPhi = new TArrayF(1);
179   fEMCALFidCutMinPhi->SetAt(80.,0); 
180   fEMCALFidCutMaxPhi = new TArrayF(1);
181   fEMCALFidCutMaxPhi->SetAt(187.,0); 
182   
183   fPHOSFidCutMinEta = new TArrayF(1);
184   fPHOSFidCutMinEta->SetAt(-0.13,0); 
185   fPHOSFidCutMaxEta = new TArrayF(1);
186   fPHOSFidCutMaxEta->SetAt( 0.13,0); 
187   
188   fPHOSFidCutMinPhi = new TArrayF(1);
189   fPHOSFidCutMinPhi->SetAt(260.,0); 
190   fPHOSFidCutMaxPhi = new TArrayF(1);
191   fPHOSFidCutMaxPhi->SetAt(320.,0); 
192   
193 }
194
195
196 //________________________________________________________________
197 void AliFiducialCut::Print(const Option_t * opt) const
198 {
199   
200   //Print some relevant parameters set for the analysis
201   if(! opt)
202     return;
203   
204   printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
205   
206   if(fCTSFiducialCut)
207   {
208     Int_t netaregions =  fCTSFidCutMaxEta->GetSize();
209     Int_t nphiregions =  fCTSFidCutMaxPhi->GetSize();
210     printf(">> CTS Fiducial regions : phi %d eta %d\n", netaregions, nphiregions) ;
211     for(Int_t ieta = 0; ieta < netaregions; ieta++)
212       printf(" region %d : %3.2f < eta < %3.2f\n", ieta, fCTSFidCutMinEta->GetAt(ieta), fCTSFidCutMaxEta->GetAt(ieta)) ;
213     for(Int_t iphi = 0; iphi < nphiregions; iphi++)
214       printf(" region %d : %3.1f < phi < %3.1f\n", iphi, fCTSFidCutMinPhi->GetAt(iphi), fCTSFidCutMaxPhi->GetAt(iphi)) ; 
215   }
216   else printf(">>No fiducial cuts in CTS\n");
217   
218   if(fEMCALFiducialCut)
219   {
220     Int_t netaregions =  fEMCALFidCutMaxEta->GetSize();
221     Int_t nphiregions =  fEMCALFidCutMaxPhi->GetSize();
222     printf(">>EMCAL Fiducial regions : phi %d eta %d\n", netaregions, nphiregions) ;
223     for(Int_t ieta = 0; ieta < netaregions; ieta++)
224       printf(" region %d : %3.2f < eta < %3.2f\n", ieta, fEMCALFidCutMinEta->GetAt(ieta), fEMCALFidCutMaxEta->GetAt(ieta)) ;
225     for(Int_t iphi = 0; iphi < nphiregions; iphi++)
226       printf(" region %d : %3.1f < phi < %3.1f\n", iphi, fEMCALFidCutMinPhi->GetAt(iphi), fEMCALFidCutMaxPhi->GetAt(iphi)) ; 
227   }
228   else printf(">>No fiducial cuts in EMCAL\n");
229   
230   if(fPHOSFiducialCut)
231   {
232     Int_t netaregions =  fPHOSFidCutMaxEta->GetSize();
233     Int_t nphiregions =  fPHOSFidCutMaxPhi->GetSize();
234     printf(">>PHOS Fiducial regions : phi %d eta %d\n", netaregions, nphiregions) ;
235     for(Int_t ieta = 0; ieta < netaregions; ieta++)
236       printf(" region %d : %3.2f < eta < %3.2f\n", ieta, fPHOSFidCutMinEta->GetAt(ieta), fPHOSFidCutMaxEta->GetAt(ieta)) ;
237     for(Int_t iphi = 0; iphi < nphiregions; iphi++)
238       printf(" region %d : %3.1f < phi < %3.1f\n", iphi, fPHOSFidCutMinPhi->GetAt(iphi), fPHOSFidCutMaxPhi->GetAt(iphi)) ; 
239   }
240   else printf(">>No fiducial cuts in PHOS\n");
241   printf("    \n") ;
242   
243
244
245 //_______________________________________________________________________________________
246 void AliFiducialCut::SetSimpleCTSFiducialCut(Float_t eta, Float_t minphi, Float_t maxphi)
247 {
248   
249   //Method to set simple acceptance cut to CTS
250   
251   fCTSFidCutMinEta->Set(1);
252   fCTSFidCutMaxEta->Set(1);
253   fCTSFidCutMinPhi->Set(1);
254   fCTSFidCutMaxPhi->Set(1);
255   
256   fCTSFidCutMinEta->SetAt(-eta,0);
257   fCTSFidCutMaxEta->SetAt( eta,0);
258   fCTSFidCutMinPhi->SetAt(minphi,0);
259   fCTSFidCutMaxPhi->SetAt(maxphi,0);
260   
261 }
262
263 //_________________________________________________________________________________________
264 void AliFiducialCut::SetSimpleEMCALFiducialCut(Float_t eta, Float_t minphi, Float_t maxphi)
265 {
266   //Method to set simple acceptance cut to EMCAL
267   
268   fEMCALFidCutMinEta->Set(1);
269   fEMCALFidCutMaxEta->Set(1);
270   fEMCALFidCutMinPhi->Set(1);
271   fEMCALFidCutMaxPhi->Set(1);
272   
273   fEMCALFidCutMinEta->SetAt(-eta,0);
274   fEMCALFidCutMaxEta->SetAt( eta,0);
275   fEMCALFidCutMinPhi->SetAt(minphi,0);
276   fEMCALFidCutMaxPhi->SetAt(maxphi,0);
277   
278 }
279
280 //________________________________________________________________________________________
281 void AliFiducialCut::SetSimplePHOSFiducialCut(Float_t eta, Float_t minphi, Float_t maxphi)
282 {
283   //Method to set simple acceptance cut to PHOS
284   
285   fPHOSFidCutMinEta->Set(1);
286   fPHOSFidCutMaxEta->Set(1);
287   fPHOSFidCutMinPhi->Set(1);
288   fPHOSFidCutMaxPhi->Set(1);
289   
290   fPHOSFidCutMinEta->SetAt(-eta,0);
291   fPHOSFidCutMaxEta->SetAt(eta,0);
292   fPHOSFidCutMinPhi->SetAt(minphi,0);
293   fPHOSFidCutMaxPhi->SetAt(maxphi,0);
294   
295 }