]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/CaloTrackCorrBase/AliFiducialCut.cxx
AddTask added
[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(const TLorentzVector momentum, 
79                                        const TString det) const
80 {
81   //Selects EMCAL or PHOS cluster or CTS track if it is inside eta-phi defined regions
82   
83   if(det == "CTS")
84   {
85           if(!fCTSFiducialCut)  
86       return kTRUE; //Fiducial cut not requested, accept all tracks  
87           else 
88       return CheckFiducialRegion(momentum, fCTSFidCutMinPhi  , fCTSFidCutMaxPhi , fCTSFidCutMinEta  , fCTSFidCutMaxEta  );
89   }
90   else   if(det == "EMCAL") 
91   {
92           if(!fEMCALFiducialCut) 
93       return kTRUE; //Fiducial cut not requested, accept all clusters  
94           else                   
95       return CheckFiducialRegion(momentum, fEMCALFidCutMinPhi, fEMCALFidCutMaxPhi, fEMCALFidCutMinEta, fEMCALFidCutMaxEta);
96   }
97   else   if(det == "PHOS")  
98   {
99           if(!fPHOSFiducialCut) 
100       return kTRUE; //Fiducial cut not requested, accept all clusters 
101           else 
102       return CheckFiducialRegion(momentum, fPHOSFidCutMinPhi , fPHOSFidCutMaxPhi , fPHOSFidCutMinEta , fPHOSFidCutMaxEta );
103   }
104   else
105   {
106           printf("AliFiducialCut::IsInFiducialCut() - Wrong detector name = %s\n", det.Data());
107           return kFALSE;
108   }
109   
110 }
111
112 //___________________________________________________________________________________________
113 Bool_t AliFiducialCut::CheckFiducialRegion(const TLorentzVector momentum, 
114                                            const TArrayF* minphi, const TArrayF* maxphi, 
115                                            const TArrayF* mineta, const TArrayF* maxeta) const 
116 {
117   //Given the selection regions in Eta and Phi, check if particle is in this region.
118   
119   Double_t phi = momentum.Phi();
120         if(phi < 0) phi+=TMath::TwoPi() ;
121         Double_t eta =  momentum.Eta();
122         //printf("IsInFiducialCut::Det: %s, phi = %f, eta = %f\n", det.Data(),phi*TMath::RadToDeg(), eta);
123   
124   Int_t netaregions = maxeta->GetSize();
125   Int_t nphiregions = maxphi->GetSize();
126   if(netaregions !=  mineta->GetSize() || nphiregions !=  minphi->GetSize())
127                 printf("AliFiducialCut::IsInFiducialCut() - Wrong number of fiducial cut regions: nmaxeta %d != nmineta %d; nmaxphi %d != nminphi %d\n",
128            netaregions, mineta->GetSize(),  nphiregions, minphi->GetSize());
129         
130         //Eta fiducial cut
131         Bool_t bInEtaFidCut = kFALSE;
132         for(Int_t ieta = 0; ieta < netaregions; ieta++)
133                 if(eta > mineta->GetAt(ieta) && eta < maxeta->GetAt(ieta)) bInEtaFidCut = kTRUE;
134   
135         if(bInEtaFidCut){
136                 //printf("Eta cut passed\n");
137                 //Phi fiducial cut
138                 Bool_t bInPhiFidCut = kFALSE;
139                 for(Int_t iphi = 0; iphi < nphiregions; iphi++)
140                         if(phi > minphi->GetAt(iphi) *TMath::DegToRad()&& phi < maxphi->GetAt(iphi)*TMath::DegToRad()) bInPhiFidCut = kTRUE ;
141           
142                 if(bInPhiFidCut) {
143                         //printf("IsInFiducialCut:: %s cluster/track accepted\n",det.Data());
144                         return kTRUE;
145                 }
146                 else return kFALSE;
147     
148         }//In eta fid cut
149         else
150     return kFALSE;
151 }
152
153
154 //_______________________________________________________________
155 void AliFiducialCut::InitParameters()
156 {
157   
158   //Initialize the parameters of the analysis.
159   
160   fEMCALFiducialCut = kTRUE ;  
161   fPHOSFiducialCut  = kTRUE ;
162   fCTSFiducialCut   = kTRUE ;
163   
164   fCTSFidCutMinEta = new TArrayF(1);
165   fCTSFidCutMinEta->SetAt(-0.9,0); 
166   fCTSFidCutMaxEta = new TArrayF(1);
167   fCTSFidCutMaxEta->SetAt( 0.9,0); 
168   
169   fCTSFidCutMinPhi = new TArrayF(1);
170   fCTSFidCutMinPhi->SetAt(0.  ,0); 
171   fCTSFidCutMaxPhi = new TArrayF(1);
172   fCTSFidCutMaxPhi->SetAt(360.,0); 
173   
174   fEMCALFidCutMinEta = new TArrayF(1);
175   fEMCALFidCutMinEta->SetAt(-0.7,0); 
176   fEMCALFidCutMaxEta = new TArrayF(1);
177   fEMCALFidCutMaxEta->SetAt( 0.7,0); 
178   
179   fEMCALFidCutMinPhi = new TArrayF(1);
180   fEMCALFidCutMinPhi->SetAt(80.,0); 
181   fEMCALFidCutMaxPhi = new TArrayF(1);
182   fEMCALFidCutMaxPhi->SetAt(187.,0); 
183   
184   fPHOSFidCutMinEta = new TArrayF(1);
185   fPHOSFidCutMinEta->SetAt(-0.13,0); 
186   fPHOSFidCutMaxEta = new TArrayF(1);
187   fPHOSFidCutMaxEta->SetAt( 0.13,0); 
188   
189   fPHOSFidCutMinPhi = new TArrayF(1);
190   fPHOSFidCutMinPhi->SetAt(260.,0); 
191   fPHOSFidCutMaxPhi = new TArrayF(1);
192   fPHOSFidCutMaxPhi->SetAt(320.,0); 
193   
194 }
195
196
197 //________________________________________________________________
198 void AliFiducialCut::Print(const Option_t * opt) const
199 {
200   
201   //Print some relevant parameters set for the analysis
202   if(! opt)
203     return;
204   
205   printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
206   
207   if(fCTSFiducialCut)
208   {
209     Int_t netaregions =  fCTSFidCutMaxEta->GetSize();
210     Int_t nphiregions =  fCTSFidCutMaxPhi->GetSize();
211     printf(">> CTS Fiducial regions : phi %d eta %d\n", netaregions, nphiregions) ;
212     for(Int_t ieta = 0; ieta < netaregions; ieta++)
213       printf(" region %d : %3.2f < eta < %3.2f\n", ieta, fCTSFidCutMinEta->GetAt(ieta), fCTSFidCutMaxEta->GetAt(ieta)) ;
214     for(Int_t iphi = 0; iphi < nphiregions; iphi++)
215       printf(" region %d : %3.1f < phi < %3.1f\n", iphi, fCTSFidCutMinPhi->GetAt(iphi), fCTSFidCutMaxPhi->GetAt(iphi)) ; 
216   }
217   else printf(">>No fiducial cuts in CTS\n");
218   
219   if(fEMCALFiducialCut)
220   {
221     Int_t netaregions =  fEMCALFidCutMaxEta->GetSize();
222     Int_t nphiregions =  fEMCALFidCutMaxPhi->GetSize();
223     printf(">>EMCAL Fiducial regions : phi %d eta %d\n", netaregions, nphiregions) ;
224     for(Int_t ieta = 0; ieta < netaregions; ieta++)
225       printf(" region %d : %3.2f < eta < %3.2f\n", ieta, fEMCALFidCutMinEta->GetAt(ieta), fEMCALFidCutMaxEta->GetAt(ieta)) ;
226     for(Int_t iphi = 0; iphi < nphiregions; iphi++)
227       printf(" region %d : %3.1f < phi < %3.1f\n", iphi, fEMCALFidCutMinPhi->GetAt(iphi), fEMCALFidCutMaxPhi->GetAt(iphi)) ; 
228   }
229   else printf(">>No fiducial cuts in EMCAL\n");
230   
231   if(fPHOSFiducialCut)
232   {
233     Int_t netaregions =  fPHOSFidCutMaxEta->GetSize();
234     Int_t nphiregions =  fPHOSFidCutMaxPhi->GetSize();
235     printf(">>PHOS Fiducial regions : phi %d eta %d\n", netaregions, nphiregions) ;
236     for(Int_t ieta = 0; ieta < netaregions; ieta++)
237       printf(" region %d : %3.2f < eta < %3.2f\n", ieta, fPHOSFidCutMinEta->GetAt(ieta), fPHOSFidCutMaxEta->GetAt(ieta)) ;
238     for(Int_t iphi = 0; iphi < nphiregions; iphi++)
239       printf(" region %d : %3.1f < phi < %3.1f\n", iphi, fPHOSFidCutMinPhi->GetAt(iphi), fPHOSFidCutMaxPhi->GetAt(iphi)) ; 
240   }
241   else printf(">>No fiducial cuts in PHOS\n");
242   printf("    \n") ;
243   
244
245
246 //_______________________________________________________________
247 void AliFiducialCut::SetSimpleCTSFiducialCut(const Float_t eta, 
248                                              const Float_t minphi, 
249                                              const Float_t maxphi)
250 {
251   
252   //Method to set simple acceptance cut to CTS
253   
254   fCTSFidCutMinEta->Set(1);
255   fCTSFidCutMaxEta->Set(1);
256   fCTSFidCutMinPhi->Set(1);
257   fCTSFidCutMaxPhi->Set(1);
258   
259   fCTSFidCutMinEta->SetAt(-eta,0);
260   fCTSFidCutMaxEta->SetAt( eta,0);
261   fCTSFidCutMinPhi->SetAt(minphi,0);
262   fCTSFidCutMaxPhi->SetAt(maxphi,0);
263   
264 }
265
266 //__________________________________________________________________
267 void AliFiducialCut::SetSimpleEMCALFiducialCut(const Float_t eta, 
268                                                const Float_t minphi, 
269                                                const Float_t maxphi)
270 {
271   //Method to set simple acceptance cut to EMCAL
272   
273   fEMCALFidCutMinEta->Set(1);
274   fEMCALFidCutMaxEta->Set(1);
275   fEMCALFidCutMinPhi->Set(1);
276   fEMCALFidCutMaxPhi->Set(1);
277   
278   fEMCALFidCutMinEta->SetAt(-eta,0);
279   fEMCALFidCutMaxEta->SetAt( eta,0);
280   fEMCALFidCutMinPhi->SetAt(minphi,0);
281   fEMCALFidCutMaxPhi->SetAt(maxphi,0);
282   
283   
284 }
285
286 //_________________________________________________________________
287 void AliFiducialCut::SetSimplePHOSFiducialCut(const Float_t eta, 
288                                               const Float_t minphi, 
289                                               const Float_t maxphi)
290 {
291   //Method to set simple acceptance cut to PHOS
292   
293   fPHOSFidCutMinEta->Set(1);
294   fPHOSFidCutMaxEta->Set(1);
295   fPHOSFidCutMinPhi->Set(1);
296   fPHOSFidCutMaxPhi->Set(1);
297   
298   fPHOSFidCutMinEta->SetAt(-eta,0);
299   fPHOSFidCutMaxEta->SetAt(eta,0);
300   fPHOSFidCutMinPhi->SetAt(minphi,0);
301   fPHOSFidCutMaxPhi->SetAt(maxphi,0);
302   
303 }