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