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