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