]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliMC.cxx
Correct streamer (Haavard)
[u/mrichter/AliRoot.git] / STEER / AliMC.cxx
CommitLineData
5d12ce38 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
d0d4a6b3 18// This class is extracted from the AliRun class
19// and contains all the MC-related functionality
20// The number of dependencies has to be reduced...
21// Author: F.Carminati
22// Federico.Carminati@cern.ch
23
c75cfb51 24#include <RVersion.h>
5d12ce38 25#include <TBrowser.h>
7ca4655f 26#include <TClonesArray.h>
27#include <TGeoManager.h>
5d12ce38 28#include <TStopwatch.h>
29#include <TSystem.h>
30#include <TVirtualMC.h>
b60e0f5e 31
21bf7095 32#include "AliLog.h"
5d12ce38 33#include "AliDetector.h"
34#include "AliGenerator.h"
35#include "AliHeader.h"
36#include "AliLego.h"
37#include "AliMC.h"
38#include "AliMCQA.h"
39#include "AliRun.h"
40#include "AliStack.h"
0561efeb 41#include "AliMagF.h"
5d12ce38 42#include "AliTrackReference.h"
43
44
45ClassImp(AliMC)
46
47//_______________________________________________________________________
48AliMC::AliMC() :
49 fGenerator(0),
50 fEventEnergy(0),
51 fSummEnergy(0),
52 fSum2Energy(0),
53 fTrRmax(1.e10),
54 fTrZmax(1.e10),
90e48c0c 55 fRDecayMax(1.e10),
38ca2ad6 56 fRDecayMin(-1.),
90e48c0c 57 fDecayPdg(0),
5d12ce38 58 fImedia(0),
59 fTransParName("\0"),
60 fMCQA(0),
61 fHitLists(0),
62 fTrackReferences(0)
63
64{
d0d4a6b3 65 //default constructor
0561efeb 66 DecayLimits();
5d12ce38 67}
68
69//_______________________________________________________________________
70AliMC::AliMC(const char *name, const char *title) :
71 TVirtualMCApplication(name, title),
72 fGenerator(0),
73 fEventEnergy(0),
74 fSummEnergy(0),
75 fSum2Energy(0),
76 fTrRmax(1.e10),
77 fTrZmax(1.e10),
90e48c0c 78 fRDecayMax(1.e10),
38ca2ad6 79 fRDecayMin(-1.),
90e48c0c 80 fDecayPdg(0),
5d12ce38 81 fImedia(new TArrayI(1000)),
82 fTransParName("\0"),
83 fMCQA(0),
84 fHitLists(new TList()),
85 fTrackReferences(new TClonesArray("AliTrackReference", 100))
86{
d0d4a6b3 87 //constructor
5d12ce38 88 // Set transport parameters
89 SetTransPar();
0561efeb 90 DecayLimits();
5d12ce38 91 // Prepare the tracking medium lists
92 for(Int_t i=0;i<1000;i++) (*fImedia)[i]=-99;
5d12ce38 93}
94
95//_______________________________________________________________________
96AliMC::AliMC(const AliMC &mc) :
97 TVirtualMCApplication(mc),
98 fGenerator(0),
99 fEventEnergy(0),
100 fSummEnergy(0),
101 fSum2Energy(0),
102 fTrRmax(1.e10),
103 fTrZmax(1.e10),
90e48c0c 104 fRDecayMax(1.e10),
38ca2ad6 105 fRDecayMin(-1.),
90e48c0c 106 fDecayPdg(0),
5d12ce38 107 fImedia(0),
108 fTransParName("\0"),
109 fMCQA(0),
110 fHitLists(0),
111 fTrackReferences(0)
112{
113 //
90e48c0c 114 // Copy constructor for AliMC
5d12ce38 115 //
116 mc.Copy(*this);
117}
118
119//_______________________________________________________________________
120AliMC::~AliMC()
121{
d0d4a6b3 122 //destructor
5d12ce38 123 delete fGenerator;
124 delete fImedia;
125 delete fMCQA;
126 delete fHitLists;
127 // Delete track references
128 if (fTrackReferences) {
129 fTrackReferences->Delete();
130 delete fTrackReferences;
131 fTrackReferences = 0;
132 }
133
134}
135
136//_______________________________________________________________________
6c4904c2 137void AliMC::Copy(TObject &) const
5d12ce38 138{
d0d4a6b3 139 //dummy Copy function
21bf7095 140 AliFatal("Not implemented!");
5d12ce38 141}
142
143//_______________________________________________________________________
144void AliMC::ConstructGeometry()
145{
146 //
4a9de4af 147 // Either load geometry from file or create it through usual
148 // loop on detectors. In the first case the method
149 // AliModule::CreateMaterials() only builds fIdtmed and is postponed
150 // at InitGeometry().
5d12ce38 151 //
152
4a9de4af 153 if(gAlice->IsRootGeometry()){
154 // Load geometry
155 const char *geomfilename = gAlice->GetGeometryFileName();
156 if(gSystem->ExpandPathName(geomfilename)){
157 AliInfo(Form("Loading geometry from file:\n %40s\n\n",geomfilename));
158 TGeoManager::Import(geomfilename);
159 }else{
160 AliInfo(Form("Geometry file %40s not found!\n",geomfilename));
161 return;
162 }
163 }else{
164 // Create modules, materials, geometry
5d12ce38 165 TStopwatch stw;
166 TIter next(gAlice->Modules());
167 AliModule *detector;
21bf7095 168 AliDebug(1, "Geometry creation:");
5d12ce38 169 while((detector = dynamic_cast<AliModule*>(next()))) {
170 stw.Start();
171 // Initialise detector materials and geometry
172 detector->CreateMaterials();
173 detector->CreateGeometry();
21bf7095 174 AliInfo(Form("%10s R:%.2fs C:%.2fs",
175 detector->GetName(),stw.RealTime(),stw.CpuTime()));
5d12ce38 176 }
4a9de4af 177 }
178
5d12ce38 179}
180
661663fa 181//_______________________________________________________________________
182void AliMC::ConstructOpGeometry()
183{
184 //
185 // Loop all detector modules and call DefineOpticalProperties() method
186 //
187
188 TIter next(gAlice->Modules());
189 AliModule *detector;
190 AliInfo("Optical properties definition");
191 while((detector = dynamic_cast<AliModule*>(next()))) {
192 // Initialise detector optical properties
193 detector->DefineOpticalProperties();
194 }
195}
196
5d12ce38 197//_______________________________________________________________________
198void AliMC::InitGeometry()
199{
200 //
27cf3cdc 201 // Initialize detectors
5d12ce38 202 //
203
4a9de4af 204 AliInfo("Initialisation:");
205 TStopwatch stw;
206 TIter next(gAlice->Modules());
207 AliModule *detector;
208 while((detector = dynamic_cast<AliModule*>(next()))) {
209 stw.Start();
27cf3cdc 210 // Initialise detector geometry
4a9de4af 211 if(gAlice->IsRootGeometry()) detector->CreateMaterials();
212 detector->Init();
4a9de4af 213 AliInfo(Form("%10s R:%.2fs C:%.2fs",
214 detector->GetName(),stw.RealTime(),stw.CpuTime()));
215 }
4787b401 216}
217
218//_______________________________________________________________________
219void AliMC::SetAllAlignableVolumes()
220{
221 //
222 // Add alignable volumes (TGeoPNEntries) looping on all
223 // active modules
224 //
225
226 AliInfo(Form("Setting entries for all alignable volumes of active detectors"));
227 AliModule *detector;
228 TIter next(gAlice->Modules());
229 while((detector = dynamic_cast<AliModule*>(next()))) {
230 detector->AddAlignableVolumes();
231 }
5d12ce38 232}
233
234//_______________________________________________________________________
2057aecb 235void AliMC::GeneratePrimaries()
5d12ce38 236{
237 //
238 // Generate primary particles and fill them in the stack.
239 //
240
241 Generator()->Generate();
242}
243
244//_______________________________________________________________________
245void AliMC::SetGenerator(AliGenerator *generator)
246{
247 //
248 // Load the event generator
249 //
250 if(!fGenerator) fGenerator = generator;
251}
252
253//_______________________________________________________________________
254void AliMC::ResetGenerator(AliGenerator *generator)
255{
256 //
257 // Load the event generator
258 //
259 if(fGenerator)
260 if(generator)
21bf7095 261 AliWarning(Form("Replacing generator %s with %s",
262 fGenerator->GetName(),generator->GetName()))
5d12ce38 263 else
21bf7095 264 AliWarning(Form("Replacing generator %s with NULL",
265 fGenerator->GetName()));
5d12ce38 266 fGenerator = generator;
267}
268
269//_______________________________________________________________________
270void AliMC::FinishRun()
271{
272 // Clean generator information
21bf7095 273 AliDebug(1, "fGenerator->FinishRun()");
5d12ce38 274 fGenerator->FinishRun();
275
276 //Output energy summary tables
21bf7095 277 AliDebug(1, "EnergySummary()");
278 ToAliDebug(1, EnergySummary());
5d12ce38 279}
280
281//_______________________________________________________________________
282void AliMC::BeginPrimary()
283{
284 //
285 // Called at the beginning of each primary track
286 //
287
288 // Reset Hits info
289 ResetHits();
290 ResetTrackReferences();
291
292}
293
294//_______________________________________________________________________
295void AliMC::PreTrack()
296{
d0d4a6b3 297 // Actions before the track's transport
5d12ce38 298 TObjArray &dets = *gAlice->Modules();
299 AliModule *module;
300
301 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
302 if((module = dynamic_cast<AliModule*>(dets[i])))
303 module->PreTrack();
304
305 fMCQA->PreTrack();
306}
307
308//_______________________________________________________________________
309void AliMC::Stepping()
310{
311 //
312 // Called at every step during transport
313 //
0561efeb 314
38ca2ad6 315 Int_t id = DetFromMate(gMC->CurrentMedium());
5d12ce38 316 if (id < 0) return;
317
0561efeb 318
319 if ( gMC->IsNewTrack() &&
320 gMC->TrackTime() == 0. &&
38ca2ad6 321 fRDecayMin >= 0. &&
0561efeb 322 fRDecayMax > fRDecayMin &&
323 gMC->TrackPid() == fDecayPdg )
324 {
325 FixParticleDecaytime();
326 }
327
328
329
5d12ce38 330 //
331 // --- If lego option, do it and leave
332 if (gAlice->Lego())
333 gAlice->Lego()->StepManager();
334 else {
335 Int_t copy;
336 //Update energy deposition tables
337 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
ab03c084 338 //
339 // write tracke reference for track which is dissapearing - MI
b60e0f5e 340 if (gMC->IsTrackDisappeared()) {
341 if (gMC->Etot()>0.05) AddTrackReference(GetCurrentTrackNumber());
342 }
5d12ce38 343
344 //Call the appropriate stepping routine;
345 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
346 if(det && det->StepManagerIsEnabled()) {
e569033b 347 if(AliLog::GetGlobalDebugLevel()>0) fMCQA->StepManager(id);
5d12ce38 348 det->StepManager();
349 }
350 }
351}
352
353//_______________________________________________________________________
354void AliMC::EnergySummary()
355{
356 //
357 // Print summary of deposited energy
358 //
359
360 Int_t ndep=0;
361 Float_t edtot=0;
362 Float_t ed, ed2;
363 Int_t kn, i, left, j, id;
364 const Float_t kzero=0;
365 Int_t ievent=gAlice->GetRunLoader()->GetHeader()->GetEvent()+1;
366 //
367 // Energy loss information
368 if(ievent) {
369 printf("***************** Energy Loss Information per event (GEV) *****************\n");
370 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
371 ed=fSummEnergy[kn];
372 if(ed>0) {
373 fEventEnergy[ndep]=kn;
374 if(ievent>1) {
375 ed=ed/ievent;
376 ed2=fSum2Energy[kn];
377 ed2=ed2/ievent;
378 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
379 } else
380 ed2=99;
381 fSummEnergy[ndep]=ed;
382 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
383 edtot+=ed;
384 ndep++;
385 }
386 }
387 for(kn=0;kn<(ndep-1)/3+1;kn++) {
388 left=ndep-kn*3;
389 for(i=0;i<(3<left?3:left);i++) {
390 j=kn*3+i;
391 id=Int_t (fEventEnergy[j]+0.1);
392 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
393 }
394 printf("\n");
395 }
396 //
397 // Relative energy loss in different detectors
398 printf("******************** Relative Energy Loss per event ********************\n");
399 printf("Total energy loss per event %10.3f GeV\n",edtot);
400 for(kn=0;kn<(ndep-1)/5+1;kn++) {
401 left=ndep-kn*5;
402 for(i=0;i<(5<left?5:left);i++) {
403 j=kn*5+i;
404 id=Int_t (fEventEnergy[j]+0.1);
405 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
406 }
407 printf("\n");
408 }
409 for(kn=0;kn<75;kn++) printf("*");
410 printf("\n");
411 }
412 //
413 // Reset the TArray's
414 // fEventEnergy.Set(0);
415 // fSummEnergy.Set(0);
416 // fSum2Energy.Set(0);
417}
418
419//_____________________________________________________________________________
420void AliMC::BeginEvent()
421{
422 //
423 // Clean-up previous event
424 // Energy scores
21bf7095 425 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
426 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
427 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
428 AliDebug(1, " BEGINNING EVENT ");
429 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
430 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
431 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
5d12ce38 432
433 AliRunLoader *runloader=gAlice->GetRunLoader();
434
435 /*******************************/
436 /* Clean after eventual */
437 /* previous event */
438 /*******************************/
439
440
441 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
442 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
443 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
21bf7095 444 AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
5d12ce38 445
446 fEventEnergy.Reset();
447 // Clean detector information
448
449 if (runloader->Stack())
450 runloader->Stack()->Reset();//clean stack -> tree is unloaded
451 else
452 runloader->MakeStack();//or make a new one
504b172d 453
454 if(gAlice->Lego() == 0x0)
455 {
21bf7095 456 AliDebug(1, "fRunLoader->MakeTree(K)");
504b172d 457 runloader->MakeTree("K");
458 }
459
21bf7095 460 AliDebug(1, "gMC->SetStack(fRunLoader->Stack())");
5d12ce38 461 gMC->SetStack(gAlice->GetRunLoader()->Stack());//Was in InitMC - but was moved here
462 //because we don't have guarantee that
463 //stack pointer is not going to change from event to event
464 //since it bellobgs to header and is obtained via RunLoader
465 //
466 // Reset all Detectors & kinematics & make/reset trees
467 //
468
469 runloader->GetHeader()->Reset(gAlice->GetRunNumber(),gAlice->GetEvNumber(),
470 gAlice->GetEventNrInRun());
471// fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
472
504b172d 473 if(gAlice->Lego())
474 {
475 gAlice->Lego()->BeginEvent();
476 return;
477 }
478
5d12ce38 479
21bf7095 480 AliDebug(1, "ResetHits()");
5d12ce38 481 ResetHits();
504b172d 482
21bf7095 483 AliDebug(1, "fRunLoader->MakeTree(H)");
5d12ce38 484 runloader->MakeTree("H");
504b172d 485
21bf7095 486 AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
504b172d 487 runloader->MakeTrackRefsContainer();//for insurance
5d12ce38 488
5d12ce38 489
490 //create new branches and SetAdresses
491 TIter next(gAlice->Modules());
492 AliModule *detector;
493 while((detector = (AliModule*)next()))
494 {
21bf7095 495 AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
5d12ce38 496 detector->MakeBranch("H");
21bf7095 497 AliDebug(2, Form("%s->MakeBranchTR()",detector->GetName()));
5d12ce38 498 detector->MakeBranchTR();
21bf7095 499 AliDebug(2, Form("%s->SetTreeAddress()",detector->GetName()));
5d12ce38 500 detector->SetTreeAddress();
501 }
ab03c084 502 // make branch for AliRun track References
503 TTree * treeTR = runloader->TreeTR();
504 if (treeTR){
505 // make branch for central track references
506 if (!fTrackReferences) fTrackReferences = new TClonesArray("AliTrackReference",0);
507 TBranch *branch;
508 branch = treeTR->Branch("AliRun",&fTrackReferences);
509 branch->SetAddress(&fTrackReferences);
510 }
511 //
5d12ce38 512}
513
514//_______________________________________________________________________
515void AliMC::ResetHits()
516{
517 //
518 // Reset all Detectors hits
519 //
520 TIter next(gAlice->Modules());
521 AliModule *detector;
522 while((detector = dynamic_cast<AliModule*>(next()))) {
523 detector->ResetHits();
524 }
525}
526
527//_______________________________________________________________________
528void AliMC::PostTrack()
529{
d0d4a6b3 530 // Posts tracks for each module
4a9de4af 531 TObjArray &dets = *gAlice->Modules();
532 AliModule *module;
533
534 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
535 if((module = dynamic_cast<AliModule*>(dets[i])))
536 module->PostTrack();
5d12ce38 537}
538
539//_______________________________________________________________________
540void AliMC::FinishPrimary()
541{
542 //
543 // Called at the end of each primary track
544 //
545 AliRunLoader *runloader=gAlice->GetRunLoader();
546 // static Int_t count=0;
547 // const Int_t times=10;
548 // This primary is finished, purify stack
c75cfb51 549#if ROOT_VERSION_CODE > 262152
05c7517d 550 if (!(gMC->SecondariesAreOrdered()))
551 runloader->Stack()->ReorderKine();
c75cfb51 552#endif
5d12ce38 553 runloader->Stack()->PurifyKine();
504b172d 554
5d12ce38 555 TIter next(gAlice->Modules());
556 AliModule *detector;
557 while((detector = dynamic_cast<AliModule*>(next()))) {
558 detector->FinishPrimary();
504b172d 559 AliLoader* loader = detector->GetLoader();
560 if(loader)
561 {
562 TTree* treeH = loader->TreeH();
563 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
564 }
5d12ce38 565 }
566
567 // Write out track references if any
568 if (runloader->TreeTR()) runloader->TreeTR()->Fill();
569}
570
571//_______________________________________________________________________
572void AliMC::FinishEvent()
573{
574 //
575 // Called at the end of the event.
576 //
577
578 //
38ca2ad6 579
5d12ce38 580 if(gAlice->Lego()) gAlice->Lego()->FinishEvent();
581
582 TIter next(gAlice->Modules());
583 AliModule *detector;
584 while((detector = dynamic_cast<AliModule*>(next()))) {
585 detector->FinishEvent();
586 }
587
588 //Update the energy deposit tables
589 Int_t i;
590 for(i=0;i<fEventEnergy.GetSize();i++)
591 {
592 fSummEnergy[i]+=fEventEnergy[i];
593 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
594 }
595
596 AliRunLoader *runloader=gAlice->GetRunLoader();
597
598 AliHeader* header = runloader->GetHeader();
599 AliStack* stack = runloader->Stack();
600 if ( (header == 0x0) || (stack == 0x0) )
601 {//check if we got header and stack. If not cry and exit aliroot
21bf7095 602 AliFatal("Can not get the stack or header from LOADER");
5d12ce38 603 return;//never reached
604 }
605 // Update Header information
606 header->SetNprimary(stack->GetNprimary());
607 header->SetNtrack(stack->GetNtrack());
608
609
610 // Write out the kinematics
38ca2ad6 611 if (!gAlice->Lego()) stack->FinishEvent();
5d12ce38 612
613 // Write out the event Header information
614 TTree* treeE = runloader->TreeE();
615 if (treeE)
616 {
617 header->SetStack(stack);
618 treeE->Fill();
619 }
620 else
621 {
21bf7095 622 AliError("Can not get TreeE from RL");
5d12ce38 623 }
624
504b172d 625 if(gAlice->Lego() == 0x0)
626 {
627 runloader->WriteKinematics("OVERWRITE");
628 runloader->WriteTrackRefs("OVERWRITE");
629 runloader->WriteHits("OVERWRITE");
630 }
631
21bf7095 632 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
633 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
634 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
635 AliDebug(1, " FINISHING EVENT ");
636 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
637 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
638 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
5d12ce38 639}
640
641//_______________________________________________________________________
642void AliMC::Field(const Double_t* x, Double_t* b) const
643{
d0d4a6b3 644 // Calculates field "b" at point "x"
5d12ce38 645 gAlice->Field(x,b);
646}
647
648//_______________________________________________________________________
649void AliMC::Init()
650{
d0d4a6b3 651 // MC initialization
5d12ce38 652
653 //=================Create Materials and geometry
654 gMC->Init();
4787b401 655 //Set alignable volumes for the whole geometry
656 SetAllAlignableVolumes();
5d12ce38 657 //Read the cuts for all materials
658 ReadTransPar();
659 //Build the special IMEDIA table
660 MediaTable();
661
662 //Compute cross-sections
663 gMC->BuildPhysics();
664
5d12ce38 665 //Initialise geometry deposition table
666 fEventEnergy.Set(gMC->NofVolumes()+1);
667 fSummEnergy.Set(gMC->NofVolumes()+1);
668 fSum2Energy.Set(gMC->NofVolumes()+1);
669
670 //
671 fMCQA = new AliMCQA(gAlice->GetNdets());
672
673 // Register MC in configuration
674 AliConfig::Instance()->Add(gMC);
675
676}
677
678//_______________________________________________________________________
679void AliMC::MediaTable()
680{
681 //
682 // Built media table to get from the media number to
683 // the detector id
684 //
685
686 Int_t kz, nz, idt, lz, i, k, ind;
687 // Int_t ibeg;
688 TObjArray &dets = *gAlice->Detectors();
689 AliModule *det;
690 Int_t ndets=gAlice->GetNdets();
691 //
692 // For all detectors
693 for (kz=0;kz<ndets;kz++) {
694 // If detector is defined
695 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
696 TArrayI &idtmed = *(det->GetIdtmed());
697 for(nz=0;nz<100;nz++) {
0561efeb 698
5d12ce38 699 // Find max and min material number
700 if((idt=idtmed[nz])) {
701 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
702 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
703 }
704 }
705 if(det->LoMedium() > det->HiMedium()) {
706 det->LoMedium() = 0;
707 det->HiMedium() = 0;
708 } else {
709 if(det->HiMedium() > fImedia->GetSize()) {
21bf7095 710 AliError(Form("Increase fImedia from %d to %d",
711 fImedia->GetSize(),det->HiMedium()));
5d12ce38 712 return;
713 }
714 // Tag all materials in rage as belonging to detector kz
715 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
716 (*fImedia)[lz]=kz;
717 }
718 }
719 }
720 }
721 //
722 // Print summary table
21bf7095 723 AliInfo("Tracking media ranges:");
724 ToAliInfo(
5d12ce38 725 for(i=0;i<(ndets-1)/6+1;i++) {
726 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
727 ind=i*6+k;
728 det=dynamic_cast<AliModule*>(dets[ind]);
729 if(det)
730 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
731 det->HiMedium());
732 else
733 printf(" %6s: %3d -> %3d;","NULL",0,0);
734 }
735 printf("\n");
736 }
21bf7095 737 )
5d12ce38 738}
739
740//_______________________________________________________________________
741void AliMC::ReadTransPar()
742{
743 //
744 // Read filename to set the transport parameters
745 //
746
747
748 const Int_t kncuts=10;
749 const Int_t knflags=11;
750 const Int_t knpars=kncuts+knflags;
751 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
752 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
753 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
754 "MULS","PAIR","PHOT","RAYL"};
755 char line[256];
756 char detName[7];
757 char* filtmp;
758 Float_t cut[kncuts];
759 Int_t flag[knflags];
760 Int_t i, itmed, iret, ktmed, kz;
761 FILE *lun;
762 //
763 // See whether the file is there
764 filtmp=gSystem->ExpandPathName(fTransParName.Data());
765 lun=fopen(filtmp,"r");
766 delete [] filtmp;
767 if(!lun) {
21bf7095 768 AliWarning(Form("File %s does not exist!",fTransParName.Data()));
5d12ce38 769 return;
770 }
771 //
5d12ce38 772 while(1) {
773 // Initialise cuts and flags
774 for(i=0;i<kncuts;i++) cut[i]=-99;
775 for(i=0;i<knflags;i++) flag[i]=-99;
776 itmed=0;
777 for(i=0;i<256;i++) line[i]='\0';
778 // Read up to the end of line excluded
779 iret=fscanf(lun,"%[^\n]",line);
780 if(iret<0) {
781 //End of file
782 fclose(lun);
5d12ce38 783 return;
784 }
785 // Read the end of line
786 fscanf(lun,"%*c");
787 if(!iret) continue;
788 if(line[0]=='*') continue;
789 // Read the numbers
790 iret=sscanf(line,"%s %d %f %f %f %f %f %f %f %f %f %f %d %d %d %d %d %d %d %d %d %d %d",
791 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
792 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
793 &flag[8],&flag[9],&flag[10]);
794 if(!iret) continue;
795 if(iret<0) {
796 //reading error
21bf7095 797 AliWarning(Form("Error reading file %s",fTransParName.Data()));
5d12ce38 798 continue;
799 }
800 // Check that the module exist
801 AliModule *mod = gAlice->GetModule(detName);
802 if(mod) {
803 // Get the array of media numbers
804 TArrayI &idtmed = *mod->GetIdtmed();
805 // Check that the tracking medium code is valid
806 if(0<=itmed && itmed < 100) {
807 ktmed=idtmed[itmed];
808 if(!ktmed) {
21bf7095 809 AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
5d12ce38 810 continue;
811 }
812 // Set energy thresholds
813 for(kz=0;kz<kncuts;kz++) {
814 if(cut[kz]>=0) {
21bf7095 815 AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
816 kpars[kz],cut[kz],itmed,mod->GetName()));
5d12ce38 817 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
818 }
819 }
820 // Set transport mechanisms
821 for(kz=0;kz<knflags;kz++) {
822 if(flag[kz]>=0) {
21bf7095 823 AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
824 kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
5d12ce38 825 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
826 }
827 }
828 } else {
21bf7095 829 AliWarning(Form("Invalid medium code %d",itmed));
5d12ce38 830 continue;
831 }
832 } else {
21bf7095 833 AliDebug(1, Form("%s not present",detName));
5d12ce38 834 continue;
835 }
836 }
837}
838
839//_______________________________________________________________________
840void AliMC::SetTransPar(const char *filename)
841{
842 //
843 // Sets the file name for transport parameters
844 //
845 fTransParName = filename;
846}
847
848//_______________________________________________________________________
0054628d 849void AliMC::Browse(TBrowser *b)
5d12ce38 850{
851 //
852 // Called when the item "Run" is clicked on the left pane
853 // of the Root browser.
854 // It displays the Root Trees and all detectors.
855 //
856 //detectors are in folders anyway
857 b->Add(fMCQA,"AliMCQA");
858}
859
5d12ce38 860//_______________________________________________________________________
861void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
862{
863 //
864 // Add a hit to detector id
865 //
866 TObjArray &dets = *gAlice->Modules();
867 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
868}
869
870//_______________________________________________________________________
871void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
872{
873 //
874 // Add digit to detector id
875 //
876 TObjArray &dets = *gAlice->Modules();
877 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
878}
879
880//_______________________________________________________________________
881Int_t AliMC::GetCurrentTrackNumber() const {
882 //
883 // Returns current track
884 //
885 return gAlice->GetRunLoader()->Stack()->GetCurrentTrackNumber();
886}
887
888//_______________________________________________________________________
889void AliMC::DumpPart (Int_t i) const
890{
891 //
892 // Dumps particle i in the stack
893 //
894 AliRunLoader * runloader = gAlice->GetRunLoader();
895 if (runloader->Stack())
896 runloader->Stack()->DumpPart(i);
897}
898
899//_______________________________________________________________________
900void AliMC::DumpPStack () const
901{
902 //
903 // Dumps the particle stack
904 //
905 AliRunLoader * runloader = gAlice->GetRunLoader();
906 if (runloader->Stack())
907 runloader->Stack()->DumpPStack();
908}
909
910//_______________________________________________________________________
911Int_t AliMC::GetNtrack() const {
912 //
913 // Returns number of tracks in stack
914 //
915 Int_t ntracks = -1;
916 AliRunLoader * runloader = gAlice->GetRunLoader();
917 if (runloader->Stack())
918 ntracks = runloader->Stack()->GetNtrack();
919 return ntracks;
920}
921
922//_______________________________________________________________________
923Int_t AliMC::GetPrimary(Int_t track) const
924{
925 //
926 // return number of primary that has generated track
927 //
928 Int_t nprimary = -999;
929 AliRunLoader * runloader = gAlice->GetRunLoader();
930 if (runloader->Stack())
931 nprimary = runloader->Stack()->GetPrimary(track);
932 return nprimary;
933}
934
935//_______________________________________________________________________
936TParticle* AliMC::Particle(Int_t i) const
937{
d0d4a6b3 938 // Returns the i-th particle from the stack taking into account
939 // the remaping done by PurifyKine
5d12ce38 940 AliRunLoader * runloader = gAlice->GetRunLoader();
941 if (runloader)
942 if (runloader->Stack())
943 return runloader->Stack()->Particle(i);
944 return 0x0;
945}
946
947//_______________________________________________________________________
948TObjArray* AliMC::Particles() const {
949 //
950 // Returns pointer to Particles array
951 //
952 AliRunLoader * runloader = gAlice->GetRunLoader();
953 if (runloader)
954 if (runloader->Stack())
955 return runloader->Stack()->Particles();
956 return 0x0;
957}
958
959//_______________________________________________________________________
960void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
961 Float_t *vpos, Float_t *polar, Float_t tof,
2057aecb 962 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
5d12ce38 963{
964// Delegate to stack
965//
966 AliRunLoader * runloader = gAlice->GetRunLoader();
967 if (runloader)
968 if (runloader->Stack())
969 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
970 mech, ntr, weight, is);
971}
972
973//_______________________________________________________________________
974void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
975 Double_t px, Double_t py, Double_t pz, Double_t e,
976 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
977 Double_t polx, Double_t poly, Double_t polz,
2057aecb 978 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
5d12ce38 979{
980 // Delegate to stack
981 //
982 AliRunLoader * runloader = gAlice->GetRunLoader();
983 if (runloader)
984 if (runloader->Stack())
985 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
986 polx, poly, polz, mech, ntr, weight, is);
987}
988
989//_______________________________________________________________________
2057aecb 990void AliMC::SetHighWaterMark(Int_t nt) const
5d12ce38 991{
992 //
993 // Set high water mark for last track in event
994 AliRunLoader * runloader = gAlice->GetRunLoader();
995 if (runloader)
996 if (runloader->Stack())
997 runloader->Stack()->SetHighWaterMark(nt);
998}
999
1000//_______________________________________________________________________
2057aecb 1001void AliMC::KeepTrack(Int_t track) const
5d12ce38 1002{
1003 //
1004 // Delegate to stack
1005 //
1006 AliRunLoader * runloader = gAlice->GetRunLoader();
1007 if (runloader)
1008 if (runloader->Stack())
1009 runloader->Stack()->KeepTrack(track);
1010}
1011
1012//_______________________________________________________________________
2057aecb 1013void AliMC::FlagTrack(Int_t track) const
5d12ce38 1014{
1015 // Delegate to stack
1016 //
1017 AliRunLoader * runloader = gAlice->GetRunLoader();
1018 if (runloader)
1019 if (runloader->Stack())
1020 runloader->Stack()->FlagTrack(track);
1021}
1022
1023//_______________________________________________________________________
2057aecb 1024void AliMC::SetCurrentTrack(Int_t track) const
5d12ce38 1025{
1026 //
1027 // Set current track number
1028 //
1029 AliRunLoader * runloader = gAlice->GetRunLoader();
1030 if (runloader)
1031 if (runloader->Stack())
1032 runloader->Stack()->SetCurrentTrack(track);
1033}
1034
1035//_______________________________________________________________________
1036void AliMC::AddTrackReference(Int_t label){
1037 //
1038 // add a trackrefernce to the list
1039 if (!fTrackReferences) {
21bf7095 1040 AliError("Container trackrefernce not active");
5d12ce38 1041 return;
1042 }
1043 Int_t nref = fTrackReferences->GetEntriesFast();
1044 TClonesArray &lref = *fTrackReferences;
1045 new(lref[nref]) AliTrackReference(label);
1046}
1047
1048
1049
1050//_______________________________________________________________________
1051void AliMC::ResetTrackReferences()
1052{
1053 //
1054 // Reset all references
1055 //
1056 if (fTrackReferences) fTrackReferences->Clear();
1057
1058 TIter next(gAlice->Modules());
1059 AliModule *detector;
1060 while((detector = dynamic_cast<AliModule*>(next()))) {
1061 detector->ResetTrackReferences();
1062 }
1063}
1064
1065void AliMC::RemapTrackReferencesIDs(Int_t *map)
1066{
1067 //
1068 // Remapping track reference
1069 // Called at finish primary
1070 //
1071 if (!fTrackReferences) return;
4e490a96 1072 Int_t nEntries = fTrackReferences->GetEntries();
1073 for (Int_t i=0; i < nEntries; i++){
1074 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
1075 if (ref) {
1076 Int_t newID = map[ref->GetTrack()];
1077 if (newID>=0) ref->SetTrack(newID);
1078 else {
1079 ref->SetBit(kNotDeleted,kFALSE);
1080 fTrackReferences->RemoveAt(i);
1081 }
1082 } // if ref
5d12ce38 1083 }
1084 fTrackReferences->Compress();
1085}
0561efeb 1086
1087void AliMC::FixParticleDecaytime()
1088{
1089 //
1090 // Fix the particle decay time according to rmin and rmax for decays
1091 //
1092
1093 TLorentzVector p;
1094 gMC->TrackMomentum(p);
1095 Double_t tmin, tmax;
1096 Double_t b;
1097
1098 // Transverse velocity
1099 Double_t vt = p.Pt() / p.E();
1100
1101 if ((b = gAlice->Field()->SolenoidField()) > 0.) { // [kG]
1102
1103 // Radius of helix
1104
1105 Double_t rho = p.Pt() / 0.0003 / b; // [cm]
1106
1107 // Revolution frequency
1108
1109 Double_t omega = vt / rho;
1110
1111 // Maximum and minimum decay time
1112 //
1113 // Check for curlers first
1114 if (fRDecayMax * fRDecayMax / rho / rho / 2. > 1.) return;
1115
1116 //
1117
1118 tmax = TMath::ACos(1. - fRDecayMax * fRDecayMax / rho / rho / 2.) / omega; // [ct]
1119 tmin = TMath::ACos(1. - fRDecayMin * fRDecayMin / rho / rho / 2.) / omega; // [ct]
1120 } else {
1121 tmax = fRDecayMax / vt; // [ct]
1122 tmin = fRDecayMin / vt; // [ct]
1123 }
1124
1125 //
1126 // Dial t using the two limits
1127 Double_t t = tmin + (tmax - tmin) * gRandom->Rndm(); // [ct]
1128 //
1129 //
1130 // Force decay time in transport code
1131 //
fa7ca4da 1132 gMC->ForceDecayTime(t / 2.99792458e10);
0561efeb 1133}