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