1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.54 2001/02/07 10:39:05 hristov
19 Remove default value for argument
21 Revision 1.53 2001/02/06 11:02:26 hristov
22 New SetTrack interface added, added check for unfilled particles in FinishEvent (I.Hrivnacova)
24 Revision 1.52 2001/02/05 16:22:25 buncic
25 Added TreeS to GetEvent().
27 Revision 1.51 2001/02/02 15:16:20 morsch
28 SetHighWaterMark method added to mark last particle in event.
30 Revision 1.50 2001/01/27 10:32:00 hristov
31 Leave the loop when primaries are filled (I.Hrivnacova)
33 Revision 1.49 2001/01/26 19:58:48 hristov
34 Major upgrade of AliRoot code
36 Revision 1.48 2001/01/17 10:50:50 hristov
37 Corrections to destructors
39 Revision 1.47 2000/12/18 10:44:01 morsch
40 Possibility to set field map by passing pointer to objet of type AliMagF via
43 gAlice->SetField(new AliMagFCM("Map2", "$(ALICE_ROOT)/data/field01.dat",2,1.,10.));
45 Revision 1.46 2000/12/14 19:29:27 fca
46 galice.cuts was not read any more
48 Revision 1.45 2000/11/30 07:12:49 alibrary
49 Introducing new Rndm and QA classes
51 Revision 1.44 2000/10/26 13:58:59 morsch
52 Add possibility to choose the lego generator (of type AliGeneratorLego or derived) when running
53 RunLego(). Default is the base class AliGeneratorLego.
55 Revision 1.43 2000/10/09 09:43:17 fca
56 Special remapping of hits for TPC and TRD. End-of-primary action introduced
58 Revision 1.42 2000/10/02 21:28:14 fca
59 Removal of useless dependecies via forward declarations
61 Revision 1.41 2000/07/13 16:19:09 fca
62 Mainly coding conventions + some small bug fixes
64 Revision 1.40 2000/07/12 08:56:25 fca
65 Coding convention correction and warning removal
67 Revision 1.39 2000/07/11 18:24:59 fca
68 Coding convention corrections + few minor bug fixes
70 Revision 1.38 2000/06/20 13:05:45 fca
71 Writing down the TREE headers before job starts
73 Revision 1.37 2000/06/09 20:05:11 morsch
74 Introduce possibility to chose magnetic field version 3: AliMagFDM + field02.dat
76 Revision 1.36 2000/06/08 14:03:58 hristov
77 Only one initializer for a default argument
79 Revision 1.35 2000/06/07 10:13:14 hristov
80 Delete only existent objects.
82 Revision 1.34 2000/05/18 10:45:38 fca
83 Delete Particle Factory properly
85 Revision 1.33 2000/05/16 13:10:40 fca
86 New method IsNewTrack and fix for a problem in Father-Daughter relations
88 Revision 1.32 2000/04/27 10:38:21 fca
89 Correct termination of Lego Run and introduce Lego getter in AliRun
91 Revision 1.31 2000/04/26 10:17:32 fca
92 Changes in Lego for G4 compatibility
94 Revision 1.30 2000/04/18 19:11:40 fca
95 Introduce variable Config.C function signature
97 Revision 1.29 2000/04/07 11:12:34 fca
98 G4 compatibility changes
100 Revision 1.28 2000/04/05 06:51:06 fca
101 Workaround for an HP compiler problem
103 Revision 1.27 2000/03/22 18:08:07 fca
104 Rationalisation of the virtual MC interfaces
106 Revision 1.26 2000/03/22 13:42:26 fca
107 SetGenerator does not replace an existing generator, ResetGenerator does
109 Revision 1.25 2000/02/23 16:25:22 fca
110 AliVMC and AliGeant3 classes introduced
111 ReadEuclid moved from AliRun to AliModule
113 Revision 1.24 2000/01/19 17:17:20 fca
114 Introducing a list of lists of hits -- more hits allowed for detector now
116 Revision 1.23 1999/12/03 11:14:31 fca
117 Fixing previous wrong checking
119 Revision 1.21 1999/11/25 10:40:08 fca
120 Fixing daughters information also in primary tracks
122 Revision 1.20 1999/10/04 18:08:49 fca
123 Adding protection against inconsistent Euclid files
125 Revision 1.19 1999/09/29 07:50:40 fca
126 Introduction of the Copyright and cvs Log
130 ///////////////////////////////////////////////////////////////////////////////
132 // Control class for Alice C++ //
133 // Only one single instance of this class exists. //
134 // The object is created in main program aliroot //
135 // and is pointed by the global gAlice. //
137 // -Supports the list of all Alice Detectors (fModules). //
138 // -Supports the list of particles (fParticles). //
139 // -Supports the Trees. //
140 // -Supports the geometry. //
141 // -Supports the event display. //
144 <img src="picts/AliRunClass.gif">
149 <img src="picts/alirun.gif">
153 ///////////////////////////////////////////////////////////////////////////////
158 #include <iostream.h>
166 #include <TObjectTable.h>
168 #include <TGeometry.h>
170 #include "TBrowser.h"
172 #include "TParticle.h"
174 #include "AliDisplay.h"
177 #include "AliMagFC.h"
178 #include "AliMagFCM.h"
179 #include "AliMagFDM.h"
181 #include "TRandom3.h"
183 #include "AliGenerator.h"
184 #include "AliLegoGenerator.h"
186 #include "AliDetector.h"
190 static AliHeader *gAliHeader;
194 //_____________________________________________________________________________
198 // Default constructor for AliRun
223 fPDGDB = 0; //Particle factory object!
225 fConfigFunction = "\0";
228 fTransParName = "\0";
229 fBaseFileName = ".\0";
231 fParticleMap = new TObjArray(10000);
234 //_____________________________________________________________________________
235 AliRun::AliRun(const char *name, const char *title)
239 // Constructor for the main processor.
240 // Creates the geometry
241 // Creates the list of Detectors.
242 // Creates the list of particles.
259 fConfigFunction = "Config();";
261 // Set random number generator
262 gRandom = fRandom = new TRandom3();
264 if (gSystem->Getenv("CONFIG_SEED")) {
265 gRandom->SetSeed((UInt_t)atoi(gSystem->Getenv("CONFIG_SEED")));
268 gROOT->GetListOfBrowsables()->Add(this,name);
270 // create the support list for the various Detectors
271 fModules = new TObjArray(77);
273 // Create the TNode geometry for the event display
275 BuildSimpleGeometry();
285 // Create the particle stack
286 fParticles = new TClonesArray("TParticle",1000);
290 // Create default mag field
295 // Prepare the tracking medium lists
296 fImedia = new TArrayI(1000);
297 for(i=0;i<1000;i++) (*fImedia)[i]=-99;
300 fPDGDB = TDatabasePDG::Instance(); //Particle factory object!
302 // Create HitLists list
303 fHitLists = new TList();
306 fBaseFileName = ".\0";
308 fParticleMap = new TObjArray(10000);
312 //_____________________________________________________________________________
316 // Default AliRun destructor
336 fParticles->Delete();
344 //_____________________________________________________________________________
345 void AliRun::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
348 // Add a hit to detector id
350 TObjArray &dets = *fModules;
351 if(dets[id]) ((AliModule*) dets[id])->AddHit(track,vol,hits);
354 //_____________________________________________________________________________
355 void AliRun::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
358 // Add digit to detector id
360 TObjArray &dets = *fModules;
361 if(dets[id]) ((AliModule*) dets[id])->AddDigit(tracks,digits);
364 //_____________________________________________________________________________
365 void AliRun::Browse(TBrowser *b)
368 // Called when the item "Run" is clicked on the left pane
369 // of the Root browser.
370 // It displays the Root Trees and all detectors.
372 if (fTreeK) b->Add(fTreeK,fTreeK->GetName());
373 if (fTreeH) b->Add(fTreeH,fTreeH->GetName());
374 if (fTreeD) b->Add(fTreeD,fTreeD->GetName());
375 if (fTreeE) b->Add(fTreeE,fTreeE->GetName());
376 if (fTreeR) b->Add(fTreeR,fTreeR->GetName());
377 if (fTreeS) b->Add(fTreeS,fTreeS->GetName());
379 TIter next(fModules);
381 while((detector = (AliModule*)next())) {
382 b->Add(detector,detector->GetName());
384 b->Add(fMCQA,"AliMCQA");
387 //_____________________________________________________________________________
391 // Initialize Alice geometry
396 //_____________________________________________________________________________
397 void AliRun::BuildSimpleGeometry()
400 // Create a simple TNode geometry used by Root display engine
402 // Initialise geometry
404 fGeometry = new TGeometry("AliceGeom","Galice Geometry for Hits");
405 new TMaterial("void","Vacuum",0,0,0); //Everything is void
406 TBRIK *brik = new TBRIK("S_alice","alice volume","void",2000,2000,3000);
407 brik->SetVisibility(0);
408 new TNode("alice","alice","S_alice");
411 //_____________________________________________________________________________
412 void AliRun::CleanDetectors()
415 // Clean Detectors at the end of event
417 TIter next(fModules);
419 while((detector = (AliModule*)next())) {
420 detector->FinishEvent();
424 //_____________________________________________________________________________
425 void AliRun::CleanParents()
428 // Clean Particles stack.
429 // Set parent/daughter relations
431 TObjArray &particles = *fParticleMap;
434 for(i=0; i<fHgwmk+1; i++) {
435 part = (TParticle *)particles.At(i);
436 if(part) if(!part->TestBit(kDaughtersBit)) {
437 part->SetFirstDaughter(-1);
438 part->SetLastDaughter(-1);
443 //_____________________________________________________________________________
444 Int_t AliRun::DistancetoPrimitive(Int_t, Int_t)
447 // Return the distance from the mouse to the AliRun object
453 //_____________________________________________________________________________
454 void AliRun::DumpPart (Int_t i) const
457 // Dumps particle i in the stack
459 ((TParticle*) (*fParticleMap)[i])->Print();
462 //_____________________________________________________________________________
463 void AliRun::DumpPStack () const
466 // Dumps the particle stack
468 TObjArray &particles = *fParticleMap;
470 "\n\n=======================================================================\n");
471 for (Int_t i=0;i<fNtrack;i++)
473 printf("-> %d ",i); ((TParticle*) particles[i])->Print();
474 printf("--------------------------------------------------------------\n");
477 "\n=======================================================================\n\n");
480 void AliRun::SetField(AliMagF* magField)
482 // Set Magnetic Field Map
487 //_____________________________________________________________________________
488 void AliRun::SetField(Int_t type, Int_t version, Float_t scale,
489 Float_t maxField, char* filename)
492 // Set magnetic field parameters
493 // type Magnetic field transport flag 0=no field, 2=helix, 3=Runge Kutta
494 // version Magnetic field map version (only 1 active now)
495 // scale Scale factor for the magnetic field
496 // maxField Maximum value for the magnetic field
499 // --- Sanity check on mag field flags
500 if(fField) delete fField;
502 fField = new AliMagFC("Map1"," ",type,scale,maxField);
503 } else if(version<=2) {
504 fField = new AliMagFCM("Map2-3",filename,type,scale,maxField);
506 } else if(version==3) {
507 fField = new AliMagFDM("Map4",filename,type,scale,maxField);
510 Warning("SetField","Invalid map %d\n",version);
514 //_____________________________________________________________________________
515 void AliRun::PreTrack()
517 TObjArray &dets = *fModules;
520 for(Int_t i=0; i<=fNdets; i++)
521 if((module = (AliModule*)dets[i]))
527 //_____________________________________________________________________________
528 void AliRun::PostTrack()
530 TObjArray &dets = *fModules;
533 for(Int_t i=0; i<=fNdets; i++)
534 if((module = (AliModule*)dets[i]))
538 //_____________________________________________________________________________
539 void AliRun::FinishPrimary()
542 // Called at the end of each primary track
545 // static Int_t count=0;
546 // const Int_t times=10;
547 // This primary is finished, purify stack
550 TIter next(fModules);
552 while((detector = (AliModule*)next())) {
553 detector->FinishPrimary();
556 // Write out hits if any
557 if (gAlice->TreeH()) {
558 gAlice->TreeH()->Fill();
565 // if(++count%times==1) gObjectTable->Print();
568 //_____________________________________________________________________________
569 void AliRun::FinishEvent()
572 // Called at the end of the event.
576 if(fLego) fLego->FinishEvent();
578 //Update the energy deposit tables
580 for(i=0;i<fEventEnergy.GetSize();i++) {
581 fSummEnergy[i]+=fEventEnergy[i];
582 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
584 fEventEnergy.Reset();
586 // Clean detector information
589 // Write out the kinematics
593 Bool_t allFilled = kFALSE;
595 for(i=0; i<fHgwmk+1; ++i) if((part=fParticleMap->At(i))) {
596 fParticleBuffer = (TParticle*) part;
597 fParticleFileMap[i]= (Int_t) fTreeK->GetEntries();
599 (*fParticleMap)[i]=0;
601 // When all primaries were filled no particle!=0
602 // should be left => to be removed later.
603 if (allFilled) printf("Why != 0 part # %d?\n",i);
606 // // printf("Why = 0 part # %d?\n",i); => We know.
608 // we don't break now in order to be sure there is no
609 // particle !=0 left.
610 // To be removed later and replaced with break.
611 if(!allFilled) allFilled = kTRUE;
615 // Write out the digits
626 // Write out reconstructed clusters
631 // Write out the event Header information
632 if (fTreeE) fTreeE->Fill();
637 // Write Tree headers
638 if (fTreeK) fTreeK->Write(0,TObject::kOverwrite);
639 if (fTreeH) fTreeH->Write(0,TObject::kOverwrite);
640 if (fTreeD) fTreeD->Write(0,TObject::kOverwrite);
641 if (fTreeR) fTreeR->Write(0,TObject::kOverwrite);
642 if (fTreeS) fTreeS->Write(0,TObject::kOverwrite);
647 //_____________________________________________________________________________
648 void AliRun::FinishRun()
651 // Called at the end of the run.
655 if(fLego) fLego->FinishRun();
657 // Clean detector information
658 TIter next(fModules);
660 while((detector = (AliModule*)next())) {
661 detector->FinishRun();
664 //Output energy summary tables
667 TFile *file = fTreeE->GetCurrentFile();
671 fTreeE->Write(0,TObject::kOverwrite);
673 // Write AliRun info and all detectors parameters
676 // Clean tree information
678 delete fTreeK; fTreeK = 0;
681 delete fTreeH; fTreeH = 0;
684 delete fTreeD; fTreeD = 0;
687 delete fTreeR; fTreeR = 0;
690 delete fTreeE; fTreeE = 0;
697 //_____________________________________________________________________________
698 void AliRun::FlagTrack(Int_t track)
701 // Flags a track and all its family tree to be kept
708 particle=(TParticle*)fParticleMap->At(curr);
710 // If the particle is flagged the three from here upward is saved already
711 if(particle->TestBit(kKeepBit)) return;
713 // Save this particle
714 particle->SetBit(kKeepBit);
716 // Move to father if any
717 if((curr=particle->GetFirstMother())==-1) return;
721 //_____________________________________________________________________________
722 void AliRun::EnergySummary()
725 // Print summary of deposited energy
731 Int_t kn, i, left, j, id;
732 const Float_t kzero=0;
733 Int_t ievent=fHeader.GetEvent()+1;
735 // Energy loss information
737 printf("***************** Energy Loss Information per event (GEV) *****************\n");
738 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
741 fEventEnergy[ndep]=kn;
746 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
749 fSummEnergy[ndep]=ed;
750 fSum2Energy[ndep]=TMath::Min((Float_t) 99.,TMath::Max(ed2,kzero));
755 for(kn=0;kn<(ndep-1)/3+1;kn++) {
757 for(i=0;i<(3<left?3:left);i++) {
759 id=Int_t (fEventEnergy[j]+0.1);
760 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
765 // Relative energy loss in different detectors
766 printf("******************** Relative Energy Loss per event ********************\n");
767 printf("Total energy loss per event %10.3f GeV\n",edtot);
768 for(kn=0;kn<(ndep-1)/5+1;kn++) {
770 for(i=0;i<(5<left?5:left);i++) {
772 id=Int_t (fEventEnergy[j]+0.1);
773 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
777 for(kn=0;kn<75;kn++) printf("*");
781 // Reset the TArray's
782 // fEventEnergy.Set(0);
783 // fSummEnergy.Set(0);
784 // fSum2Energy.Set(0);
787 //_____________________________________________________________________________
788 AliModule *AliRun::GetModule(const char *name) const
791 // Return pointer to detector from name
793 return (AliModule*)fModules->FindObject(name);
796 //_____________________________________________________________________________
797 AliDetector *AliRun::GetDetector(const char *name) const
800 // Return pointer to detector from name
802 return (AliDetector*)fModules->FindObject(name);
805 //_____________________________________________________________________________
806 Int_t AliRun::GetModuleID(const char *name) const
809 // Return galice internal detector identifier from name
812 TObject *mod=fModules->FindObject(name);
813 if(mod) i=fModules->IndexOf(mod);
817 //_____________________________________________________________________________
818 Int_t AliRun::GetEvent(Int_t event)
821 // Connect the Trees Kinematics and Hits for event # event
822 // Set branch addresses
825 // Reset existing structures
831 // Delete Trees already connected
832 if (fTreeK) delete fTreeK;
833 if (fTreeH) delete fTreeH;
834 if (fTreeD) delete fTreeD;
835 if (fTreeR) delete fTreeR;
836 if (fTreeS) delete fTreeS;
838 // Get header from file
839 if(fTreeE) fTreeE->GetEntry(event);
840 else Error("GetEvent","Cannot file Header Tree\n");
841 TFile *file = fTreeE->GetCurrentFile();
845 // Get Kine Tree from file
847 sprintf(treeName,"TreeK%d",event);
848 fTreeK = (TTree*)gDirectory->Get(treeName);
849 if (fTreeK) fTreeK->SetBranchAddress("Particles", &fParticleBuffer);
850 else Error("GetEvent","cannot find Kine Tree for event:%d\n",event);
851 // Create the particle stack
852 if(!fParticles) fParticles = new TClonesArray("TParticle",1000);
853 // Build the pointer list
855 fParticleMap->Clear();
856 fParticleMap->Expand(fTreeK->GetEntries());
858 fParticleMap = new TObjArray(fTreeK->GetEntries());
862 // Get Hits Tree header from file
863 sprintf(treeName,"TreeH%d",event);
864 fTreeH = (TTree*)gDirectory->Get(treeName);
866 Error("GetEvent","cannot find Hits Tree for event:%d\n",event);
871 // Get Digits Tree header from file
872 sprintf(treeName,"TreeD%d",event);
873 fTreeD = (TTree*)gDirectory->Get(treeName);
875 // Warning("GetEvent","cannot find Digits Tree for event:%d\n",event);
880 // Get SDigits Tree header from file
881 sprintf(treeName,"TreeS%d",event);
882 fTreeS = (TTree*)gDirectory->Get(treeName);
884 // Warning("GetEvent","cannot find SDigits Tree for event:%d\n",event);
889 // Get Reconstruct Tree header from file
890 sprintf(treeName,"TreeR%d",event);
891 fTreeR = (TTree*)gDirectory->Get(treeName);
893 // printf("WARNING: cannot find Reconstructed Tree for event:%d\n",event);
898 // Set Trees branch addresses
899 TIter next(fModules);
901 while((detector = (AliModule*)next())) {
902 detector->SetTreeAddress();
905 fNtrack = Int_t (fTreeK->GetEntries());
909 //_____________________________________________________________________________
910 TGeometry *AliRun::GetGeometry()
913 // Import Alice geometry from current file
914 // Return pointer to geometry object
916 if (!fGeometry) fGeometry = (TGeometry*)gDirectory->Get("AliceGeom");
918 // Unlink and relink nodes in detectors
919 // This is bad and there must be a better way...
922 TIter next(fModules);
924 while((detector = (AliModule*)next())) {
925 TList *dnodes=detector->Nodes();
928 for ( j=0; j<dnodes->GetSize(); j++) {
929 node = (TNode*) dnodes->At(j);
930 node1 = fGeometry->GetNode(node->GetName());
931 dnodes->Remove(node);
932 dnodes->AddAt(node1,j);
938 //_____________________________________________________________________________
939 void AliRun::GetNextTrack(Int_t &mtrack, Int_t &ipart, Float_t *pmom,
940 Float_t &e, Float_t *vpos, Float_t *polar,
944 // Return next track from stack of particles
949 for(Int_t i=fNtrack-1; i>=0; i--) {
950 track=(TParticle*) fParticleMap->At(i);
951 if(track) if(!track->TestBit(kDoneBit)) {
953 // The track exists and has not yet been processed
955 ipart=track->GetPdgCode();
963 track->GetPolarisation(pol);
968 track->SetBit(kDoneBit);
974 // stop and start timer when we start a primary track
975 Int_t nprimaries = fHeader.GetNprimary();
976 if (fCurrent >= nprimaries) return;
977 if (fCurrent < nprimaries-1) {
979 track=(TParticle*) fParticleMap->At(fCurrent+1);
980 // track->SetProcessTime(fTimer.CpuTime());
985 //_____________________________________________________________________________
986 Int_t AliRun::GetPrimary(Int_t track) const
989 // return number of primary that has generated track
997 part = (TParticle *)fParticleMap->At(current);
998 parent=part->GetFirstMother();
999 if(parent<0) return current;
1003 //_____________________________________________________________________________
1004 void AliRun::InitMC(const char *setup)
1007 // Initialize the Alice setup
1011 Warning("Init","Cannot initialise AliRun twice!\n");
1015 gROOT->LoadMacro(setup);
1016 gInterpreter->ProcessLine(fConfigFunction.Data());
1019 gMC->DefineParticles(); //Create standard MC particles
1021 TObject *objfirst, *objlast;
1023 fNdets = fModules->GetLast()+1;
1026 //=================Create Materials and geometry
1029 // Added also after in case of interactive initialisation of modules
1030 fNdets = fModules->GetLast()+1;
1032 TIter next(fModules);
1033 AliModule *detector;
1034 while((detector = (AliModule*)next())) {
1035 detector->SetTreeAddress();
1036 objlast = gDirectory->GetList()->Last();
1038 // Add Detector histograms in Detector list of histograms
1039 if (objlast) objfirst = gDirectory->GetList()->After(objlast);
1040 else objfirst = gDirectory->GetList()->First();
1042 detector->Histograms()->Add(objfirst);
1043 objfirst = gDirectory->GetList()->After(objfirst);
1046 ReadTransPar(); //Read the cuts for all materials
1048 MediaTable(); //Build the special IMEDIA table
1050 //Initialise geometry deposition table
1051 fEventEnergy.Set(gMC->NofVolumes()+1);
1052 fSummEnergy.Set(gMC->NofVolumes()+1);
1053 fSum2Energy.Set(gMC->NofVolumes()+1);
1055 //Compute cross-sections
1056 gMC->BuildPhysics();
1058 //Write Geometry object to current file.
1063 fMCQA = new AliMCQA(fNdets);
1066 // Save stuff at the beginning of the file to avoid file corruption
1070 //_____________________________________________________________________________
1071 void AliRun::MediaTable()
1074 // Built media table to get from the media number to
1077 Int_t kz, nz, idt, lz, i, k, ind;
1079 TObjArray &dets = *gAlice->Detectors();
1082 // For all detectors
1083 for (kz=0;kz<fNdets;kz++) {
1084 // If detector is defined
1085 if((det=(AliModule*) dets[kz])) {
1086 TArrayI &idtmed = *(det->GetIdtmed());
1087 for(nz=0;nz<100;nz++) {
1088 // Find max and min material number
1089 if((idt=idtmed[nz])) {
1090 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
1091 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
1094 if(det->LoMedium() > det->HiMedium()) {
1095 det->LoMedium() = 0;
1096 det->HiMedium() = 0;
1098 if(det->HiMedium() > fImedia->GetSize()) {
1099 Error("MediaTable","Increase fImedia from %d to %d",
1100 fImedia->GetSize(),det->HiMedium());
1103 // Tag all materials in rage as belonging to detector kz
1104 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
1111 // Print summary table
1112 printf(" Traking media ranges:\n");
1113 for(i=0;i<(fNdets-1)/6+1;i++) {
1114 for(k=0;k< (6<fNdets-i*6?6:fNdets-i*6);k++) {
1116 det=(AliModule*)dets[ind];
1118 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
1121 printf(" %6s: %3d -> %3d;","NULL",0,0);
1127 //____________________________________________________________________________
1128 void AliRun::SetGenerator(AliGenerator *generator)
1131 // Load the event generator
1133 if(!fGenerator) fGenerator = generator;
1136 //____________________________________________________________________________
1137 void AliRun::ResetGenerator(AliGenerator *generator)
1140 // Load the event generator
1144 Warning("ResetGenerator","Replacing generator %s with %s\n",
1145 fGenerator->GetName(),generator->GetName());
1147 Warning("ResetGenerator","Replacing generator %s with NULL\n",
1148 fGenerator->GetName());
1149 fGenerator = generator;
1152 //____________________________________________________________________________
1153 void AliRun::SetTransPar(char *filename)
1155 fTransParName = filename;
1158 //____________________________________________________________________________
1159 void AliRun::SetBaseFile(char *filename)
1161 fBaseFileName = filename;
1164 //____________________________________________________________________________
1165 void AliRun::ReadTransPar()
1168 // Read filename to set the transport parameters
1172 const Int_t kncuts=10;
1173 const Int_t knflags=11;
1174 const Int_t knpars=kncuts+knflags;
1175 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
1176 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
1177 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
1178 "MULS","PAIR","PHOT","RAYL"};
1182 Float_t cut[kncuts];
1183 Int_t flag[knflags];
1184 Int_t i, itmed, iret, ktmed, kz;
1187 // See whether the file is there
1188 filtmp=gSystem->ExpandPathName(fTransParName.Data());
1189 lun=fopen(filtmp,"r");
1192 Warning("ReadTransPar","File %s does not exist!\n",fTransParName.Data());
1196 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
1197 printf(" *%59s\n","*");
1198 printf(" * Please check carefully what you are doing!%10s\n","*");
1199 printf(" *%59s\n","*");
1202 // Initialise cuts and flags
1203 for(i=0;i<kncuts;i++) cut[i]=-99;
1204 for(i=0;i<knflags;i++) flag[i]=-99;
1206 for(i=0;i<256;i++) line[i]='\0';
1207 // Read up to the end of line excluded
1208 iret=fscanf(lun,"%[^\n]",line);
1212 printf(" *%59s\n","*");
1213 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
1216 // Read the end of line
1219 if(line[0]=='*') continue;
1221 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",
1222 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
1223 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
1224 &flag[8],&flag[9],&flag[10]);
1228 Warning("ReadTransPar","Error reading file %s\n",fTransParName.Data());
1231 // Check that the module exist
1232 AliModule *mod = GetModule(detName);
1234 // Get the array of media numbers
1235 TArrayI &idtmed = *mod->GetIdtmed();
1236 // Check that the tracking medium code is valid
1237 if(0<=itmed && itmed < 100) {
1238 ktmed=idtmed[itmed];
1240 Warning("ReadTransPar","Invalid tracking medium code %d for %s\n",itmed,mod->GetName());
1243 // Set energy thresholds
1244 for(kz=0;kz<kncuts;kz++) {
1246 printf(" * %-6s set to %10.3E for tracking medium code %4d for %s\n",
1247 kpars[kz],cut[kz],itmed,mod->GetName());
1248 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
1251 // Set transport mechanisms
1252 for(kz=0;kz<knflags;kz++) {
1254 printf(" * %-6s set to %10d for tracking medium code %4d for %s\n",
1255 kpars[kncuts+kz],flag[kz],itmed,mod->GetName());
1256 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
1260 Warning("ReadTransPar","Invalid medium code %d *\n",itmed);
1264 Warning("ReadTransPar","Module %s not present\n",detName);
1270 //_____________________________________________________________________________
1271 void AliRun::MakeBranchInTree(TTree *tree, const char* name, void* address, Int_t size, char *file)
1274 printf("* MakeBranch * Making Branch %s \n",name);
1276 TBranch *branch = tree->Branch(name,address,size);
1279 char * outFile = new char[strlen(gAlice->GetBaseFile())+strlen(file)+2];
1280 sprintf(outFile,"%s/%s",gAlice->GetBaseFile(),file);
1281 TDirectory *cwd = gDirectory;
1282 branch->SetFile(outFile);
1283 TIter next( branch->GetListOfBranches());
1284 while ((branch=(TBranch*)next())) {
1285 branch->SetFile(outFile);
1288 printf("* MakeBranch * Diverting Branch %s to file %s\n",name,file);
1294 //_____________________________________________________________________________
1295 void AliRun::MakeBranchInTree(TTree *tree, const char* name, const char *classname, void* address, Int_t size, Int_t splitlevel, char *file)
1297 TDirectory *cwd = gDirectory;
1298 TBranch *branch = tree->Branch(name,classname,address,size,splitlevel);
1300 printf("* MakeBranch * Making Branch %s \n",name);
1302 char * outFile = new char[strlen(gAlice->GetBaseFile())+strlen(file)+2];
1303 sprintf(outFile,"%s/%s",gAlice->GetBaseFile(),file);
1304 branch->SetFile(outFile);
1305 TIter next( branch->GetListOfBranches());
1306 while ((branch=(TBranch*)next())) {
1307 branch->SetFile(outFile);
1310 printf("* MakeBranch * Diverting Branch %s to file %s\n",name,file);
1315 //_____________________________________________________________________________
1316 void AliRun::MakeTree(Option_t *option, char *file)
1319 // Create the ROOT trees
1320 // Loop on all detectors to create the Root branch (if any)
1326 char *oK = strstr(option,"K");
1327 char *oH = strstr(option,"H");
1328 char *oE = strstr(option,"E");
1329 char *oD = strstr(option,"D");
1330 char *oR = strstr(option,"R");
1331 char *oS = strstr(option,"S");
1334 if (oK && !fTreeK) {
1335 sprintf(hname,"TreeK%d",fEvent);
1336 fTreeK = new TTree(hname,"Kinematics");
1337 // Create a branch for particles
1338 MakeBranchInTree(fTreeK,
1339 "Particles", "TParticle", &fParticleBuffer, 4000, 1, file) ;
1342 if (oH && !fTreeH) {
1343 sprintf(hname,"TreeH%d",fEvent);
1344 fTreeH = new TTree(hname,"Hits");
1345 fTreeH->SetAutoSave(1000000000); //no autosave
1348 if (oD && !fTreeD) {
1349 sprintf(hname,"TreeD%d",fEvent);
1350 fTreeD = new TTree(hname,"Digits");
1353 if (oS && !fTreeS) {
1354 sprintf(hname,"TreeS%d",fEvent);
1355 fTreeS = new TTree(hname,"SDigits");
1358 if (oR && !fTreeR) {
1359 sprintf(hname,"TreeR%d",fEvent);
1360 fTreeR = new TTree(hname,"Reconstruction");
1363 if (oE && !fTreeE) {
1364 fTreeE = new TTree("TE","Header");
1365 // Create a branch for Header
1366 MakeBranchInTree(fTreeE,
1367 "Header", "AliHeader", &gAliHeader, 4000, 1, file) ;
1372 // Create a branch for hits/digits for each detector
1373 // Each branch is a TClonesArray. Each data member of the Hits classes
1374 // will be in turn a subbranch of the detector master branch
1375 TIter next(fModules);
1376 AliModule *detector;
1377 while((detector = (AliModule*)next())) {
1378 if (oH || oR) detector->MakeBranch(option,file);
1382 //_____________________________________________________________________________
1383 Int_t AliRun::PurifyKine(Int_t lastSavedTrack, Int_t nofTracks)
1386 // PurifyKine with external parameters
1388 fHgwmk = lastSavedTrack;
1389 fNtrack = nofTracks;
1394 //_____________________________________________________________________________
1395 TParticle* AliRun::Particle(Int_t i)
1397 if(!(*fParticleMap)[i]) {
1398 Int_t nentries = fParticles->GetEntries();
1399 fTreeK->GetEntry(fParticleFileMap[i]);
1400 new ((*fParticles)[nentries]) TParticle(*fParticleBuffer);
1401 fParticleMap->AddAt((*fParticles)[nentries],i);
1403 return (TParticle *) (*fParticleMap)[i];
1406 //_____________________________________________________________________________
1407 void AliRun::PurifyKine()
1410 // Compress kinematic tree keeping only flagged particles
1411 // and renaming the particle id's in all the hits
1413 // TClonesArray &particles = *fParticles;
1414 TObjArray &particles = *fParticleMap;
1415 int nkeep=fHgwmk+1, parent, i;
1416 TParticle *part, *father;
1417 TArrayI map(particles.GetLast()+1);
1419 // Save in Header total number of tracks before compression
1420 fHeader.SetNtrack(fHeader.GetNtrack()+fNtrack-fHgwmk);
1422 // If no tracks generated return now
1423 if(fHgwmk+1 == fNtrack) return;
1425 Int_t toshrink = fNtrack-fHgwmk-1;
1427 // First pass, invalid Daughter information
1428 for(i=0; i<fNtrack; i++) {
1429 // Preset map, to be removed later
1430 if(i<=fHgwmk) map[i]=i ;
1433 // particles.UncheckedAt(i)->ResetBit(kDaughtersBit);
1434 if((part=(TParticle*) particles.At(i))) part->ResetBit(kDaughtersBit);
1437 // Invalid daughter information for the parent of the first particle
1438 // generated. This may or may not be the current primary according to
1439 // whether decays have been recorded among the primaries
1440 part = (TParticle *)particles.At(fHgwmk+1);
1441 particles.At(part->GetFirstMother())->ResetBit(kDaughtersBit);
1442 // Second pass, build map between old and new numbering
1443 for(i=fHgwmk+1; i<fNtrack; i++) {
1444 if(particles.At(i)->TestBit(kKeepBit)) {
1446 // This particle has to be kept
1448 // If old and new are different, have to move the pointer
1449 if(i!=nkeep) particles[nkeep]=particles.At(i);
1450 part = (TParticle*) particles.At(nkeep);
1452 // as the parent is always *before*, it must be already
1453 // in place. This is what we are checking anyway!
1454 if((parent=part->GetFirstMother())>fHgwmk)
1455 if(map[parent]==-99) Fatal("PurifyKine","map[%d] = -99!\n",parent);
1456 else part->SetFirstMother(map[parent]);
1462 // Fix daughters information
1463 for (i=fHgwmk+1; i<nkeep; i++) {
1464 part = (TParticle *)particles.At(i);
1465 parent = part->GetFirstMother();
1467 father = (TParticle *)particles.At(parent);
1468 if(father->TestBit(kDaughtersBit)) {
1470 if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
1471 if(i>father->GetLastDaughter()) father->SetLastDaughter(i);
1473 // Initialise daughters info for first pass
1474 father->SetFirstDaughter(i);
1475 father->SetLastDaughter(i);
1476 father->SetBit(kDaughtersBit);
1481 // Now loop on all registered hit lists
1482 TIter next(fHitLists);
1483 TCollection *hitList;
1484 while((hitList = (TCollection*)next())) {
1485 TIter nexthit(hitList);
1487 while((hit = (AliHit*)nexthit())) {
1488 hit->SetTrack(map[hit->GetTrack()]);
1493 // This for detectors which have a special mapping mechanism
1494 // for hits, such as TPC and TRD
1497 TIter nextmod(fModules);
1498 AliModule *detector;
1499 while((detector = (AliModule*)nextmod())) {
1500 detector->RemapTrackHitIDs(map.GetArray());
1503 // Now the output bit, from fHgwmk to nkeep we write everything and we erase
1504 if(nkeep>fParticleFileMap.GetSize()) fParticleFileMap.Set(Int_t (nkeep*1.5));
1507 for (i=fHgwmk+1; i<nkeep; ++i) {
1508 fParticleBuffer = (TParticle*) particles.At(i);
1509 fParticleFileMap[i]=(Int_t) fTreeK->GetEntries();
1514 for (i=nkeep; i<fNtrack; ++i) particles[i]=0;
1516 fLoadPoint-=toshrink;
1517 for(i=fLoadPoint; i<fLoadPoint+toshrink; ++i) fParticles->RemoveAt(i);
1524 //_____________________________________________________________________________
1525 void AliRun::BeginEvent()
1528 // Reset all Detectors & kinematics & trees
1535 fLego->BeginEvent();
1545 // Initialise event header
1546 fHeader.Reset(fRun,fEvent);
1550 sprintf(hname,"TreeK%d",fEvent);
1551 fTreeK->SetName(hname);
1555 sprintf(hname,"TreeH%d",fEvent);
1556 fTreeH->SetName(hname);
1560 sprintf(hname,"TreeD%d",fEvent);
1561 fTreeD->SetName(hname);
1565 sprintf(hname,"TreeS%d",fEvent);
1566 fTreeS->SetName(hname);
1570 sprintf(hname,"TreeR%d",fEvent);
1571 fTreeR->SetName(hname);
1574 //_____________________________________________________________________________
1575 void AliRun::ResetDigits()
1578 // Reset all Detectors digits
1580 TIter next(fModules);
1581 AliModule *detector;
1582 while((detector = (AliModule*)next())) {
1583 detector->ResetDigits();
1587 //_____________________________________________________________________________
1588 void AliRun::ResetSDigits()
1591 // Reset all Detectors digits
1593 TIter next(fModules);
1594 AliModule *detector;
1595 while((detector = (AliModule*)next())) {
1596 detector->ResetSDigits();
1600 //_____________________________________________________________________________
1601 void AliRun::ResetHits()
1604 // Reset all Detectors hits
1606 TIter next(fModules);
1607 AliModule *detector;
1608 while((detector = (AliModule*)next())) {
1609 detector->ResetHits();
1613 //_____________________________________________________________________________
1614 void AliRun::ResetPoints()
1617 // Reset all Detectors points
1619 TIter next(fModules);
1620 AliModule *detector;
1621 while((detector = (AliModule*)next())) {
1622 detector->ResetPoints();
1626 //_____________________________________________________________________________
1627 void AliRun::RunMC(Int_t nevent, const char *setup)
1630 // Main function to be called to process a galice run
1632 // Root > gAlice.Run();
1633 // a positive number of events will cause the finish routine
1637 // check if initialisation has been done
1638 if (!fInitDone) InitMC(setup);
1640 // Create the Root Tree with one branch per detector
1642 if (gSystem->Getenv("CONFIG_SPLIT_FILE")) {
1644 MakeTree("K","Kine.root");
1645 MakeTree("H","Hits.root");
1646 MakeTree("R","Reco.root");
1651 gMC->ProcessRun(nevent);
1653 // End of this run, close files
1654 if(nevent>0) FinishRun();
1657 //_____________________________________________________________________________
1659 void AliRun::Hits2Digits(const char *selected)
1661 Hits2SDigits(selected);
1662 SDigits2Digits(selected);
1665 //_____________________________________________________________________________
1667 void AliRun::Hits2SDigits(const char *selected)
1670 // Main function to be called to convert hits to digits.
1672 gAlice->GetEvent(0);
1674 TObjArray *detectors = gAlice->Detectors();
1676 TIter next(detectors);
1678 AliDetector *detector;
1680 TDirectory *cwd = gDirectory;
1684 while((detector = (AliDetector*)next())) {
1686 if (strcmp(detector->GetName(),selected)) continue;
1688 if (detector->IsActive()){
1689 if (gSystem->Getenv("CONFIG_SPLIT_FILE")) {
1691 cout << "Processing " << detector->GetName() << "..." << endl;
1692 char * outFile = new char[strlen (detector->GetName())+18];
1693 sprintf(outFile,"SDigits.%s.root",detector->GetName());
1694 detector->MakeBranch("S",outFile);
1697 detector->MakeBranch("S");
1700 detector->Hits2SDigits();
1705 //_____________________________________________________________________________
1707 void AliRun::SDigits2Digits(const char *selected)
1710 // Main function to be called to convert hits to digits.
1712 gAlice->GetEvent(0);
1714 TObjArray *detectors = gAlice->Detectors();
1716 TIter next(detectors);
1718 AliDetector *detector;
1720 TDirectory *cwd = gDirectory;
1724 while((detector = (AliDetector*)next())) {
1726 if (strcmp(detector->GetName(),selected)) continue;
1728 if (detector->IsActive()){
1729 if (gSystem->Getenv("CONFIG_SPLIT_FILE")) {
1731 cout << "Processing " << detector->GetName() << "..." << endl;
1732 char * outFile = new char[strlen (detector->GetName())+16];
1733 sprintf(outFile,"Digits.%s.root",detector->GetName());
1734 detector->MakeBranch("D",outFile);
1737 detector->MakeBranch("D");
1740 detector->SDigits2Digits();
1745 //_____________________________________________________________________________
1746 void AliRun::RunLego(const char *setup, Int_t nc1, Float_t c1min,
1747 Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
1748 Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener)
1751 // Generates lego plots of:
1752 // - radiation length map phi vs theta
1753 // - radiation length map phi vs eta
1754 // - interaction length map
1755 // - g/cm2 length map
1757 // ntheta bins in theta, eta
1758 // themin minimum angle in theta (degrees)
1759 // themax maximum angle in theta (degrees)
1761 // phimin minimum angle in phi (degrees)
1762 // phimax maximum angle in phi (degrees)
1763 // rmin minimum radius
1764 // rmax maximum radius
1767 // The number of events generated = ntheta*nphi
1768 // run input parameters in macro setup (default="Config.C")
1770 // Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
1773 <img src="picts/AliRunLego1.gif">
1778 <img src="picts/AliRunLego2.gif">
1783 <img src="picts/AliRunLego3.gif">
1788 // check if initialisation has been done
1789 if (!fInitDone) InitMC(setup);
1790 //Save current generator
1791 AliGenerator *gen=Generator();
1793 // Set new generator
1794 if (!gener) gener = new AliLegoGenerator();
1795 ResetGenerator(gener);
1797 // Configure Generator
1798 gener->SetRadiusRange(rmin, rmax);
1799 gener->SetZMax(zmax);
1800 gener->SetCoor1Range(nc1, c1min, c1max);
1801 gener->SetCoor2Range(nc2, c2min, c2max);
1804 //Create Lego object
1805 fLego = new AliLego("lego",gener);
1807 //Prepare MC for Lego Run
1812 gMC->ProcessRun(nc1*nc2+1);
1814 // Create only the Root event Tree
1817 // End of this run, close files
1819 // Restore current generator
1820 ResetGenerator(gen);
1821 // Delete Lego Object
1822 delete fLego; fLego=0;
1825 //_____________________________________________________________________________
1826 void AliRun::SetConfigFunction(const char * config)
1829 // Set the signature of the function contained in Config.C to configure
1832 fConfigFunction=config;
1835 //_____________________________________________________________________________
1836 void AliRun::SetCurrentTrack(Int_t track)
1839 // Set current track number
1844 //_____________________________________________________________________________
1845 void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
1846 Float_t *vpos, Float_t *polar, Float_t tof,
1847 AliMCProcess mech, Int_t &ntr, Float_t weight)
1850 // Load a track on the stack
1852 // done 0 if the track has to be transported
1854 // parent identifier of the parent track. -1 for a primary
1855 // pdg particle code
1856 // pmom momentum GeV/c
1858 // polar polarisation
1859 // tof time of flight in seconds
1860 // mecha production mechanism
1861 // ntr on output the number of the track stored
1863 TClonesArray &particles = *fParticles;
1864 TParticle *particle;
1866 const Int_t kfirstdaughter=-1;
1867 const Int_t klastdaughter=-1;
1869 // const Float_t tlife=0;
1872 // Here we get the static mass
1873 // For MC is ok, but a more sophisticated method could be necessary
1874 // if the calculated mass is required
1875 // also, this method is potentially dangerous if the mass
1876 // used in the MC is not the same of the PDG database
1878 mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
1879 Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
1880 pmom[1]*pmom[1]+pmom[2]*pmom[2]);
1882 //printf("Loading particle %s mass %f ene %f No %d ip %d pos %f %f %f mom %f %f %f kS %d m %s\n",
1883 //pname,mass,e,fNtrack,pdg,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],kS,mecha);
1885 particle=new(particles[fLoadPoint++]) TParticle(pdg,kS,parent,-1,kfirstdaughter,
1886 klastdaughter,pmom[0],pmom[1],pmom[2],
1887 e,vpos[0],vpos[1],vpos[2],tof);
1888 particle->SetPolarisation(TVector3(polar[0],polar[1],polar[2]));
1889 particle->SetWeight(weight);
1890 particle->SetUniqueID(mech);
1891 if(!done) particle->SetBit(kDoneBit);
1892 // Declare that the daughter information is valid
1893 particle->SetBit(kDaughtersBit);
1894 // Add the particle to the stack
1895 fParticleMap->AddAtAndExpand(particle,fNtrack);
1898 particle=(TParticle*) fParticleMap->At(parent);
1899 particle->SetLastDaughter(fNtrack);
1900 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
1903 // This is a primary track. Set high water mark for this event
1906 // Set also number if primary tracks
1907 fHeader.SetNprimary(fHgwmk+1);
1908 fHeader.SetNtrack(fHgwmk+1);
1914 // Here we get the static mass
1915 // For MC is ok, but a more sophisticated method could be necessary
1916 // if the calculated mass is required
1917 // also, this method is potentially dangerous if the mass
1918 // used in the MC is not the same of the PDG database
1920 Float_t mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
1921 Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
1922 pmom[1]*pmom[1]+pmom[2]*pmom[2]);
1924 SetTrack(done, parent, pdg, pmom[0], pmom[1], pmom[2], e,
1925 vpos[0], vpos[1], vpos[2], tof, polar[0],polar[1],polar[2],
1930 //_____________________________________________________________________________
1931 void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg,
1932 Double_t px, Double_t py, Double_t pz, Double_t e,
1933 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
1934 Double_t polx, Double_t poly, Double_t polz,
1935 AliMCProcess mech, Int_t &ntr, Float_t weight)
1938 // Load a track on the stack
1940 // done 0 if the track has to be transported
1942 // parent identifier of the parent track. -1 for a primary
1943 // pdg particle code
1944 // kS generation status code
1945 // px, py, pz momentum GeV/c
1946 // vx, vy, vz position
1947 // polar polarisation
1948 // tof time of flight in seconds
1949 // mech production mechanism
1950 // ntr on output the number of the track stored
1952 // New method interface:
1953 // arguments were changed to be in correspondence with TParticle
1955 // Note: the energy is not calculated from the static mass but
1956 // it is passed by argument e.
1958 TClonesArray &particles = *fParticles;
1961 const Int_t kFirstDaughter=-1;
1962 const Int_t kLastDaughter=-1;
1965 = new(particles[fLoadPoint++]) TParticle(pdg, kS, parent, -1,
1966 kFirstDaughter, kLastDaughter,
1967 px, py, pz, e, vx, vy, vz, tof);
1969 particle->SetPolarisation(polx, poly, polz);
1970 particle->SetWeight(weight);
1971 particle->SetUniqueID(mech);
1973 if(!done) particle->SetBit(kDoneBit);
1975 // Declare that the daughter information is valid
1976 particle->SetBit(kDaughtersBit);
1977 // Add the particle to the stack
1978 fParticleMap->AddAtAndExpand(particle,fNtrack);
1981 particle=(TParticle*) fParticleMap->At(parent);
1982 particle->SetLastDaughter(fNtrack);
1983 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
1986 // This is a primary track. Set high water mark for this event
1989 // Set also number if primary tracks
1990 fHeader.SetNprimary(fHgwmk+1);
1991 fHeader.SetNtrack(fHgwmk+1);
1996 //_____________________________________________________________________________
1997 void AliRun::SetHighWaterMark(const Int_t nt)
2000 // Set high water mark for last track in event
2003 // Set also number if primary tracks
2004 fHeader.SetNprimary(fHgwmk+1);
2005 fHeader.SetNtrack(fHgwmk+1);
2008 //_____________________________________________________________________________
2009 void AliRun::KeepTrack(const Int_t track)
2012 // flags a track to be kept
2014 fParticleMap->At(track)->SetBit(kKeepBit);
2017 //_____________________________________________________________________________
2018 void AliRun::StepManager(Int_t id)
2021 // Called at every step during transport
2025 // --- If lego option, do it and leave
2027 fLego->StepManager();
2030 //Update energy deposition tables
2031 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
2033 //Call the appropriate stepping routine;
2034 AliModule *det = (AliModule*)fModules->At(id);
2036 fMCQA->StepManager(id);
2042 //_____________________________________________________________________________
2043 void AliRun::Streamer(TBuffer &R__b)
2045 // Stream an object of class AliRun.
2047 if (R__b.IsReading()) {
2048 if (!gAlice) gAlice = this;
2050 AliRun::Class()->ReadBuffer(R__b, this);
2052 gROOT->GetListOfBrowsables()->Add(this,"Run");
2054 fTreeE = (TTree*)gDirectory->Get("TE");
2055 if (fTreeE) fTreeE->SetBranchAddress("Header", &gAliHeader);
2056 else Error("Streamer","cannot find Header Tree\n");
2057 fTreeE->GetEntry(0);
2061 AliRun::Class()->WriteBuffer(R__b, this);