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