]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALPIDv1.cxx
Added the position in global coordinate system of the particle in PHOS
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALPIDv1.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 /* $Id$ */
17
18 //_________________________________________________________________________
19 // Implementation version v1 of the EMCAL particle identifier 
20
21 //*-- Author: Yves Schutz (SUBATECH) 
22 //
23 // --- ROOT system ---
24 //#include "TROOT.h"
25 #include "TTree.h"
26 #include "TBenchmark.h"
27 #include "TSystem.h"
28   
29   // --- Standard library ---
30   
31   // --- AliRoot header files ---
32 #include "AliGenerator.h"
33 #include "AliEMCALPIDv1.h"
34 #include "AliEMCALRecParticle.h"
35 #include "AliEMCALGetter.h"
36   
37   ClassImp( AliEMCALPIDv1) 
38
39 //____________________________________________________________________________
40 AliEMCALPIDv1::AliEMCALPIDv1():AliEMCALPID()
41
42   // default ctor
43   
44   InitParameters() ; 
45   fDefaultInit = kTRUE ; 
46   
47 }
48
49 //____________________________________________________________________________
50 AliEMCALPIDv1::AliEMCALPIDv1(const AliEMCALPIDv1 & pid ):AliEMCALPID(pid)
51
52   // ctor
53   InitParameters() ; 
54   Init() ;
55   
56 }
57
58 //____________________________________________________________________________
59 AliEMCALPIDv1::AliEMCALPIDv1(const TString alirunFileName, const TString eventFolderName):AliEMCALPID(alirunFileName, eventFolderName)
60
61   //ctor with the indication on where to look for the track segments
62  
63   InitParameters() ; 
64   Init() ;
65   fDefaultInit = kFALSE ; 
66
67 }
68
69 //____________________________________________________________________________
70 AliEMCALPIDv1::~AliEMCALPIDv1()
71
72   // dtor
73 }
74
75 //____________________________________________________________________________
76 const TString AliEMCALPIDv1::BranchName() const 
77 {  
78
79   return GetName() ;
80 }
81  
82 //____________________________________________________________________________
83 Float_t  AliEMCALPIDv1::GetCalibratedEnergy(Float_t e) const
84 {
85 //      It calibrates Energy depending on the recpoint energy.
86 //      The energy of the reconstructed cluster is corrected with 
87 //      the formula A + B* E  + C* E^2, whose parameters where obtained 
88 //      through the study of the reconstructed energy distribution of 
89 //      monoenergetic photons.
90  
91   //Float_t p[]={0.,0.,0.};
92   //for (Int_t i=0; i<3; i++) p[i] = GetParameterCalibration(i);
93   Float_t enerec = e ; // p[0] +  p[1]*e + p[2]*e*e;
94   return enerec ;
95
96 }
97 //____________________________________________________________________________
98 TVector3 AliEMCALPIDv1::GetMomentumDirection(AliEMCALTowerRecPoint * emc)const 
99
100   // Calculates the momentum direction:
101   // direction is given by IP and this RecPoint
102   
103
104   TVector3 dir(0,0,0) ; 
105   TVector3 emcglobalpos ;
106   TMatrix  dummy ;
107   
108   emc->GetGlobalPosition(emcglobalpos, dummy) ;
109   
110
111   dir = emcglobalpos ;  
112   dir.SetZ( -dir.Z() ) ;   // why ?  
113   // dir.SetMag(1.) ; Removed to avoid warings !!!!!!!!!!!!!! TO BE REVISED
114
115   //account correction to the position of IP
116   Float_t xo,yo,zo ; //Coordinates of the origin
117   gAlice->Generator()->GetOrigin(xo,yo,zo) ;
118   TVector3 origin(xo,yo,zo);
119   dir = dir - origin ;
120
121   return dir ;  
122 }
123
124 //____________________________________________________________________________
125 void AliEMCALPIDv1::Init()
126 {
127   // Make all memory allocations that are not possible in default constructor
128   // Add the PID task to the list of EMCAL tasks
129
130
131   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName.Data()) ; 
132
133   if ( !gime->PID() ) 
134     gime->PostPID(this) ;
135 }
136
137 //____________________________________________________________________________
138 void AliEMCALPIDv1::InitParameters()
139 {
140   // Initialize the parameters
141   fRecParticlesInRun = 0 ; 
142   fNEvent            = 0 ;            
143   fRecParticlesInRun = 0 ;
144 }
145
146 //____________________________________________________________________________
147
148 void  AliEMCALPIDv1::Exec(Option_t * option) 
149 {
150   //Steering method
151
152   if(strstr(option,"tim"))
153     gBenchmark->Start("EMCALPID");
154   
155   if(strstr(option,"print")) {
156     Print("") ; 
157     return ; 
158   }
159   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle()) ; 
160
161   if (fLastEvent == -1) 
162     fLastEvent = gime->MaxEvent() - 1 ;
163   else 
164     fLastEvent = TMath::Min(fLastEvent,gime->MaxEvent());
165   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
166
167   Int_t ievent ;
168
169
170   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
171     gime->Event(ievent,"R") ;
172     MakeRecParticles() ;
173     WriteRecParticles();
174     if(strstr(option,"deb"))
175       PrintRecParticles(option) ;
176     //increment the total number of rec particles per run 
177     fRecParticlesInRun += gime->RecParticles()->GetEntriesFast() ; 
178     
179   }
180   if(strstr(option,"tim")){
181     gBenchmark->Stop("EMCALPID");
182     printf("Exec: took %f seconds for PID %f seconds per event", 
183          gBenchmark->GetCpuTime("EMCALPID"),  
184          gBenchmark->GetCpuTime("EMCALPID")/nEvents) ;
185   } 
186
187   Unload();
188 }
189
190 //____________________________________________________________________________
191 void  AliEMCALPIDv1::MakeRecParticles(){
192
193   // Makes a RecParticle out of a TrackSegment
194   
195   AliEMCALGetter * gime = AliEMCALGetter::Instance() ; 
196   TObjArray * aECARecPoints = gime->ECARecPoints() ; 
197   if ( !aECARecPoints ) {
198     Fatal("MakeRecParticles", "RecPoints or TrackSegments not found !") ;  
199   }
200   TClonesArray * recParticles  = gime->RecParticles() ; 
201   recParticles->Clear();
202
203   TIter next(aECARecPoints) ; 
204   AliEMCALTowerRecPoint * eca ; 
205   Int_t index = 0 ; 
206   AliEMCALRecParticle * rp ; 
207   while ( (eca = (AliEMCALTowerRecPoint *)next()) ) {
208     
209     new( (*recParticles)[index] ) AliEMCALRecParticle() ;
210     rp = (AliEMCALRecParticle *)recParticles->At(index) ; 
211     rp->SetRecPoint(index) ;
212     rp->SetIndexInList(index) ;
213         
214     // Now set type (reconstructed) of the particle
215
216     // Choose the cluster energy range
217     
218     Float_t  lambda[2] ;
219     eca->GetElipsAxis(lambda) ;
220     
221     if((lambda[0]>0.01) && (lambda[1]>0.01)){
222       // Looking PCA. Define and calculate the data (X),
223       // introduce in the function X2P that gives the components (P).  
224
225       Float_t  spher = 0. ;
226       Float_t  emaxdtotal = 0. ; 
227       
228       if((lambda[0]+lambda[1])!=0) 
229         spher=fabs(lambda[0]-lambda[1])/(lambda[0]+lambda[1]); 
230       
231       emaxdtotal=eca->GetMaximalEnergy()/eca->GetEnergy(); 
232     }
233     
234     //    Float_t time = eca->GetTime() ;
235       
236     //Set momentum, energy and other parameters 
237     Float_t  encal = GetCalibratedEnergy(eca->GetEnergy());
238     TVector3 dir   = GetMomentumDirection(eca) ; 
239     // dir.SetMag(encal) ;Removed to avoid warings !!!!!!!!!!!!!! TO BE REVISED
240     rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),encal) ;
241     rp->SetCalcMass(0);
242     rp->Name(); //If photon sets the particle pdg name to gamma
243     rp->SetProductionVertex(0,0,0,0);
244     rp->SetFirstMother(-1);
245     rp->SetLastMother(-1);
246     rp->SetFirstDaughter(-1);
247     rp->SetLastDaughter(-1);
248     rp->SetPolarisation(0,0,0);
249     //Set the position in global coordinate system from the RecPoint
250     //AliEMCALGeometry * geom = gime->EMCALGeometry() ; 
251     //AliEMCALTowerRecPoint  * erp = gime->ECARecPoint(rp->GetEMCALRPIndex()) ; 
252     TVector3 pos ; 
253     //geom->GetGlobal(erp, pos) ; !!!!!!!!!! to check 
254     rp->SetPos(pos);
255
256
257     index++ ; 
258   }
259   
260 }
261
262 //____________________________________________________________________________
263 void  AliEMCALPIDv1:: Print(Option_t * /*option*/) const
264 {
265   // Print the parameters used for the particle type identification
266
267     printf("Print: =============== AliEMCALPID1 ================") ;
268     printf("Making PID\n");
269     printf("    Pricipal analysis file from 0.5 to 100 %s\n", fFileName.Data() ) ; 
270     printf("    Name of parameters file     %s\n", fFileNamePar.Data() )  ;
271 }
272
273 //____________________________________________________________________________
274 void  AliEMCALPIDv1::Print() const
275 {
276   // Print the parameters used for the particle type identification
277
278     Info("Print", "=============== AliEMCALPIDv1 ================") ;
279 }
280
281 //____________________________________________________________________________
282 void AliEMCALPIDv1::PrintRecParticles(Option_t * option)
283 {
284   // Print table of reconstructed particles
285
286   AliEMCALGetter *gime = AliEMCALGetter::Instance() ; 
287
288   TClonesArray * recParticles = gime->RecParticles() ; 
289
290   printf("\nevent %i", gAlice->GetEvNumber()); 
291   printf("       found %i", recParticles->GetEntriesFast()); 
292   printf(" RecParticles\n"); 
293
294   if(strstr(option,"all")) {  // printing found TS
295     printf("\n  PARTICLE         Index    \n"); 
296     
297     Int_t index ;
298     for (index = 0 ; index < recParticles->GetEntries() ; index++) {
299       AliEMCALRecParticle * rp = (AliEMCALRecParticle * ) recParticles->At(index) ;       
300       printf("\n");
301       printf(rp->Name().Data());  
302       printf(" %i", rp->GetIndexInList());  
303       printf(" %i", rp->GetType());
304     }
305   }
306 }
307
308 //____________________________________________________________________________
309 void AliEMCALPIDv1::Unload() 
310 {
311   // Unloads RecPoints, TrackSegments and RecParticles from the folder 
312   AliEMCALGetter * gime = AliEMCALGetter::Instance() ;  
313   gime->EmcalLoader()->UnloadRecPoints() ;
314   gime->EmcalLoader()->UnloadTracks() ;
315   gime->EmcalLoader()->UnloadRecParticles() ;
316 }
317
318 //____________________________________________________________________________
319 void  AliEMCALPIDv1::WriteRecParticles()
320 {
321   // Write RecParticles array to file
322   AliEMCALGetter *gime = AliEMCALGetter::Instance() ; 
323
324   TClonesArray * recParticles = gime->RecParticles() ; 
325   recParticles->Expand(recParticles->GetEntriesFast() ) ;
326   TTree * treeP =  gime->TreeP() ;
327
328   
329   
330   //First rp
331   Int_t bufferSize = 32000 ;    
332   TBranch * rpBranch = treeP->Branch("EMCALRP",&recParticles,bufferSize);
333   rpBranch->SetTitle(BranchName());
334
335   rpBranch->Fill() ;
336  
337   gime->WriteRecParticles("OVERWRITE");
338   gime->WritePID("OVERWRITE");
339
340 }
341