]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSv4.cxx
Timing added to reconstruction
[u/mrichter/AliRoot.git] / PHOS / AliPHOSv4.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 of the PHOS manager class for fast simulations     
20 // Tracks particles until the reach a grossly designed PHOS module
21 // Modify the particles property (momentum, energy, type) according to
22 //  the PHOS response function. The result is called a virtual reconstructed
23 //  particle.                
24 //
25 //*-- Author: Yves Schutz (SUBATECH)
26
27 // --- ROOT system ---
28
29 #include "TBRIK.h"
30 #include "TNode.h"
31 #include "TParticle.h"
32 #include "TTree.h"
33 #include "TGeometry.h"
34 #include "TFile.h"
35
36 // --- Standard library ---
37
38 #include <stdio.h>
39
40 // --- AliRoot header files ---
41 #include "AliPHOSFastRecParticle.h"
42 #include "AliPHOSGeometry.h"
43 #include "AliPHOSv4.h"
44 #include "AliPHOSGetter.h"
45 #include "AliRun.h"
46 #include "AliConst.h"
47 #include "AliMC.h"
48
49 ClassImp(AliPHOSv4)
50
51 //____________________________________________________________________________
52 AliPHOSv4::AliPHOSv4() : AliPHOS()
53 {
54   // default ctor : initialize data member
55    fBigBoxX = 0. ;                      
56    fBigBoxY = 0. ;                      
57    fBigBoxZ = 0. ;                       
58    fFastRecParticles = 0 ;        
59    fNRecParticles = 0 ;                
60    fRan = 0 ;                            
61    fResPara1 = 0. ;                       
62    fResPara2 = 0. ;                        
63    fResPara3 = 0. ;                      
64    fPosParaA0 = 0. ;                      
65    fPosParaA1 = 0. ;
66    fPosParaB0 = 0. ;     
67    fPosParaB1 = 0. ;    
68    fPosParaB2 = 0. ;    
69
70 }
71
72 //____________________________________________________________________________
73 AliPHOSv4::AliPHOSv4(const char *name, const char *title):
74   AliPHOS(name,title)
75 {
76   // ctor
77
78   
79   // create the geometry parameters object  
80   // and post it to a folder (Post retrieves the correct geometry)
81   AliPHOSGetter::GetInstance(gDirectory->GetName(), 0)->PostGeometry() ;    
82     
83     SetBigBox(0, GetGeometry()->GetOuterBoxSize(0) ) ;
84     SetBigBox(1, GetGeometry()->GetOuterBoxSize(1) + GetGeometry()->GetCPVBoxSize(1) ) ; 
85     SetBigBox(2, GetGeometry()->GetOuterBoxSize(0) ); 
86     
87     fNRecParticles = 0 ; 
88     fFastRecParticles = new AliPHOSFastRecParticle::FastRecParticlesList("AliPHOSFastRecParticle", 100) ;
89
90     fResPara1 = 0.030 ;    // GeV
91     fResPara2 = 0.00003 ; 
92     fResPara3 = 0.00001 ; 
93     
94     fPosParaA0 = 2.87 ;    // mm
95     fPosParaA1 = -0.0975 ;  
96     fPosParaB0 = 0.257 ;   
97     fPosParaB1 = 0.137 ; 
98     fPosParaB2 = 0.00619 ; 
99 }
100
101 //____________________________________________________________________________
102 AliPHOSv4::~AliPHOSv4()
103 {
104   // dtor
105  
106   fFastRecParticles->Delete() ; 
107   delete fFastRecParticles ;
108   fFastRecParticles = 0 ; 
109
110 }
111
112 //____________________________________________________________________________
113 void AliPHOSv4::AddRecParticle(const AliPHOSFastRecParticle & rp)
114 {  
115   // Add a virtually reconstructed particle to the list 
116
117   new( (*fFastRecParticles)[fNRecParticles] ) AliPHOSFastRecParticle(rp) ;
118   fNRecParticles++ ; 
119 }
120
121 //____________________________________________________________________________
122 void AliPHOSv4::BuildGeometry()
123 {
124   // Build the PHOS geometry for the ROOT display
125    //BEGIN_HTML
126   /*
127     <H2>
128      PHOS FAST in ALICE displayed by root
129     </H2>
130     <H4> All Views </H4>
131     <P>
132     <CENTER>
133     <IMG Align=BOTTOM ALT="Fast All Views" SRC="../images/AliPHOSv4AllViews.gif"> 
134     </CENTER></P>
135     <H4> Front View </H4>
136     <P>
137     <CENTER>
138     <IMG Align=BOTTOM ALT="Fast Front View" SRC="../images/AliPHOSv4FrontView.gif"> 
139     </CENTER></P>
140   */
141   //END_HTML  
142  
143   const Int_t kColorPHOS = kRed ;
144   
145   Double_t const kRADDEG = 180.0 / kPI ;
146   
147   new TBRIK( "BigBox", "PHOS box", "void", GetBigBox(0)/2, 
148              GetBigBox(1)/2, 
149              GetBigBox(2)/2 );
150   
151   // position PHOS into ALICE
152
153   Float_t r = GetGeometry()->GetIPtoCrystalSurface() + GetBigBox(1) / 2.0 ;
154   Int_t number = 988 ; 
155   Float_t pphi =  TMath::ATan( GetBigBox(0)  / ( 2.0 * GetGeometry()->GetIPtoCrystalSurface() ) ) ;
156   pphi *= kRADDEG ;
157   TNode * top = gAlice->GetGeometry()->GetNode("alice") ;
158  
159   char * nodename = new char[20] ;  
160   char * rotname  = new char[20] ; 
161
162   for( Int_t i = 1; i <= GetGeometry()->GetNModules(); i++ ) { 
163    Float_t angle = pphi * 2 * ( i - GetGeometry()->GetNModules() / 2.0 - 0.5 ) ;
164    sprintf(rotname, "%s%d", "rot", number++) ;
165    new TRotMatrix(rotname, rotname, 90, angle, 90, 90 + angle, 0, 0);
166    top->cd();
167    sprintf(nodename,"%s%d", "Module", i) ;    
168    Float_t x =  r * TMath::Sin( angle / kRADDEG ) ;
169    Float_t y = -r * TMath::Cos( angle / kRADDEG ) ;
170    TNode * bigboxnode = new TNode(nodename, nodename, "BigBox", x, y, 0, rotname ) ;
171    bigboxnode->SetLineColor(kColorPHOS) ;
172    fNodes->Add(bigboxnode) ;
173   }
174   delete[] nodename ; 
175   delete[] rotname ; 
176 }
177
178 //____________________________________________________________________________
179 void AliPHOSv4::CreateGeometry()
180 {
181   // Create the geometry for GEANT
182   
183   AliPHOSv4 *phostmp = (AliPHOSv4*)gAlice->GetModule("PHOS") ;
184   
185   if ( phostmp == NULL ) {
186     
187     fprintf(stderr, "PHOS detector not found!\n") ;
188     return ;
189     
190   }
191
192   // Get pointer to the array containing media indeces
193   Int_t *idtmed = fIdtmed->GetArray() - 699 ;
194   
195   Float_t bigbox[3] ; 
196   bigbox[0] =   GetBigBox(0) / 2.0 ;
197   bigbox[1] =   GetBigBox(1) / 2.0 ;
198   bigbox[2] =   GetBigBox(2) / 2.0 ;
199   
200   gMC->Gsvolu("PHOS", "BOX ", idtmed[798], bigbox, 3) ;
201   
202   // --- Position  PHOS mdules in ALICE setup ---
203   
204   Int_t idrotm[99] ;
205   Double_t const kRADDEG = 180.0 / kPI ;
206   
207   for( Int_t i = 1; i <= GetGeometry()->GetNModules(); i++ ) {
208     
209     Float_t angle = GetGeometry()->GetPHOSAngle(i) ;
210     AliMatrix(idrotm[i-1], 90.0, angle, 90.0, 90.0+angle, 0.0, 0.0) ;
211  
212     Float_t r = GetGeometry()->GetIPtoCrystalSurface() + GetBigBox(1) / 2.0 ;
213
214     Float_t xP1 = r * TMath::Sin( angle / kRADDEG ) ;
215     Float_t yP1 = -r * TMath::Cos( angle / kRADDEG ) ;
216     gMC->Gspos("PHOS", i, "ALIC", xP1, yP1, 0.0, idrotm[i-1], "ONLY") ;
217  
218   } // for GetNModules
219
220 }
221
222
223 //____________________________________________________________________________
224 void AliPHOSv4::Init(void)
225 {
226   // Prints out an information message
227   
228   Int_t i;
229
230   printf("\n");
231   for(i=0;i<35;i++) printf("*");
232   printf(" FAST PHOS_INIT ");
233   for(i=0;i<35;i++) printf("*");
234   printf("\n");
235
236   // Here the PHOS initialisation code (if any!)
237
238   for(i=0;i<80;i++) printf("*");
239   printf("\n");
240   
241 }
242
243 //___________________________________________________________________________
244 Float_t AliPHOSv4::GetBigBox(Int_t index) const
245 {
246   // Get the X, Y or Z dimension of the box describing a PHOS module
247   
248   Float_t rv = 0 ; 
249
250   switch (index) {
251   case 0:
252     rv = fBigBoxX ; 
253     break ; 
254   case 1:
255      rv = fBigBoxY ; 
256     break ; 
257   case 2:
258      rv = fBigBoxZ ; 
259     break ; 
260  }
261   return rv ; 
262 }
263
264 //___________________________________________________________________________
265 void AliPHOSv4::MakeBranch(Option_t* opt, const char *file)
266 {  
267   // Create new branch in the current reconstructed Root Tree
268  
269   AliDetector::MakeBranch(opt,file) ;
270   
271   char branchname[10];
272   sprintf(branchname,"%s",GetName());
273   const char *cd = strstr(opt,"R");
274   
275   if (fFastRecParticles && gAlice->TreeR() && cd) {
276     MakeBranchInTree(gAlice->TreeR(), 
277                      branchname, &fFastRecParticles, fBufferSize, file);
278   }
279 }
280
281 //____________________________________________________________________________
282 Double_t AliPHOSv4::MakeEnergy(const Double_t energy)
283 {  
284   // Smears the energy according to the energy dependent energy resolution.
285   // A gaussian distribution is assumed
286
287   Double_t sigma  = SigmaE(energy) ; 
288   return  fRan.Gaus(energy, sigma) ;   
289 }
290
291 //____________________________________________________________________________
292 TVector3 AliPHOSv4::MakePosition(const Double_t energy, const TVector3 pos, const Double_t theta, const Double_t phi)
293 {
294   // Smears the impact position according to the energy dependent position resolution
295   // A gaussian position distribution is assumed
296
297   TVector3 newpos ;
298   Double_t sigma = SigmaP( energy, theta*180./TMath::Pi() ) ;
299   Double_t x = fRan.Gaus( pos.X(), sigma ) ;
300   sigma = SigmaP( energy, phi*180./TMath::Pi() ) ;
301   Double_t z = fRan.Gaus( pos.Z(), sigma ) ;
302   Double_t y = pos.Y() ; 
303   
304   newpos.SetX(x) ; 
305   newpos.SetY(y) ; 
306   newpos.SetZ(z) ; 
307               
308   return newpos ; 
309 }
310
311 //____________________________________________________________________________
312 void AliPHOSv4::MakeRecParticle(const Int_t modid, const TVector3 pos, AliPHOSFastRecParticle & rp)
313 {
314   // Modify the primary particle properties according
315   //  1. the response function of PHOS
316   //  2. the performance of the EMC+PPSD setup
317   
318   Int_t type = MakeType( rp ) ;
319   rp.SetType(type) ;
320
321   
322   // get the detected energy
323
324   TLorentzVector momentum ;  
325   rp.Momentum(momentum) ; 
326   Double_t kineticenergy = TMath::Sqrt( TMath::Power(momentum.E(), 2) - TMath::Power(rp.GetMass(), 2) ) ; 
327   Double_t modifiedkineticenergy = MakeEnergy(kineticenergy ) ;
328   Double_t modifiedenergy = TMath::Sqrt( TMath::Power(modifiedkineticenergy, 2)  
329                                          + TMath::Power( rp.GetMass(), 2) ) ;
330  
331   // get the angle of incidence 
332   
333   Double_t incidencetheta = 90. * TMath::Pi() /180 - rp.Theta() ; 
334   Double_t incidencephi   = ( 270 + GetGeometry()->GetPHOSAngle(modid) ) * TMath::Pi() / 180. - rp.Phi() ;   
335
336   // get the detected direction
337   
338   TVector3 modifiedposition = MakePosition(kineticenergy, pos, incidencetheta, incidencephi) ; 
339   modifiedposition *= modifiedkineticenergy / modifiedposition.Mag() ; 
340
341   // Set the modified 4-momentum of the reconstructed particle
342
343   rp.SetMomentum(modifiedposition.X(), modifiedposition.Y(), modifiedposition.Z(), modifiedenergy) ; 
344
345  }
346
347 //____________________________________________________________________________
348 Int_t AliPHOSv4::MakeType(AliPHOSFastRecParticle & rp )
349 {
350   // Generate a particle type using the performance of the EMC+PPSD setup
351
352   Int_t rv =   AliPHOSFastRecParticle::kUNDEFINED ;
353   Int_t charge = (Int_t)rp.GetPDG()->Charge() ;
354   Int_t test ; 
355   Float_t ran ; 
356   if ( charge != 0 && ( TMath::Abs(rp.GetPdgCode()) != 11 ) ) 
357     test = - 1 ;
358   else
359     test = rp.GetPdgCode() ; 
360
361   switch (test) { 
362
363   case 22:    // it's a photon
364     ran = fRan.Rndm() ; 
365     if ( ran <= 0.5 )  // 50 % 
366       rv =  AliPHOSFastRecParticle::kGAMMA ; 
367     else {
368       ran = fRan.Rndm() ; 
369       if( ran <= 0.9498 )
370         rv =  AliPHOSFastRecParticle::kNEUTRALEM ; 
371       else
372         rv =  AliPHOSFastRecParticle::kNEUTRALHA ; 
373     }     
374     break ; 
375
376   case 2112:  // it's a neutron
377     ran = fRan.Rndm() ; 
378     if ( ran <= 0.9998 )
379       rv =  AliPHOSFastRecParticle::kNEUTRALHA ; 
380     else 
381       rv = AliPHOSFastRecParticle::kNEUTRALEM ; 
382     break ; 
383
384   case -2112: // it's a anti-neutron
385     ran = fRan.Rndm() ; 
386     if ( ran <= 0.9984 )
387       rv =  AliPHOSFastRecParticle::kNEUTRALHA ; 
388     else 
389       rv =  AliPHOSFastRecParticle::kNEUTRALEM ; 
390     break ; 
391     
392   case 11:    // it's a electron
393     ran = fRan.Rndm() ; 
394     if ( ran <= 0.9996 )
395       rv =  AliPHOSFastRecParticle::kELECTRON ; 
396     else 
397       rv =  AliPHOSFastRecParticle::kCHARGEDHA ; 
398     break; 
399
400   case -11:   // it's a positon
401     ran = fRan.Rndm() ; 
402     if ( ran <= 0.9996 )
403       rv =  AliPHOSFastRecParticle::kELECTRON ; 
404     else 
405       rv =  AliPHOSFastRecParticle::kCHARGEDHA ; 
406     break; 
407
408   case -1:    // it's a charged
409     ran = fRan.Rndm() ; 
410     if ( ran <= 0.9996 )
411       rv =  AliPHOSFastRecParticle::kCHARGEDHA ; 
412     else 
413       rv =  AliPHOSFastRecParticle::kGAMMA ; 
414
415     break ; 
416   }
417     
418   
419   return rv ;
420 }
421
422 //___________________________________________________________________________
423 void AliPHOSv4::ResetPoints()
424 {
425   // This overloads the method in AliDetector
426   
427   ResetFastRecParticles() ; 
428 }
429
430 //___________________________________________________________________________
431 void AliPHOSv4::ResetFastRecParticles()
432 {
433   // Resets the list of virtual reconstructed particles
434  
435   if (fFastRecParticles) 
436     fFastRecParticles->Clear() ;
437   fNRecParticles = 0 ; 
438 }
439
440 //___________________________________________________________________________
441 void AliPHOSv4::SetBigBox(Int_t index, Float_t value)
442 {
443   // Set the size of the Box describing a PHOS module
444   
445   switch (index) {
446   case 0:
447     fBigBoxX = value ; 
448     break ; 
449   case 1:
450     fBigBoxY = value ; 
451     break ; 
452   case 2:
453     fBigBoxZ = value ; 
454     break ; 
455  }
456
457 }
458
459 //____________________________________________________________________________
460 Double_t AliPHOSv4::SigmaE(Double_t energy)
461 {
462   // Calculates the energy dependent energy resolution
463   
464   Double_t rv = -1 ; 
465   
466   rv = TMath::Sqrt( TMath::Power(fResPara1/energy, 2) 
467                + TMath::Power(fResPara2/TMath::Sqrt(energy), 2) 
468                + TMath::Power(fResPara3, 2) ) ;  
469
470   return rv * energy ; 
471 }
472
473 //____________________________________________________________________________
474 Double_t AliPHOSv4::SigmaP(Double_t energy, Double_t incidence)
475 {
476   // Calculates the energy dependent position resolution 
477
478   Double_t paraA = fPosParaA0 + fPosParaA1 * incidence ; 
479   Double_t paraB = fPosParaB0 + fPosParaB1 * incidence + fPosParaB2 * incidence * incidence ; 
480
481   return ( paraA / TMath::Sqrt(energy) + paraB ) * 0.1   ; // in cm  
482 }
483
484 //____________________________________________________________________________
485 void AliPHOSv4::StepManager(void)
486 {
487   // Only verifies if the particle reaches PHOS and stops the tracking 
488
489   Int_t primary =  gAlice->GetPrimary( gAlice->CurrentTrack() ); 
490
491   TLorentzVector lv ; 
492   gMC->TrackPosition(lv) ;
493   TVector3 pos = lv.Vect() ; 
494   Int_t modid  ; 
495   gMC->CurrentVolID(modid);
496   
497   Float_t energy = gMC->Etot() ; //Total energy of current track
498
499   //Calculating mass of current particle
500   TDatabasePDG * pdg = TDatabasePDG::Instance() ;
501   TParticlePDG * partPDG = pdg->GetParticle(gMC->TrackPid()) ;
502   Float_t mass = partPDG->Mass() ;
503
504   if(energy > mass){
505     pos.SetMag(TMath::Sqrt(energy*energy-mass*mass)) ;
506     TLorentzVector pTrack(pos, energy) ;  
507
508     TParticle * part = new TParticle(gMC->TrackPid(), 0,-1,-1,-1,-1, pTrack, lv)  ;
509         
510     AliPHOSFastRecParticle rp(*part) ;
511     rp.SetPrimary(primary) ; 
512
513     // Adds the response of PHOS to the particle
514     MakeRecParticle(modid, pos, rp) ;
515     
516     // add the `track' particle to the FastRecParticles list
517   
518     AddRecParticle(rp) ;
519
520     part->Delete() ;
521   }
522   // stop the track as soon PHOS is reached
523   
524   gMC->StopTrack() ; 
525
526 }
527