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.55 2001/02/12 15:52:54 buncic
19 Removed OpenBaseFile().
21 Revision 1.54 2001/02/07 10:39:05 hristov
22 Remove default value for argument
24 Revision 1.53 2001/02/06 11:02:26 hristov
25 New SetTrack interface added, added check for unfilled particles in FinishEvent (I.Hrivnacova)
27 Revision 1.52 2001/02/05 16:22:25 buncic
28 Added TreeS to GetEvent().
30 Revision 1.51 2001/02/02 15:16:20 morsch
31 SetHighWaterMark method added to mark last particle in event.
33 Revision 1.50 2001/01/27 10:32:00 hristov
34 Leave the loop when primaries are filled (I.Hrivnacova)
36 Revision 1.49 2001/01/26 19:58:48 hristov
37 Major upgrade of AliRoot code
39 Revision 1.48 2001/01/17 10:50:50 hristov
40 Corrections to destructors
42 Revision 1.47 2000/12/18 10:44:01 morsch
43 Possibility to set field map by passing pointer to objet of type AliMagF via
46 gAlice->SetField(new AliMagFCM("Map2", "$(ALICE_ROOT)/data/field01.dat",2,1.,10.));
48 Revision 1.46 2000/12/14 19:29:27 fca
49 galice.cuts was not read any more
51 Revision 1.45 2000/11/30 07:12:49 alibrary
52 Introducing new Rndm and QA classes
54 Revision 1.44 2000/10/26 13:58:59 morsch
55 Add possibility to choose the lego generator (of type AliGeneratorLego or derived) when running
56 RunLego(). Default is the base class AliGeneratorLego.
58 Revision 1.43 2000/10/09 09:43:17 fca
59 Special remapping of hits for TPC and TRD. End-of-primary action introduced
61 Revision 1.42 2000/10/02 21:28:14 fca
62 Removal of useless dependecies via forward declarations
64 Revision 1.41 2000/07/13 16:19:09 fca
65 Mainly coding conventions + some small bug fixes
67 Revision 1.40 2000/07/12 08:56:25 fca
68 Coding convention correction and warning removal
70 Revision 1.39 2000/07/11 18:24:59 fca
71 Coding convention corrections + few minor bug fixes
73 Revision 1.38 2000/06/20 13:05:45 fca
74 Writing down the TREE headers before job starts
76 Revision 1.37 2000/06/09 20:05:11 morsch
77 Introduce possibility to chose magnetic field version 3: AliMagFDM + field02.dat
79 Revision 1.36 2000/06/08 14:03:58 hristov
80 Only one initializer for a default argument
82 Revision 1.35 2000/06/07 10:13:14 hristov
83 Delete only existent objects.
85 Revision 1.34 2000/05/18 10:45:38 fca
86 Delete Particle Factory properly
88 Revision 1.33 2000/05/16 13:10:40 fca
89 New method IsNewTrack and fix for a problem in Father-Daughter relations
91 Revision 1.32 2000/04/27 10:38:21 fca
92 Correct termination of Lego Run and introduce Lego getter in AliRun
94 Revision 1.31 2000/04/26 10:17:32 fca
95 Changes in Lego for G4 compatibility
97 Revision 1.30 2000/04/18 19:11:40 fca
98 Introduce variable Config.C function signature
100 Revision 1.29 2000/04/07 11:12:34 fca
101 G4 compatibility changes
103 Revision 1.28 2000/04/05 06:51:06 fca
104 Workaround for an HP compiler problem
106 Revision 1.27 2000/03/22 18:08:07 fca
107 Rationalisation of the virtual MC interfaces
109 Revision 1.26 2000/03/22 13:42:26 fca
110 SetGenerator does not replace an existing generator, ResetGenerator does
112 Revision 1.25 2000/02/23 16:25:22 fca
113 AliVMC and AliGeant3 classes introduced
114 ReadEuclid moved from AliRun to AliModule
116 Revision 1.24 2000/01/19 17:17:20 fca
117 Introducing a list of lists of hits -- more hits allowed for detector now
119 Revision 1.23 1999/12/03 11:14:31 fca
120 Fixing previous wrong checking
122 Revision 1.21 1999/11/25 10:40:08 fca
123 Fixing daughters information also in primary tracks
125 Revision 1.20 1999/10/04 18:08:49 fca
126 Adding protection against inconsistent Euclid files
128 Revision 1.19 1999/09/29 07:50:40 fca
129 Introduction of the Copyright and cvs Log
133 ///////////////////////////////////////////////////////////////////////////////
135 // Control class for Alice C++ //
136 // Only one single instance of this class exists. //
137 // The object is created in main program aliroot //
138 // and is pointed by the global gAlice. //
140 // -Supports the list of all Alice Detectors (fModules). //
141 // -Supports the list of particles (fParticles). //
142 // -Supports the Trees. //
143 // -Supports the geometry. //
144 // -Supports the event display. //
147 <img src="picts/AliRunClass.gif">
152 <img src="picts/alirun.gif">
156 ///////////////////////////////////////////////////////////////////////////////
161 #include <iostream.h>
169 #include <TObjectTable.h>
171 #include <TGeometry.h>
173 #include "TBrowser.h"
175 #include "TParticle.h"
177 #include "AliDisplay.h"
180 #include "AliMagFC.h"
181 #include "AliMagFCM.h"
182 #include "AliMagFDM.h"
184 #include "TRandom3.h"
186 #include "AliGenerator.h"
187 #include "AliLegoGenerator.h"
189 #include "AliDetector.h"
193 static AliHeader *gAliHeader;
197 //_____________________________________________________________________________
201 // Default constructor for AliRun
226 fPDGDB = 0; //Particle factory object!
228 fConfigFunction = "\0";
231 fTransParName = "\0";
232 fBaseFileName = ".\0";
234 fParticleMap = new TObjArray(10000);
237 //_____________________________________________________________________________
238 AliRun::AliRun(const char *name, const char *title)
242 // Constructor for the main processor.
243 // Creates the geometry
244 // Creates the list of Detectors.
245 // Creates the list of particles.
262 fConfigFunction = "Config();";
264 // Set random number generator
265 gRandom = fRandom = new TRandom3();
267 if (gSystem->Getenv("CONFIG_SEED")) {
268 gRandom->SetSeed((UInt_t)atoi(gSystem->Getenv("CONFIG_SEED")));
271 gROOT->GetListOfBrowsables()->Add(this,name);
273 // create the support list for the various Detectors
274 fModules = new TObjArray(77);
276 // Create the TNode geometry for the event display
278 BuildSimpleGeometry();
288 // Create the particle stack
289 fParticles = new TClonesArray("TParticle",1000);
293 // Create default mag field
298 // Prepare the tracking medium lists
299 fImedia = new TArrayI(1000);
300 for(i=0;i<1000;i++) (*fImedia)[i]=-99;
303 fPDGDB = TDatabasePDG::Instance(); //Particle factory object!
305 // Create HitLists list
306 fHitLists = new TList();
309 fBaseFileName = ".\0";
311 fParticleMap = new TObjArray(10000);
315 //_____________________________________________________________________________
319 // Default AliRun destructor
339 fParticles->Delete();
347 //_____________________________________________________________________________
348 void AliRun::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
351 // Add a hit to detector id
353 TObjArray &dets = *fModules;
354 if(dets[id]) ((AliModule*) dets[id])->AddHit(track,vol,hits);
357 //_____________________________________________________________________________
358 void AliRun::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
361 // Add digit to detector id
363 TObjArray &dets = *fModules;
364 if(dets[id]) ((AliModule*) dets[id])->AddDigit(tracks,digits);
367 //_____________________________________________________________________________
368 void AliRun::Browse(TBrowser *b)
371 // Called when the item "Run" is clicked on the left pane
372 // of the Root browser.
373 // It displays the Root Trees and all detectors.
375 if (fTreeK) b->Add(fTreeK,fTreeK->GetName());
376 if (fTreeH) b->Add(fTreeH,fTreeH->GetName());
377 if (fTreeD) b->Add(fTreeD,fTreeD->GetName());
378 if (fTreeE) b->Add(fTreeE,fTreeE->GetName());
379 if (fTreeR) b->Add(fTreeR,fTreeR->GetName());
380 if (fTreeS) b->Add(fTreeS,fTreeS->GetName());
382 TIter next(fModules);
384 while((detector = (AliModule*)next())) {
385 b->Add(detector,detector->GetName());
387 b->Add(fMCQA,"AliMCQA");
390 //_____________________________________________________________________________
394 // Initialize Alice geometry
399 //_____________________________________________________________________________
400 void AliRun::BuildSimpleGeometry()
403 // Create a simple TNode geometry used by Root display engine
405 // Initialise geometry
407 fGeometry = new TGeometry("AliceGeom","Galice Geometry for Hits");
408 new TMaterial("void","Vacuum",0,0,0); //Everything is void
409 TBRIK *brik = new TBRIK("S_alice","alice volume","void",2000,2000,3000);
410 brik->SetVisibility(0);
411 new TNode("alice","alice","S_alice");
414 //_____________________________________________________________________________
415 void AliRun::CleanDetectors()
418 // Clean Detectors at the end of event
420 TIter next(fModules);
422 while((detector = (AliModule*)next())) {
423 detector->FinishEvent();
427 //_____________________________________________________________________________
428 void AliRun::CleanParents()
431 // Clean Particles stack.
432 // Set parent/daughter relations
434 TObjArray &particles = *fParticleMap;
437 for(i=0; i<fHgwmk+1; i++) {
438 part = (TParticle *)particles.At(i);
439 if(part) if(!part->TestBit(kDaughtersBit)) {
440 part->SetFirstDaughter(-1);
441 part->SetLastDaughter(-1);
446 //_____________________________________________________________________________
447 Int_t AliRun::DistancetoPrimitive(Int_t, Int_t)
450 // Return the distance from the mouse to the AliRun object
456 //_____________________________________________________________________________
457 void AliRun::DumpPart (Int_t i) const
460 // Dumps particle i in the stack
462 ((TParticle*) (*fParticleMap)[i])->Print();
465 //_____________________________________________________________________________
466 void AliRun::DumpPStack () const
469 // Dumps the particle stack
471 TObjArray &particles = *fParticleMap;
473 "\n\n=======================================================================\n");
474 for (Int_t i=0;i<fNtrack;i++)
476 printf("-> %d ",i); ((TParticle*) particles[i])->Print();
477 printf("--------------------------------------------------------------\n");
480 "\n=======================================================================\n\n");
483 void AliRun::SetField(AliMagF* magField)
485 // Set Magnetic Field Map
490 //_____________________________________________________________________________
491 void AliRun::SetField(Int_t type, Int_t version, Float_t scale,
492 Float_t maxField, char* filename)
495 // Set magnetic field parameters
496 // type Magnetic field transport flag 0=no field, 2=helix, 3=Runge Kutta
497 // version Magnetic field map version (only 1 active now)
498 // scale Scale factor for the magnetic field
499 // maxField Maximum value for the magnetic field
502 // --- Sanity check on mag field flags
503 if(fField) delete fField;
505 fField = new AliMagFC("Map1"," ",type,scale,maxField);
506 } else if(version<=2) {
507 fField = new AliMagFCM("Map2-3",filename,type,scale,maxField);
509 } else if(version==3) {
510 fField = new AliMagFDM("Map4",filename,type,scale,maxField);
513 Warning("SetField","Invalid map %d\n",version);
517 //_____________________________________________________________________________
518 void AliRun::PreTrack()
520 TObjArray &dets = *fModules;
523 for(Int_t i=0; i<=fNdets; i++)
524 if((module = (AliModule*)dets[i]))
530 //_____________________________________________________________________________
531 void AliRun::PostTrack()
533 TObjArray &dets = *fModules;
536 for(Int_t i=0; i<=fNdets; i++)
537 if((module = (AliModule*)dets[i]))
541 //_____________________________________________________________________________
542 void AliRun::FinishPrimary()
545 // Called at the end of each primary track
548 // static Int_t count=0;
549 // const Int_t times=10;
550 // This primary is finished, purify stack
553 TIter next(fModules);
555 while((detector = (AliModule*)next())) {
556 detector->FinishPrimary();
559 // Write out hits if any
560 if (gAlice->TreeH()) {
561 gAlice->TreeH()->Fill();
568 // if(++count%times==1) gObjectTable->Print();
571 //_____________________________________________________________________________
572 void AliRun::FinishEvent()
575 // Called at the end of the event.
579 if(fLego) fLego->FinishEvent();
581 //Update the energy deposit tables
583 for(i=0;i<fEventEnergy.GetSize();i++) {
584 fSummEnergy[i]+=fEventEnergy[i];
585 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
587 fEventEnergy.Reset();
589 // Clean detector information
592 // Write out the kinematics
595 if(fTreeK->GetEntries() ==0) {
596 // set the fParticleFileMap size for the first time
597 fParticleFileMap.Set(fHgwmk+1);
600 Bool_t allFilled = kFALSE;
602 for(i=0; i<fHgwmk+1; ++i) if((part=fParticleMap->At(i))) {
603 fParticleBuffer = (TParticle*) part;
604 fParticleFileMap[i]= (Int_t) fTreeK->GetEntries();
606 (*fParticleMap)[i]=0;
608 // When all primaries were filled no particle!=0
609 // should be left => to be removed later.
610 if (allFilled) printf("Why != 0 part # %d?\n",i);
613 // // printf("Why = 0 part # %d?\n",i); => We know.
615 // we don't break now in order to be sure there is no
616 // particle !=0 left.
617 // To be removed later and replaced with break.
618 if(!allFilled) allFilled = kTRUE;
622 // Set number of tracks to event header
623 fHeader.SetNtrack(fNtrack);
625 // Write out the digits
636 // Write out reconstructed clusters
641 // Write out the event Header information
642 if (fTreeE) fTreeE->Fill();
647 // Write Tree headers
648 if (fTreeK) fTreeK->Write(0,TObject::kOverwrite);
649 if (fTreeH) fTreeH->Write(0,TObject::kOverwrite);
650 if (fTreeD) fTreeD->Write(0,TObject::kOverwrite);
651 if (fTreeR) fTreeR->Write(0,TObject::kOverwrite);
652 if (fTreeS) fTreeS->Write(0,TObject::kOverwrite);
657 //_____________________________________________________________________________
658 void AliRun::FinishRun()
661 // Called at the end of the run.
665 if(fLego) fLego->FinishRun();
667 // Clean detector information
668 TIter next(fModules);
670 while((detector = (AliModule*)next())) {
671 detector->FinishRun();
674 //Output energy summary tables
677 TFile *file = fTreeE->GetCurrentFile();
681 fTreeE->Write(0,TObject::kOverwrite);
683 // Write AliRun info and all detectors parameters
686 // Clean tree information
688 delete fTreeK; fTreeK = 0;
691 delete fTreeH; fTreeH = 0;
694 delete fTreeD; fTreeD = 0;
697 delete fTreeR; fTreeR = 0;
700 delete fTreeE; fTreeE = 0;
707 //_____________________________________________________________________________
708 void AliRun::FlagTrack(Int_t track)
711 // Flags a track and all its family tree to be kept
718 particle=(TParticle*)fParticleMap->At(curr);
720 // If the particle is flagged the three from here upward is saved already
721 if(particle->TestBit(kKeepBit)) return;
723 // Save this particle
724 particle->SetBit(kKeepBit);
726 // Move to father if any
727 if((curr=particle->GetFirstMother())==-1) return;
731 //_____________________________________________________________________________
732 void AliRun::EnergySummary()
735 // Print summary of deposited energy
741 Int_t kn, i, left, j, id;
742 const Float_t kzero=0;
743 Int_t ievent=fHeader.GetEvent()+1;
745 // Energy loss information
747 printf("***************** Energy Loss Information per event (GEV) *****************\n");
748 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
751 fEventEnergy[ndep]=kn;
756 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
759 fSummEnergy[ndep]=ed;
760 fSum2Energy[ndep]=TMath::Min((Float_t) 99.,TMath::Max(ed2,kzero));
765 for(kn=0;kn<(ndep-1)/3+1;kn++) {
767 for(i=0;i<(3<left?3:left);i++) {
769 id=Int_t (fEventEnergy[j]+0.1);
770 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
775 // Relative energy loss in different detectors
776 printf("******************** Relative Energy Loss per event ********************\n");
777 printf("Total energy loss per event %10.3f GeV\n",edtot);
778 for(kn=0;kn<(ndep-1)/5+1;kn++) {
780 for(i=0;i<(5<left?5:left);i++) {
782 id=Int_t (fEventEnergy[j]+0.1);
783 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
787 for(kn=0;kn<75;kn++) printf("*");
791 // Reset the TArray's
792 // fEventEnergy.Set(0);
793 // fSummEnergy.Set(0);
794 // fSum2Energy.Set(0);
797 //_____________________________________________________________________________
798 AliModule *AliRun::GetModule(const char *name) const
801 // Return pointer to detector from name
803 return (AliModule*)fModules->FindObject(name);
806 //_____________________________________________________________________________
807 AliDetector *AliRun::GetDetector(const char *name) const
810 // Return pointer to detector from name
812 return (AliDetector*)fModules->FindObject(name);
815 //_____________________________________________________________________________
816 Int_t AliRun::GetModuleID(const char *name) const
819 // Return galice internal detector identifier from name
822 TObject *mod=fModules->FindObject(name);
823 if(mod) i=fModules->IndexOf(mod);
827 //_____________________________________________________________________________
828 Int_t AliRun::GetEvent(Int_t event)
831 // Connect the Trees Kinematics and Hits for event # event
832 // Set branch addresses
835 // Reset existing structures
841 // Delete Trees already connected
842 if (fTreeK) delete fTreeK;
843 if (fTreeH) delete fTreeH;
844 if (fTreeD) delete fTreeD;
845 if (fTreeR) delete fTreeR;
846 if (fTreeS) delete fTreeS;
848 // Get header from file
849 if(fTreeE) fTreeE->GetEntry(event);
850 else Error("GetEvent","Cannot file Header Tree\n");
851 TFile *file = fTreeE->GetCurrentFile();
855 // Get Kine Tree from file
857 sprintf(treeName,"TreeK%d",event);
858 fTreeK = (TTree*)gDirectory->Get(treeName);
859 if (fTreeK) fTreeK->SetBranchAddress("Particles", &fParticleBuffer);
860 else Error("GetEvent","cannot find Kine Tree for event:%d\n",event);
861 // Create the particle stack
862 if(!fParticles) fParticles = new TClonesArray("TParticle",1000);
863 // Build the pointer list
865 fParticleMap->Clear();
866 fParticleMap->Expand(fTreeK->GetEntries());
868 fParticleMap = new TObjArray(fTreeK->GetEntries());
872 // Get Hits Tree header from file
873 sprintf(treeName,"TreeH%d",event);
874 fTreeH = (TTree*)gDirectory->Get(treeName);
876 Error("GetEvent","cannot find Hits Tree for event:%d\n",event);
881 // Get Digits Tree header from file
882 sprintf(treeName,"TreeD%d",event);
883 fTreeD = (TTree*)gDirectory->Get(treeName);
885 // Warning("GetEvent","cannot find Digits Tree for event:%d\n",event);
890 // Get SDigits Tree header from file
891 sprintf(treeName,"TreeS%d",event);
892 fTreeS = (TTree*)gDirectory->Get(treeName);
894 // Warning("GetEvent","cannot find SDigits Tree for event:%d\n",event);
899 // Get Reconstruct Tree header from file
900 sprintf(treeName,"TreeR%d",event);
901 fTreeR = (TTree*)gDirectory->Get(treeName);
903 // printf("WARNING: cannot find Reconstructed Tree for event:%d\n",event);
908 // Set Trees branch addresses
909 TIter next(fModules);
911 while((detector = (AliModule*)next())) {
912 detector->SetTreeAddress();
915 fNtrack = Int_t (fTreeK->GetEntries());
919 //_____________________________________________________________________________
920 TGeometry *AliRun::GetGeometry()
923 // Import Alice geometry from current file
924 // Return pointer to geometry object
926 if (!fGeometry) fGeometry = (TGeometry*)gDirectory->Get("AliceGeom");
928 // Unlink and relink nodes in detectors
929 // This is bad and there must be a better way...
932 TIter next(fModules);
934 while((detector = (AliModule*)next())) {
935 TList *dnodes=detector->Nodes();
938 for ( j=0; j<dnodes->GetSize(); j++) {
939 node = (TNode*) dnodes->At(j);
940 node1 = fGeometry->GetNode(node->GetName());
941 dnodes->Remove(node);
942 dnodes->AddAt(node1,j);
948 //_____________________________________________________________________________
949 void AliRun::GetNextTrack(Int_t &mtrack, Int_t &ipart, Float_t *pmom,
950 Float_t &e, Float_t *vpos, Float_t *polar,
954 // Return next track from stack of particles
959 for(Int_t i=fNtrack-1; i>=0; i--) {
960 track=(TParticle*) fParticleMap->At(i);
961 if(track) if(!track->TestBit(kDoneBit)) {
963 // The track exists and has not yet been processed
965 ipart=track->GetPdgCode();
973 track->GetPolarisation(pol);
978 track->SetBit(kDoneBit);
984 // stop and start timer when we start a primary track
985 Int_t nprimaries = fHeader.GetNprimary();
986 if (fCurrent >= nprimaries) return;
987 if (fCurrent < nprimaries-1) {
989 track=(TParticle*) fParticleMap->At(fCurrent+1);
990 // track->SetProcessTime(fTimer.CpuTime());
995 //_____________________________________________________________________________
996 Int_t AliRun::GetPrimary(Int_t track) const
999 // return number of primary that has generated track
1001 int current, parent;
1007 part = (TParticle *)fParticleMap->At(current);
1008 parent=part->GetFirstMother();
1009 if(parent<0) return current;
1013 //_____________________________________________________________________________
1014 void AliRun::InitMC(const char *setup)
1017 // Initialize the Alice setup
1021 Warning("Init","Cannot initialise AliRun twice!\n");
1025 gROOT->LoadMacro(setup);
1026 gInterpreter->ProcessLine(fConfigFunction.Data());
1029 gMC->DefineParticles(); //Create standard MC particles
1031 TObject *objfirst, *objlast;
1033 fNdets = fModules->GetLast()+1;
1036 //=================Create Materials and geometry
1039 // Added also after in case of interactive initialisation of modules
1040 fNdets = fModules->GetLast()+1;
1042 TIter next(fModules);
1043 AliModule *detector;
1044 while((detector = (AliModule*)next())) {
1045 detector->SetTreeAddress();
1046 objlast = gDirectory->GetList()->Last();
1048 // Add Detector histograms in Detector list of histograms
1049 if (objlast) objfirst = gDirectory->GetList()->After(objlast);
1050 else objfirst = gDirectory->GetList()->First();
1052 detector->Histograms()->Add(objfirst);
1053 objfirst = gDirectory->GetList()->After(objfirst);
1056 ReadTransPar(); //Read the cuts for all materials
1058 MediaTable(); //Build the special IMEDIA table
1060 //Initialise geometry deposition table
1061 fEventEnergy.Set(gMC->NofVolumes()+1);
1062 fSummEnergy.Set(gMC->NofVolumes()+1);
1063 fSum2Energy.Set(gMC->NofVolumes()+1);
1065 //Compute cross-sections
1066 gMC->BuildPhysics();
1068 //Write Geometry object to current file.
1073 fMCQA = new AliMCQA(fNdets);
1076 // Save stuff at the beginning of the file to avoid file corruption
1080 //_____________________________________________________________________________
1081 void AliRun::MediaTable()
1084 // Built media table to get from the media number to
1087 Int_t kz, nz, idt, lz, i, k, ind;
1089 TObjArray &dets = *gAlice->Detectors();
1092 // For all detectors
1093 for (kz=0;kz<fNdets;kz++) {
1094 // If detector is defined
1095 if((det=(AliModule*) dets[kz])) {
1096 TArrayI &idtmed = *(det->GetIdtmed());
1097 for(nz=0;nz<100;nz++) {
1098 // Find max and min material number
1099 if((idt=idtmed[nz])) {
1100 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
1101 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
1104 if(det->LoMedium() > det->HiMedium()) {
1105 det->LoMedium() = 0;
1106 det->HiMedium() = 0;
1108 if(det->HiMedium() > fImedia->GetSize()) {
1109 Error("MediaTable","Increase fImedia from %d to %d",
1110 fImedia->GetSize(),det->HiMedium());
1113 // Tag all materials in rage as belonging to detector kz
1114 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
1121 // Print summary table
1122 printf(" Traking media ranges:\n");
1123 for(i=0;i<(fNdets-1)/6+1;i++) {
1124 for(k=0;k< (6<fNdets-i*6?6:fNdets-i*6);k++) {
1126 det=(AliModule*)dets[ind];
1128 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
1131 printf(" %6s: %3d -> %3d;","NULL",0,0);
1137 //____________________________________________________________________________
1138 void AliRun::SetGenerator(AliGenerator *generator)
1141 // Load the event generator
1143 if(!fGenerator) fGenerator = generator;
1146 //____________________________________________________________________________
1147 void AliRun::ResetGenerator(AliGenerator *generator)
1150 // Load the event generator
1154 Warning("ResetGenerator","Replacing generator %s with %s\n",
1155 fGenerator->GetName(),generator->GetName());
1157 Warning("ResetGenerator","Replacing generator %s with NULL\n",
1158 fGenerator->GetName());
1159 fGenerator = generator;
1162 //____________________________________________________________________________
1163 void AliRun::SetTransPar(char *filename)
1165 fTransParName = filename;
1168 //____________________________________________________________________________
1169 void AliRun::SetBaseFile(char *filename)
1171 fBaseFileName = filename;
1174 //____________________________________________________________________________
1175 void AliRun::ReadTransPar()
1178 // Read filename to set the transport parameters
1182 const Int_t kncuts=10;
1183 const Int_t knflags=11;
1184 const Int_t knpars=kncuts+knflags;
1185 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
1186 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
1187 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
1188 "MULS","PAIR","PHOT","RAYL"};
1192 Float_t cut[kncuts];
1193 Int_t flag[knflags];
1194 Int_t i, itmed, iret, ktmed, kz;
1197 // See whether the file is there
1198 filtmp=gSystem->ExpandPathName(fTransParName.Data());
1199 lun=fopen(filtmp,"r");
1202 Warning("ReadTransPar","File %s does not exist!\n",fTransParName.Data());
1206 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
1207 printf(" *%59s\n","*");
1208 printf(" * Please check carefully what you are doing!%10s\n","*");
1209 printf(" *%59s\n","*");
1212 // Initialise cuts and flags
1213 for(i=0;i<kncuts;i++) cut[i]=-99;
1214 for(i=0;i<knflags;i++) flag[i]=-99;
1216 for(i=0;i<256;i++) line[i]='\0';
1217 // Read up to the end of line excluded
1218 iret=fscanf(lun,"%[^\n]",line);
1222 printf(" *%59s\n","*");
1223 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
1226 // Read the end of line
1229 if(line[0]=='*') continue;
1231 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",
1232 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
1233 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
1234 &flag[8],&flag[9],&flag[10]);
1238 Warning("ReadTransPar","Error reading file %s\n",fTransParName.Data());
1241 // Check that the module exist
1242 AliModule *mod = GetModule(detName);
1244 // Get the array of media numbers
1245 TArrayI &idtmed = *mod->GetIdtmed();
1246 // Check that the tracking medium code is valid
1247 if(0<=itmed && itmed < 100) {
1248 ktmed=idtmed[itmed];
1250 Warning("ReadTransPar","Invalid tracking medium code %d for %s\n",itmed,mod->GetName());
1253 // Set energy thresholds
1254 for(kz=0;kz<kncuts;kz++) {
1256 printf(" * %-6s set to %10.3E for tracking medium code %4d for %s\n",
1257 kpars[kz],cut[kz],itmed,mod->GetName());
1258 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
1261 // Set transport mechanisms
1262 for(kz=0;kz<knflags;kz++) {
1264 printf(" * %-6s set to %10d for tracking medium code %4d for %s\n",
1265 kpars[kncuts+kz],flag[kz],itmed,mod->GetName());
1266 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
1270 Warning("ReadTransPar","Invalid medium code %d *\n",itmed);
1274 Warning("ReadTransPar","Module %s not present\n",detName);
1280 //_____________________________________________________________________________
1281 void AliRun::MakeBranchInTree(TTree *tree, const char* name, void* address, Int_t size, char *file)
1284 printf("* MakeBranch * Making Branch %s \n",name);
1286 TBranch *branch = tree->Branch(name,address,size);
1289 char * outFile = new char[strlen(gAlice->GetBaseFile())+strlen(file)+2];
1290 sprintf(outFile,"%s/%s",gAlice->GetBaseFile(),file);
1291 TDirectory *cwd = gDirectory;
1292 branch->SetFile(outFile);
1293 TIter next( branch->GetListOfBranches());
1294 while ((branch=(TBranch*)next())) {
1295 branch->SetFile(outFile);
1298 printf("* MakeBranch * Diverting Branch %s to file %s\n",name,file);
1304 //_____________________________________________________________________________
1305 void AliRun::MakeBranchInTree(TTree *tree, const char* name, const char *classname, void* address, Int_t size, Int_t splitlevel, char *file)
1307 TDirectory *cwd = gDirectory;
1308 TBranch *branch = tree->Branch(name,classname,address,size,splitlevel);
1310 printf("* MakeBranch * Making Branch %s \n",name);
1312 char * outFile = new char[strlen(gAlice->GetBaseFile())+strlen(file)+2];
1313 sprintf(outFile,"%s/%s",gAlice->GetBaseFile(),file);
1314 branch->SetFile(outFile);
1315 TIter next( branch->GetListOfBranches());
1316 while ((branch=(TBranch*)next())) {
1317 branch->SetFile(outFile);
1320 printf("* MakeBranch * Diverting Branch %s to file %s\n",name,file);
1325 //_____________________________________________________________________________
1326 void AliRun::MakeTree(Option_t *option, char *file)
1329 // Create the ROOT trees
1330 // Loop on all detectors to create the Root branch (if any)
1336 char *oK = strstr(option,"K");
1337 char *oH = strstr(option,"H");
1338 char *oE = strstr(option,"E");
1339 char *oD = strstr(option,"D");
1340 char *oR = strstr(option,"R");
1341 char *oS = strstr(option,"S");
1344 if (oK && !fTreeK) {
1345 sprintf(hname,"TreeK%d",fEvent);
1346 fTreeK = new TTree(hname,"Kinematics");
1347 // Create a branch for particles
1348 MakeBranchInTree(fTreeK,
1349 "Particles", "TParticle", &fParticleBuffer, 4000, 1, file) ;
1352 if (oH && !fTreeH) {
1353 sprintf(hname,"TreeH%d",fEvent);
1354 fTreeH = new TTree(hname,"Hits");
1355 fTreeH->SetAutoSave(1000000000); //no autosave
1358 if (oD && !fTreeD) {
1359 sprintf(hname,"TreeD%d",fEvent);
1360 fTreeD = new TTree(hname,"Digits");
1363 if (oS && !fTreeS) {
1364 sprintf(hname,"TreeS%d",fEvent);
1365 fTreeS = new TTree(hname,"SDigits");
1368 if (oR && !fTreeR) {
1369 sprintf(hname,"TreeR%d",fEvent);
1370 fTreeR = new TTree(hname,"Reconstruction");
1373 if (oE && !fTreeE) {
1374 fTreeE = new TTree("TE","Header");
1375 // Create a branch for Header
1376 MakeBranchInTree(fTreeE,
1377 "Header", "AliHeader", &gAliHeader, 4000, 1, file) ;
1382 // Create a branch for hits/digits for each detector
1383 // Each branch is a TClonesArray. Each data member of the Hits classes
1384 // will be in turn a subbranch of the detector master branch
1385 TIter next(fModules);
1386 AliModule *detector;
1387 while((detector = (AliModule*)next())) {
1388 if (oH || oR) detector->MakeBranch(option,file);
1392 //_____________________________________________________________________________
1393 Int_t AliRun::PurifyKine(Int_t lastSavedTrack, Int_t nofTracks)
1396 // PurifyKine with external parameters
1398 fHgwmk = lastSavedTrack;
1399 fNtrack = nofTracks;
1404 //_____________________________________________________________________________
1405 TParticle* AliRun::Particle(Int_t i)
1407 if(!(*fParticleMap)[i]) {
1408 Int_t nentries = fParticles->GetEntries();
1410 // algorithmic way of getting entry index
1411 // (primary particles are filled after secondaries)
1413 if (i<fHeader.GetNprimary())
1414 entry = i+fHeader.GetNsecondary();
1416 entry = i-fHeader.GetNprimary();
1418 // only check the algorithmic way and give
1419 // the fatal error if it is wrong
1420 if (entry != fParticleFileMap[i]) {
1422 "!!!! The algorithmic way is WRONG: !!!\n entry: %d map: %d",
1423 entry, fParticleFileMap[i]);
1426 fTreeK->GetEntry(fParticleFileMap[i]);
1427 new ((*fParticles)[nentries]) TParticle(*fParticleBuffer);
1428 fParticleMap->AddAt((*fParticles)[nentries],i);
1430 return (TParticle *) (*fParticleMap)[i];
1433 //_____________________________________________________________________________
1434 void AliRun::PurifyKine()
1437 // Compress kinematic tree keeping only flagged particles
1438 // and renaming the particle id's in all the hits
1440 // TClonesArray &particles = *fParticles;
1441 TObjArray &particles = *fParticleMap;
1442 int nkeep=fHgwmk+1, parent, i;
1443 TParticle *part, *father;
1444 TArrayI map(particles.GetLast()+1);
1446 // Save in Header total number of tracks before compression
1447 fHeader.SetNtrack(fHeader.GetNtrack()+fNtrack-fHgwmk);
1449 // If no tracks generated return now
1450 if(fHgwmk+1 == fNtrack) return;
1452 Int_t toshrink = fNtrack-fHgwmk-1;
1454 // First pass, invalid Daughter information
1455 for(i=0; i<fNtrack; i++) {
1456 // Preset map, to be removed later
1457 if(i<=fHgwmk) map[i]=i ;
1460 // particles.UncheckedAt(i)->ResetBit(kDaughtersBit);
1461 if((part=(TParticle*) particles.At(i))) part->ResetBit(kDaughtersBit);
1464 // Invalid daughter information for the parent of the first particle
1465 // generated. This may or may not be the current primary according to
1466 // whether decays have been recorded among the primaries
1467 part = (TParticle *)particles.At(fHgwmk+1);
1468 particles.At(part->GetFirstMother())->ResetBit(kDaughtersBit);
1469 // Second pass, build map between old and new numbering
1470 for(i=fHgwmk+1; i<fNtrack; i++) {
1471 if(particles.At(i)->TestBit(kKeepBit)) {
1473 // This particle has to be kept
1475 // If old and new are different, have to move the pointer
1476 if(i!=nkeep) particles[nkeep]=particles.At(i);
1477 part = (TParticle*) particles.At(nkeep);
1479 // as the parent is always *before*, it must be already
1480 // in place. This is what we are checking anyway!
1481 if((parent=part->GetFirstMother())>fHgwmk)
1482 if(map[parent]==-99) Fatal("PurifyKine","map[%d] = -99!\n",parent);
1483 else part->SetFirstMother(map[parent]);
1489 // Fix daughters information
1490 for (i=fHgwmk+1; i<nkeep; i++) {
1491 part = (TParticle *)particles.At(i);
1492 parent = part->GetFirstMother();
1494 father = (TParticle *)particles.At(parent);
1495 if(father->TestBit(kDaughtersBit)) {
1497 if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
1498 if(i>father->GetLastDaughter()) father->SetLastDaughter(i);
1500 // Initialise daughters info for first pass
1501 father->SetFirstDaughter(i);
1502 father->SetLastDaughter(i);
1503 father->SetBit(kDaughtersBit);
1508 // Now loop on all registered hit lists
1509 TIter next(fHitLists);
1510 TCollection *hitList;
1511 while((hitList = (TCollection*)next())) {
1512 TIter nexthit(hitList);
1514 while((hit = (AliHit*)nexthit())) {
1515 hit->SetTrack(map[hit->GetTrack()]);
1520 // This for detectors which have a special mapping mechanism
1521 // for hits, such as TPC and TRD
1524 TIter nextmod(fModules);
1525 AliModule *detector;
1526 while((detector = (AliModule*)nextmod())) {
1527 detector->RemapTrackHitIDs(map.GetArray());
1530 // Now the output bit, from fHgwmk to nkeep we write everything and we erase
1531 if(nkeep>fParticleFileMap.GetSize()) fParticleFileMap.Set(Int_t (nkeep*1.5));
1534 for (i=fHgwmk+1; i<nkeep; ++i) {
1535 fParticleBuffer = (TParticle*) particles.At(i);
1536 fParticleFileMap[i]=(Int_t) fTreeK->GetEntries();
1541 for (i=nkeep; i<fNtrack; ++i) particles[i]=0;
1543 fLoadPoint-=toshrink;
1544 for(i=fLoadPoint; i<fLoadPoint+toshrink; ++i) fParticles->RemoveAt(i);
1551 //_____________________________________________________________________________
1552 void AliRun::BeginEvent()
1555 // Reset all Detectors & kinematics & trees
1562 fLego->BeginEvent();
1572 // Initialise event header
1573 fHeader.Reset(fRun,fEvent);
1577 sprintf(hname,"TreeK%d",fEvent);
1578 fTreeK->SetName(hname);
1582 sprintf(hname,"TreeH%d",fEvent);
1583 fTreeH->SetName(hname);
1587 sprintf(hname,"TreeD%d",fEvent);
1588 fTreeD->SetName(hname);
1592 sprintf(hname,"TreeS%d",fEvent);
1593 fTreeS->SetName(hname);
1597 sprintf(hname,"TreeR%d",fEvent);
1598 fTreeR->SetName(hname);
1601 //_____________________________________________________________________________
1602 void AliRun::ResetDigits()
1605 // Reset all Detectors digits
1607 TIter next(fModules);
1608 AliModule *detector;
1609 while((detector = (AliModule*)next())) {
1610 detector->ResetDigits();
1614 //_____________________________________________________________________________
1615 void AliRun::ResetSDigits()
1618 // Reset all Detectors digits
1620 TIter next(fModules);
1621 AliModule *detector;
1622 while((detector = (AliModule*)next())) {
1623 detector->ResetSDigits();
1627 //_____________________________________________________________________________
1628 void AliRun::ResetHits()
1631 // Reset all Detectors hits
1633 TIter next(fModules);
1634 AliModule *detector;
1635 while((detector = (AliModule*)next())) {
1636 detector->ResetHits();
1640 //_____________________________________________________________________________
1641 void AliRun::ResetPoints()
1644 // Reset all Detectors points
1646 TIter next(fModules);
1647 AliModule *detector;
1648 while((detector = (AliModule*)next())) {
1649 detector->ResetPoints();
1653 //_____________________________________________________________________________
1654 void AliRun::RunMC(Int_t nevent, const char *setup)
1657 // Main function to be called to process a galice run
1659 // Root > gAlice.Run();
1660 // a positive number of events will cause the finish routine
1664 // check if initialisation has been done
1665 if (!fInitDone) InitMC(setup);
1667 // Create the Root Tree with one branch per detector
1669 if (gSystem->Getenv("CONFIG_SPLIT_FILE")) {
1671 MakeTree("K","Kine.root");
1672 MakeTree("H","Hits.root");
1673 MakeTree("R","Reco.root");
1678 gMC->ProcessRun(nevent);
1680 // End of this run, close files
1681 if(nevent>0) FinishRun();
1684 //_____________________________________________________________________________
1686 void AliRun::Hits2Digits(const char *selected)
1688 Hits2SDigits(selected);
1689 SDigits2Digits(selected);
1692 //_____________________________________________________________________________
1694 void AliRun::Hits2SDigits(const char *selected)
1697 // Main function to be called to convert hits to digits.
1699 gAlice->GetEvent(0);
1701 TObjArray *detectors = gAlice->Detectors();
1703 TIter next(detectors);
1705 AliDetector *detector;
1707 TDirectory *cwd = gDirectory;
1711 while((detector = (AliDetector*)next())) {
1713 if (strcmp(detector->GetName(),selected)) continue;
1715 if (detector->IsActive()){
1716 if (gSystem->Getenv("CONFIG_SPLIT_FILE")) {
1718 cout << "Processing " << detector->GetName() << "..." << endl;
1719 char * outFile = new char[strlen (detector->GetName())+18];
1720 sprintf(outFile,"SDigits.%s.root",detector->GetName());
1721 detector->MakeBranch("S",outFile);
1724 detector->MakeBranch("S");
1727 detector->Hits2SDigits();
1732 //_____________________________________________________________________________
1734 void AliRun::SDigits2Digits(const char *selected)
1737 // Main function to be called to convert hits to digits.
1739 gAlice->GetEvent(0);
1741 TObjArray *detectors = gAlice->Detectors();
1743 TIter next(detectors);
1745 AliDetector *detector;
1747 TDirectory *cwd = gDirectory;
1751 while((detector = (AliDetector*)next())) {
1753 if (strcmp(detector->GetName(),selected)) continue;
1755 if (detector->IsActive()){
1756 if (gSystem->Getenv("CONFIG_SPLIT_FILE")) {
1758 cout << "Processing " << detector->GetName() << "..." << endl;
1759 char * outFile = new char[strlen (detector->GetName())+16];
1760 sprintf(outFile,"Digits.%s.root",detector->GetName());
1761 detector->MakeBranch("D",outFile);
1764 detector->MakeBranch("D");
1767 detector->SDigits2Digits();
1772 //_____________________________________________________________________________
1773 void AliRun::RunLego(const char *setup, Int_t nc1, Float_t c1min,
1774 Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
1775 Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener)
1778 // Generates lego plots of:
1779 // - radiation length map phi vs theta
1780 // - radiation length map phi vs eta
1781 // - interaction length map
1782 // - g/cm2 length map
1784 // ntheta bins in theta, eta
1785 // themin minimum angle in theta (degrees)
1786 // themax maximum angle in theta (degrees)
1788 // phimin minimum angle in phi (degrees)
1789 // phimax maximum angle in phi (degrees)
1790 // rmin minimum radius
1791 // rmax maximum radius
1794 // The number of events generated = ntheta*nphi
1795 // run input parameters in macro setup (default="Config.C")
1797 // Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
1800 <img src="picts/AliRunLego1.gif">
1805 <img src="picts/AliRunLego2.gif">
1810 <img src="picts/AliRunLego3.gif">
1815 // check if initialisation has been done
1816 if (!fInitDone) InitMC(setup);
1817 //Save current generator
1818 AliGenerator *gen=Generator();
1820 // Set new generator
1821 if (!gener) gener = new AliLegoGenerator();
1822 ResetGenerator(gener);
1824 // Configure Generator
1825 gener->SetRadiusRange(rmin, rmax);
1826 gener->SetZMax(zmax);
1827 gener->SetCoor1Range(nc1, c1min, c1max);
1828 gener->SetCoor2Range(nc2, c2min, c2max);
1831 //Create Lego object
1832 fLego = new AliLego("lego",gener);
1834 //Prepare MC for Lego Run
1839 gMC->ProcessRun(nc1*nc2+1);
1841 // Create only the Root event Tree
1844 // End of this run, close files
1846 // Restore current generator
1847 ResetGenerator(gen);
1848 // Delete Lego Object
1849 delete fLego; fLego=0;
1852 //_____________________________________________________________________________
1853 void AliRun::SetConfigFunction(const char * config)
1856 // Set the signature of the function contained in Config.C to configure
1859 fConfigFunction=config;
1862 //_____________________________________________________________________________
1863 void AliRun::SetCurrentTrack(Int_t track)
1866 // Set current track number
1871 //_____________________________________________________________________________
1872 void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
1873 Float_t *vpos, Float_t *polar, Float_t tof,
1874 AliMCProcess mech, Int_t &ntr, Float_t weight)
1877 // Load a track on the stack
1879 // done 0 if the track has to be transported
1881 // parent identifier of the parent track. -1 for a primary
1882 // pdg particle code
1883 // pmom momentum GeV/c
1885 // polar polarisation
1886 // tof time of flight in seconds
1887 // mecha production mechanism
1888 // ntr on output the number of the track stored
1890 TClonesArray &particles = *fParticles;
1891 TParticle *particle;
1893 const Int_t kfirstdaughter=-1;
1894 const Int_t klastdaughter=-1;
1896 // const Float_t tlife=0;
1899 // Here we get the static mass
1900 // For MC is ok, but a more sophisticated method could be necessary
1901 // if the calculated mass is required
1902 // also, this method is potentially dangerous if the mass
1903 // used in the MC is not the same of the PDG database
1905 mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
1906 Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
1907 pmom[1]*pmom[1]+pmom[2]*pmom[2]);
1909 //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",
1910 //pname,mass,e,fNtrack,pdg,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],kS,mecha);
1912 particle=new(particles[fLoadPoint++]) TParticle(pdg,kS,parent,-1,kfirstdaughter,
1913 klastdaughter,pmom[0],pmom[1],pmom[2],
1914 e,vpos[0],vpos[1],vpos[2],tof);
1915 particle->SetPolarisation(TVector3(polar[0],polar[1],polar[2]));
1916 particle->SetWeight(weight);
1917 particle->SetUniqueID(mech);
1918 if(!done) particle->SetBit(kDoneBit);
1919 // Declare that the daughter information is valid
1920 particle->SetBit(kDaughtersBit);
1921 // Add the particle to the stack
1922 fParticleMap->AddAtAndExpand(particle,fNtrack);
1925 particle=(TParticle*) fParticleMap->At(parent);
1926 particle->SetLastDaughter(fNtrack);
1927 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
1930 // This is a primary track. Set high water mark for this event
1933 // Set also number if primary tracks
1934 fHeader.SetNprimary(fHgwmk+1);
1935 fHeader.SetNtrack(fHgwmk+1);
1941 // Here we get the static mass
1942 // For MC is ok, but a more sophisticated method could be necessary
1943 // if the calculated mass is required
1944 // also, this method is potentially dangerous if the mass
1945 // used in the MC is not the same of the PDG database
1947 Float_t mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
1948 Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
1949 pmom[1]*pmom[1]+pmom[2]*pmom[2]);
1951 SetTrack(done, parent, pdg, pmom[0], pmom[1], pmom[2], e,
1952 vpos[0], vpos[1], vpos[2], tof, polar[0],polar[1],polar[2],
1957 //_____________________________________________________________________________
1958 void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg,
1959 Double_t px, Double_t py, Double_t pz, Double_t e,
1960 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
1961 Double_t polx, Double_t poly, Double_t polz,
1962 AliMCProcess mech, Int_t &ntr, Float_t weight)
1965 // Load a track on the stack
1967 // done 0 if the track has to be transported
1969 // parent identifier of the parent track. -1 for a primary
1970 // pdg particle code
1971 // kS generation status code
1972 // px, py, pz momentum GeV/c
1973 // vx, vy, vz position
1974 // polar polarisation
1975 // tof time of flight in seconds
1976 // mech production mechanism
1977 // ntr on output the number of the track stored
1979 // New method interface:
1980 // arguments were changed to be in correspondence with TParticle
1982 // Note: the energy is not calculated from the static mass but
1983 // it is passed by argument e.
1985 TClonesArray &particles = *fParticles;
1988 const Int_t kFirstDaughter=-1;
1989 const Int_t kLastDaughter=-1;
1992 = new(particles[fLoadPoint++]) TParticle(pdg, kS, parent, -1,
1993 kFirstDaughter, kLastDaughter,
1994 px, py, pz, e, vx, vy, vz, tof);
1996 particle->SetPolarisation(polx, poly, polz);
1997 particle->SetWeight(weight);
1998 particle->SetUniqueID(mech);
2000 if(!done) particle->SetBit(kDoneBit);
2002 // Declare that the daughter information is valid
2003 particle->SetBit(kDaughtersBit);
2004 // Add the particle to the stack
2005 fParticleMap->AddAtAndExpand(particle,fNtrack);
2008 particle=(TParticle*) fParticleMap->At(parent);
2009 particle->SetLastDaughter(fNtrack);
2010 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
2013 // This is a primary track. Set high water mark for this event
2016 // Set also number if primary tracks
2017 fHeader.SetNprimary(fHgwmk+1);
2018 fHeader.SetNtrack(fHgwmk+1);
2023 //_____________________________________________________________________________
2024 void AliRun::SetHighWaterMark(const Int_t nt)
2027 // Set high water mark for last track in event
2030 // Set also number if primary tracks
2031 fHeader.SetNprimary(fHgwmk+1);
2032 fHeader.SetNtrack(fHgwmk+1);
2035 //_____________________________________________________________________________
2036 void AliRun::KeepTrack(const Int_t track)
2039 // flags a track to be kept
2041 fParticleMap->At(track)->SetBit(kKeepBit);
2044 //_____________________________________________________________________________
2045 void AliRun::StepManager(Int_t id)
2048 // Called at every step during transport
2052 // --- If lego option, do it and leave
2054 fLego->StepManager();
2057 //Update energy deposition tables
2058 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
2060 //Call the appropriate stepping routine;
2061 AliModule *det = (AliModule*)fModules->At(id);
2063 fMCQA->StepManager(id);
2069 //_____________________________________________________________________________
2070 void AliRun::Streamer(TBuffer &R__b)
2072 // Stream an object of class AliRun.
2074 if (R__b.IsReading()) {
2075 if (!gAlice) gAlice = this;
2077 AliRun::Class()->ReadBuffer(R__b, this);
2079 gROOT->GetListOfBrowsables()->Add(this,"Run");
2081 fTreeE = (TTree*)gDirectory->Get("TE");
2082 if (fTreeE) fTreeE->SetBranchAddress("Header", &gAliHeader);
2083 else Error("Streamer","cannot find Header Tree\n");
2084 fTreeE->GetEntry(0);
2088 AliRun::Class()->WriteBuffer(R__b, this);