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