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