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