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.60 2001/03/21 18:22:30 hristov
19 fParticleFileMap fix (I.Hrivnacova)
21 Revision 1.59 2001/03/12 17:47:03 hristov
22 Changes needed on Sun with CC 5.0
24 Revision 1.58 2001/03/09 14:27:26 morsch
25 Fix for multiple events per file: inhibit decrease of size of fParticleFileMap.
27 Revision 1.57 2001/02/23 17:40:23 buncic
28 All trees needed for simulation created in RunMC(). TreeR and its branches
29 are now created in new RunReco() method.
31 Revision 1.56 2001/02/14 15:45:20 hristov
32 Algorithmic way of getting entry index in fParticleMap. Protection of fParticleFileMap (I.Hrivnacova)
34 Revision 1.55 2001/02/12 15:52:54 buncic
35 Removed OpenBaseFile().
37 Revision 1.54 2001/02/07 10:39:05 hristov
38 Remove default value for argument
40 Revision 1.53 2001/02/06 11:02:26 hristov
41 New SetTrack interface added, added check for unfilled particles in FinishEvent (I.Hrivnacova)
43 Revision 1.52 2001/02/05 16:22:25 buncic
44 Added TreeS to GetEvent().
46 Revision 1.51 2001/02/02 15:16:20 morsch
47 SetHighWaterMark method added to mark last particle in event.
49 Revision 1.50 2001/01/27 10:32:00 hristov
50 Leave the loop when primaries are filled (I.Hrivnacova)
52 Revision 1.49 2001/01/26 19:58:48 hristov
53 Major upgrade of AliRoot code
55 Revision 1.48 2001/01/17 10:50:50 hristov
56 Corrections to destructors
58 Revision 1.47 2000/12/18 10:44:01 morsch
59 Possibility to set field map by passing pointer to objet of type AliMagF via
62 gAlice->SetField(new AliMagFCM("Map2", "$(ALICE_ROOT)/data/field01.dat",2,1.,10.));
64 Revision 1.46 2000/12/14 19:29:27 fca
65 galice.cuts was not read any more
67 Revision 1.45 2000/11/30 07:12:49 alibrary
68 Introducing new Rndm and QA classes
70 Revision 1.44 2000/10/26 13:58:59 morsch
71 Add possibility to choose the lego generator (of type AliGeneratorLego or derived) when running
72 RunLego(). Default is the base class AliGeneratorLego.
74 Revision 1.43 2000/10/09 09:43:17 fca
75 Special remapping of hits for TPC and TRD. End-of-primary action introduced
77 Revision 1.42 2000/10/02 21:28:14 fca
78 Removal of useless dependecies via forward declarations
80 Revision 1.41 2000/07/13 16:19:09 fca
81 Mainly coding conventions + some small bug fixes
83 Revision 1.40 2000/07/12 08:56:25 fca
84 Coding convention correction and warning removal
86 Revision 1.39 2000/07/11 18:24:59 fca
87 Coding convention corrections + few minor bug fixes
89 Revision 1.38 2000/06/20 13:05:45 fca
90 Writing down the TREE headers before job starts
92 Revision 1.37 2000/06/09 20:05:11 morsch
93 Introduce possibility to chose magnetic field version 3: AliMagFDM + field02.dat
95 Revision 1.36 2000/06/08 14:03:58 hristov
96 Only one initializer for a default argument
98 Revision 1.35 2000/06/07 10:13:14 hristov
99 Delete only existent objects.
101 Revision 1.34 2000/05/18 10:45:38 fca
102 Delete Particle Factory properly
104 Revision 1.33 2000/05/16 13:10:40 fca
105 New method IsNewTrack and fix for a problem in Father-Daughter relations
107 Revision 1.32 2000/04/27 10:38:21 fca
108 Correct termination of Lego Run and introduce Lego getter in AliRun
110 Revision 1.31 2000/04/26 10:17:32 fca
111 Changes in Lego for G4 compatibility
113 Revision 1.30 2000/04/18 19:11:40 fca
114 Introduce variable Config.C function signature
116 Revision 1.29 2000/04/07 11:12:34 fca
117 G4 compatibility changes
119 Revision 1.28 2000/04/05 06:51:06 fca
120 Workaround for an HP compiler problem
122 Revision 1.27 2000/03/22 18:08:07 fca
123 Rationalisation of the virtual MC interfaces
125 Revision 1.26 2000/03/22 13:42:26 fca
126 SetGenerator does not replace an existing generator, ResetGenerator does
128 Revision 1.25 2000/02/23 16:25:22 fca
129 AliVMC and AliGeant3 classes introduced
130 ReadEuclid moved from AliRun to AliModule
132 Revision 1.24 2000/01/19 17:17:20 fca
133 Introducing a list of lists of hits -- more hits allowed for detector now
135 Revision 1.23 1999/12/03 11:14:31 fca
136 Fixing previous wrong checking
138 Revision 1.21 1999/11/25 10:40:08 fca
139 Fixing daughters information also in primary tracks
141 Revision 1.20 1999/10/04 18:08:49 fca
142 Adding protection against inconsistent Euclid files
144 Revision 1.19 1999/09/29 07:50:40 fca
145 Introduction of the Copyright and cvs Log
149 ///////////////////////////////////////////////////////////////////////////////
151 // Control class for Alice C++ //
152 // Only one single instance of this class exists. //
153 // The object is created in main program aliroot //
154 // and is pointed by the global gAlice. //
156 // -Supports the list of all Alice Detectors (fModules). //
157 // -Supports the list of particles (fParticles). //
158 // -Supports the Trees. //
159 // -Supports the geometry. //
160 // -Supports the event display. //
163 <img src="picts/AliRunClass.gif">
168 <img src="picts/alirun.gif">
172 ///////////////////////////////////////////////////////////////////////////////
177 #include <iostream.h>
185 #include <TObjectTable.h>
187 #include <TGeometry.h>
189 #include "TBrowser.h"
191 #include "TParticle.h"
193 #include "AliDisplay.h"
196 #include "AliMagFC.h"
197 #include "AliMagFCM.h"
198 #include "AliMagFDM.h"
200 #include "TRandom3.h"
202 #include "AliGenerator.h"
203 #include "AliLegoGenerator.h"
205 #include "AliDetector.h"
209 static AliHeader *gAliHeader;
213 //_____________________________________________________________________________
215 : fParticleFileMap(fHeader.GetParticleFileMap())
218 // Default constructor for AliRun
243 fPDGDB = 0; //Particle factory object!
245 fConfigFunction = "\0";
248 fTransParName = "\0";
249 fBaseFileName = ".\0";
251 fParticleMap = new TObjArray(10000);
254 //_____________________________________________________________________________
255 AliRun::AliRun(const char *name, const char *title)
256 : TNamed(name,title),
257 fParticleFileMap(fHeader.GetParticleFileMap())
261 // Constructor for the main processor.
262 // Creates the geometry
263 // Creates the list of Detectors.
264 // Creates the list of particles.
281 fConfigFunction = "Config();";
283 // Set random number generator
284 gRandom = fRandom = new TRandom3();
286 if (gSystem->Getenv("CONFIG_SEED")) {
287 gRandom->SetSeed((UInt_t)atoi(gSystem->Getenv("CONFIG_SEED")));
290 gROOT->GetListOfBrowsables()->Add(this,name);
292 // create the support list for the various Detectors
293 fModules = new TObjArray(77);
295 // Create the TNode geometry for the event display
297 BuildSimpleGeometry();
307 // Create the particle stack
308 fParticles = new TClonesArray("TParticle",1000);
312 // Create default mag field
317 // Prepare the tracking medium lists
318 fImedia = new TArrayI(1000);
319 for(i=0;i<1000;i++) (*fImedia)[i]=-99;
322 fPDGDB = TDatabasePDG::Instance(); //Particle factory object!
324 // Create HitLists list
325 fHitLists = new TList();
328 fBaseFileName = ".\0";
330 fParticleMap = new TObjArray(10000);
334 //_____________________________________________________________________________
338 // Default AliRun destructor
358 fParticles->Delete();
366 //_____________________________________________________________________________
367 void AliRun::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
370 // Add a hit to detector id
372 TObjArray &dets = *fModules;
373 if(dets[id]) ((AliModule*) dets[id])->AddHit(track,vol,hits);
376 //_____________________________________________________________________________
377 void AliRun::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
380 // Add digit to detector id
382 TObjArray &dets = *fModules;
383 if(dets[id]) ((AliModule*) dets[id])->AddDigit(tracks,digits);
386 //_____________________________________________________________________________
387 void AliRun::Browse(TBrowser *b)
390 // Called when the item "Run" is clicked on the left pane
391 // of the Root browser.
392 // It displays the Root Trees and all detectors.
394 if (fTreeK) b->Add(fTreeK,fTreeK->GetName());
395 if (fTreeH) b->Add(fTreeH,fTreeH->GetName());
396 if (fTreeD) b->Add(fTreeD,fTreeD->GetName());
397 if (fTreeE) b->Add(fTreeE,fTreeE->GetName());
398 if (fTreeR) b->Add(fTreeR,fTreeR->GetName());
399 if (fTreeS) b->Add(fTreeS,fTreeS->GetName());
401 TIter next(fModules);
403 while((detector = (AliModule*)next())) {
404 b->Add(detector,detector->GetName());
406 b->Add(fMCQA,"AliMCQA");
409 //_____________________________________________________________________________
413 // Initialize Alice geometry
418 //_____________________________________________________________________________
419 void AliRun::BuildSimpleGeometry()
422 // Create a simple TNode geometry used by Root display engine
424 // Initialise geometry
426 fGeometry = new TGeometry("AliceGeom","Galice Geometry for Hits");
427 new TMaterial("void","Vacuum",0,0,0); //Everything is void
428 TBRIK *brik = new TBRIK("S_alice","alice volume","void",2000,2000,3000);
429 brik->SetVisibility(0);
430 new TNode("alice","alice","S_alice");
433 //_____________________________________________________________________________
434 void AliRun::CleanDetectors()
437 // Clean Detectors at the end of event
439 TIter next(fModules);
441 while((detector = (AliModule*)next())) {
442 detector->FinishEvent();
446 //_____________________________________________________________________________
447 void AliRun::CleanParents()
450 // Clean Particles stack.
451 // Set parent/daughter relations
453 TObjArray &particles = *fParticleMap;
456 for(i=0; i<fHgwmk+1; i++) {
457 part = (TParticle *)particles.At(i);
458 if(part) if(!part->TestBit(kDaughtersBit)) {
459 part->SetFirstDaughter(-1);
460 part->SetLastDaughter(-1);
465 //_____________________________________________________________________________
466 Int_t AliRun::DistancetoPrimitive(Int_t, Int_t)
469 // Return the distance from the mouse to the AliRun object
475 //_____________________________________________________________________________
476 void AliRun::DumpPart (Int_t i) const
479 // Dumps particle i in the stack
481 ((TParticle*) (*fParticleMap)[i])->Print();
484 //_____________________________________________________________________________
485 void AliRun::DumpPStack () const
488 // Dumps the particle stack
490 TObjArray &particles = *fParticleMap;
492 "\n\n=======================================================================\n");
493 for (Int_t i=0;i<fNtrack;i++)
495 printf("-> %d ",i); ((TParticle*) particles[i])->Print();
496 printf("--------------------------------------------------------------\n");
499 "\n=======================================================================\n\n");
502 void AliRun::SetField(AliMagF* magField)
504 // Set Magnetic Field Map
509 //_____________________________________________________________________________
510 void AliRun::SetField(Int_t type, Int_t version, Float_t scale,
511 Float_t maxField, char* filename)
514 // Set magnetic field parameters
515 // type Magnetic field transport flag 0=no field, 2=helix, 3=Runge Kutta
516 // version Magnetic field map version (only 1 active now)
517 // scale Scale factor for the magnetic field
518 // maxField Maximum value for the magnetic field
521 // --- Sanity check on mag field flags
522 if(fField) delete fField;
524 fField = new AliMagFC("Map1"," ",type,scale,maxField);
525 } else if(version<=2) {
526 fField = new AliMagFCM("Map2-3",filename,type,scale,maxField);
528 } else if(version==3) {
529 fField = new AliMagFDM("Map4",filename,type,scale,maxField);
532 Warning("SetField","Invalid map %d\n",version);
536 //_____________________________________________________________________________
537 void AliRun::PreTrack()
539 TObjArray &dets = *fModules;
542 for(Int_t i=0; i<=fNdets; i++)
543 if((module = (AliModule*)dets[i]))
549 //_____________________________________________________________________________
550 void AliRun::PostTrack()
552 TObjArray &dets = *fModules;
555 for(Int_t i=0; i<=fNdets; i++)
556 if((module = (AliModule*)dets[i]))
560 //_____________________________________________________________________________
561 void AliRun::FinishPrimary()
564 // Called at the end of each primary track
567 // static Int_t count=0;
568 // const Int_t times=10;
569 // This primary is finished, purify stack
572 TIter next(fModules);
574 while((detector = (AliModule*)next())) {
575 detector->FinishPrimary();
578 // Write out hits if any
579 if (gAlice->TreeH()) {
580 gAlice->TreeH()->Fill();
587 // if(++count%times==1) gObjectTable->Print();
590 //_____________________________________________________________________________
591 void AliRun::FinishEvent()
594 // Called at the end of the event.
598 if(fLego) fLego->FinishEvent();
600 //Update the energy deposit tables
602 for(i=0;i<fEventEnergy.GetSize();i++) {
603 fSummEnergy[i]+=fEventEnergy[i];
604 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
606 fEventEnergy.Reset();
608 // Clean detector information
611 // Write out the kinematics
614 if(fTreeK->GetEntries() ==0) {
615 // set the fParticleFileMap size for the first time
616 if (fHgwmk+1 > fParticleFileMap.GetSize())
617 fParticleFileMap.Set(fHgwmk+1);
620 Bool_t allFilled = kFALSE;
622 for(i=0; i<fHgwmk+1; ++i) if((part=fParticleMap->At(i))) {
623 fParticleBuffer = (TParticle*) part;
624 fParticleFileMap[i]= (Int_t) fTreeK->GetEntries();
626 (*fParticleMap)[i]=0;
628 // When all primaries were filled no particle!=0
629 // should be left => to be removed later.
630 if (allFilled) printf("Why != 0 part # %d?\n",i);
633 // // printf("Why = 0 part # %d?\n",i); => We know.
635 // we don't break now in order to be sure there is no
636 // particle !=0 left.
637 // To be removed later and replaced with break.
638 if(!allFilled) allFilled = kTRUE;
642 // Set number of tracks to event header
643 fHeader.SetNtrack(fNtrack);
645 // Write out the digits
656 // Write out reconstructed clusters
661 // Write out the event Header information
662 if (fTreeE) fTreeE->Fill();
667 // Write Tree headers
668 if (fTreeK) fTreeK->Write(0,TObject::kOverwrite);
669 if (fTreeH) fTreeH->Write(0,TObject::kOverwrite);
670 if (fTreeD) fTreeD->Write(0,TObject::kOverwrite);
671 if (fTreeR) fTreeR->Write(0,TObject::kOverwrite);
672 if (fTreeS) fTreeS->Write(0,TObject::kOverwrite);
677 //_____________________________________________________________________________
678 void AliRun::FinishRun()
681 // Called at the end of the run.
685 if(fLego) fLego->FinishRun();
687 if(fGenerator) fGenerator->FinishRun();
689 // Clean detector information
690 TIter next(fModules);
692 while((detector = (AliModule*)next())) {
693 detector->FinishRun();
696 //Output energy summary tables
699 TFile *file = fTreeE->GetCurrentFile();
703 fTreeE->Write(0,TObject::kOverwrite);
705 // Write AliRun info and all detectors parameters
706 Write(0,TObject::kOverwrite);
708 // Clean tree information
710 delete fTreeK; fTreeK = 0;
713 delete fTreeH; fTreeH = 0;
716 delete fTreeD; fTreeD = 0;
719 delete fTreeR; fTreeR = 0;
722 delete fTreeE; fTreeE = 0;
729 //_____________________________________________________________________________
730 void AliRun::FlagTrack(Int_t track)
733 // Flags a track and all its family tree to be kept
740 particle=(TParticle*)fParticleMap->At(curr);
742 // If the particle is flagged the three from here upward is saved already
743 if(particle->TestBit(kKeepBit)) return;
745 // Save this particle
746 particle->SetBit(kKeepBit);
748 // Move to father if any
749 if((curr=particle->GetFirstMother())==-1) return;
753 //_____________________________________________________________________________
754 void AliRun::EnergySummary()
757 // Print summary of deposited energy
763 Int_t kn, i, left, j, id;
764 const Float_t kzero=0;
765 Int_t ievent=fHeader.GetEvent()+1;
767 // Energy loss information
769 printf("***************** Energy Loss Information per event (GEV) *****************\n");
770 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
773 fEventEnergy[ndep]=kn;
778 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
781 fSummEnergy[ndep]=ed;
782 fSum2Energy[ndep]=TMath::Min((Float_t) 99.,TMath::Max(ed2,kzero));
787 for(kn=0;kn<(ndep-1)/3+1;kn++) {
789 for(i=0;i<(3<left?3:left);i++) {
791 id=Int_t (fEventEnergy[j]+0.1);
792 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
797 // Relative energy loss in different detectors
798 printf("******************** Relative Energy Loss per event ********************\n");
799 printf("Total energy loss per event %10.3f GeV\n",edtot);
800 for(kn=0;kn<(ndep-1)/5+1;kn++) {
802 for(i=0;i<(5<left?5:left);i++) {
804 id=Int_t (fEventEnergy[j]+0.1);
805 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
809 for(kn=0;kn<75;kn++) printf("*");
813 // Reset the TArray's
814 // fEventEnergy.Set(0);
815 // fSummEnergy.Set(0);
816 // fSum2Energy.Set(0);
819 //_____________________________________________________________________________
820 AliModule *AliRun::GetModule(const char *name) const
823 // Return pointer to detector from name
825 return (AliModule*)fModules->FindObject(name);
828 //_____________________________________________________________________________
829 AliDetector *AliRun::GetDetector(const char *name) const
832 // Return pointer to detector from name
834 return (AliDetector*)fModules->FindObject(name);
837 //_____________________________________________________________________________
838 Int_t AliRun::GetModuleID(const char *name) const
841 // Return galice internal detector identifier from name
844 TObject *mod=fModules->FindObject(name);
845 if(mod) i=fModules->IndexOf(mod);
849 //_____________________________________________________________________________
850 Int_t AliRun::GetEvent(Int_t event)
853 // Connect the Trees Kinematics and Hits for event # event
854 // Set branch addresses
857 // Reset existing structures
863 // Delete Trees already connected
864 if (fTreeK) delete fTreeK;
865 if (fTreeH) delete fTreeH;
866 if (fTreeD) delete fTreeD;
867 if (fTreeR) delete fTreeR;
868 if (fTreeS) delete fTreeS;
870 // Get header from file
871 if(fTreeE) fTreeE->GetEntry(event);
872 else Error("GetEvent","Cannot file Header Tree\n");
873 TFile *file = fTreeE->GetCurrentFile();
877 // Get Kine Tree from file
879 sprintf(treeName,"TreeK%d",event);
880 fTreeK = (TTree*)gDirectory->Get(treeName);
881 if (fTreeK) fTreeK->SetBranchAddress("Particles", &fParticleBuffer);
882 else Error("GetEvent","cannot find Kine Tree for event:%d\n",event);
883 // Create the particle stack
884 if(!fParticles) fParticles = new TClonesArray("TParticle",1000);
885 // Build the pointer list
887 fParticleMap->Clear();
888 fParticleMap->Expand(fTreeK->GetEntries());
890 fParticleMap = new TObjArray(fTreeK->GetEntries());
894 // Get Hits Tree header from file
895 sprintf(treeName,"TreeH%d",event);
896 fTreeH = (TTree*)gDirectory->Get(treeName);
898 Error("GetEvent","cannot find Hits Tree for event:%d\n",event);
903 // Get Digits Tree header from file
904 sprintf(treeName,"TreeD%d",event);
905 fTreeD = (TTree*)gDirectory->Get(treeName);
907 // Warning("GetEvent","cannot find Digits Tree for event:%d\n",event);
912 // Get SDigits Tree header from file
913 sprintf(treeName,"TreeS%d",event);
914 fTreeS = (TTree*)gDirectory->Get(treeName);
916 // Warning("GetEvent","cannot find SDigits Tree for event:%d\n",event);
921 // Get Reconstruct Tree header from file
922 sprintf(treeName,"TreeR%d",event);
923 fTreeR = (TTree*)gDirectory->Get(treeName);
925 // printf("WARNING: cannot find Reconstructed Tree for event:%d\n",event);
930 // Set Trees branch addresses
931 TIter next(fModules);
933 while((detector = (AliModule*)next())) {
934 detector->SetTreeAddress();
937 fNtrack = Int_t (fTreeK->GetEntries());
941 //_____________________________________________________________________________
942 TGeometry *AliRun::GetGeometry()
945 // Import Alice geometry from current file
946 // Return pointer to geometry object
948 if (!fGeometry) fGeometry = (TGeometry*)gDirectory->Get("AliceGeom");
950 // Unlink and relink nodes in detectors
951 // This is bad and there must be a better way...
954 TIter next(fModules);
956 while((detector = (AliModule*)next())) {
957 TList *dnodes=detector->Nodes();
960 for ( j=0; j<dnodes->GetSize(); j++) {
961 node = (TNode*) dnodes->At(j);
962 node1 = fGeometry->GetNode(node->GetName());
963 dnodes->Remove(node);
964 dnodes->AddAt(node1,j);
970 //_____________________________________________________________________________
971 void AliRun::GetNextTrack(Int_t &mtrack, Int_t &ipart, Float_t *pmom,
972 Float_t &e, Float_t *vpos, Float_t *polar,
976 // Return next track from stack of particles
981 for(Int_t i=fNtrack-1; i>=0; i--) {
982 track=(TParticle*) fParticleMap->At(i);
983 if(track) if(!track->TestBit(kDoneBit)) {
985 // The track exists and has not yet been processed
987 ipart=track->GetPdgCode();
995 track->GetPolarisation(pol);
1000 track->SetBit(kDoneBit);
1006 // stop and start timer when we start a primary track
1007 Int_t nprimaries = fHeader.GetNprimary();
1008 if (fCurrent >= nprimaries) return;
1009 if (fCurrent < nprimaries-1) {
1011 track=(TParticle*) fParticleMap->At(fCurrent+1);
1012 // track->SetProcessTime(fTimer.CpuTime());
1017 //_____________________________________________________________________________
1018 Int_t AliRun::GetPrimary(Int_t track) const
1021 // return number of primary that has generated track
1023 int current, parent;
1029 part = (TParticle *)fParticleMap->At(current);
1030 parent=part->GetFirstMother();
1031 if(parent<0) return current;
1035 //_____________________________________________________________________________
1036 void AliRun::InitMC(const char *setup)
1039 // Initialize the Alice setup
1043 Warning("Init","Cannot initialise AliRun twice!\n");
1047 gROOT->LoadMacro(setup);
1048 gInterpreter->ProcessLine(fConfigFunction.Data());
1051 gMC->DefineParticles(); //Create standard MC particles
1053 TObject *objfirst, *objlast;
1055 fNdets = fModules->GetLast()+1;
1058 //=================Create Materials and geometry
1061 // Added also after in case of interactive initialisation of modules
1062 fNdets = fModules->GetLast()+1;
1064 TIter next(fModules);
1065 AliModule *detector;
1066 while((detector = (AliModule*)next())) {
1067 detector->SetTreeAddress();
1068 objlast = gDirectory->GetList()->Last();
1070 // Add Detector histograms in Detector list of histograms
1071 if (objlast) objfirst = gDirectory->GetList()->After(objlast);
1072 else objfirst = gDirectory->GetList()->First();
1074 detector->Histograms()->Add(objfirst);
1075 objfirst = gDirectory->GetList()->After(objfirst);
1078 ReadTransPar(); //Read the cuts for all materials
1080 MediaTable(); //Build the special IMEDIA table
1082 //Initialise geometry deposition table
1083 fEventEnergy.Set(gMC->NofVolumes()+1);
1084 fSummEnergy.Set(gMC->NofVolumes()+1);
1085 fSum2Energy.Set(gMC->NofVolumes()+1);
1087 //Compute cross-sections
1088 gMC->BuildPhysics();
1090 //Write Geometry object to current file.
1095 fMCQA = new AliMCQA(fNdets);
1098 // Save stuff at the beginning of the file to avoid file corruption
1102 //_____________________________________________________________________________
1103 void AliRun::MediaTable()
1106 // Built media table to get from the media number to
1109 Int_t kz, nz, idt, lz, i, k, ind;
1111 TObjArray &dets = *gAlice->Detectors();
1114 // For all detectors
1115 for (kz=0;kz<fNdets;kz++) {
1116 // If detector is defined
1117 if((det=(AliModule*) dets[kz])) {
1118 TArrayI &idtmed = *(det->GetIdtmed());
1119 for(nz=0;nz<100;nz++) {
1120 // Find max and min material number
1121 if((idt=idtmed[nz])) {
1122 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
1123 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
1126 if(det->LoMedium() > det->HiMedium()) {
1127 det->LoMedium() = 0;
1128 det->HiMedium() = 0;
1130 if(det->HiMedium() > fImedia->GetSize()) {
1131 Error("MediaTable","Increase fImedia from %d to %d",
1132 fImedia->GetSize(),det->HiMedium());
1135 // Tag all materials in rage as belonging to detector kz
1136 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
1143 // Print summary table
1144 printf(" Traking media ranges:\n");
1145 for(i=0;i<(fNdets-1)/6+1;i++) {
1146 for(k=0;k< (6<fNdets-i*6?6:fNdets-i*6);k++) {
1148 det=(AliModule*)dets[ind];
1150 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
1153 printf(" %6s: %3d -> %3d;","NULL",0,0);
1159 //____________________________________________________________________________
1160 void AliRun::SetGenerator(AliGenerator *generator)
1163 // Load the event generator
1165 if(!fGenerator) fGenerator = generator;
1168 //____________________________________________________________________________
1169 void AliRun::ResetGenerator(AliGenerator *generator)
1172 // Load the event generator
1176 Warning("ResetGenerator","Replacing generator %s with %s\n",
1177 fGenerator->GetName(),generator->GetName());
1179 Warning("ResetGenerator","Replacing generator %s with NULL\n",
1180 fGenerator->GetName());
1181 fGenerator = generator;
1184 //____________________________________________________________________________
1185 void AliRun::SetTransPar(char *filename)
1187 fTransParName = filename;
1190 //____________________________________________________________________________
1191 void AliRun::SetBaseFile(char *filename)
1193 fBaseFileName = filename;
1196 //____________________________________________________________________________
1197 void AliRun::ReadTransPar()
1200 // Read filename to set the transport parameters
1204 const Int_t kncuts=10;
1205 const Int_t knflags=11;
1206 const Int_t knpars=kncuts+knflags;
1207 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
1208 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
1209 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
1210 "MULS","PAIR","PHOT","RAYL"};
1214 Float_t cut[kncuts];
1215 Int_t flag[knflags];
1216 Int_t i, itmed, iret, ktmed, kz;
1219 // See whether the file is there
1220 filtmp=gSystem->ExpandPathName(fTransParName.Data());
1221 lun=fopen(filtmp,"r");
1224 Warning("ReadTransPar","File %s does not exist!\n",fTransParName.Data());
1228 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
1229 printf(" *%59s\n","*");
1230 printf(" * Please check carefully what you are doing!%10s\n","*");
1231 printf(" *%59s\n","*");
1234 // Initialise cuts and flags
1235 for(i=0;i<kncuts;i++) cut[i]=-99;
1236 for(i=0;i<knflags;i++) flag[i]=-99;
1238 for(i=0;i<256;i++) line[i]='\0';
1239 // Read up to the end of line excluded
1240 iret=fscanf(lun,"%[^\n]",line);
1244 printf(" *%59s\n","*");
1245 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
1248 // Read the end of line
1251 if(line[0]=='*') continue;
1253 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",
1254 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
1255 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
1256 &flag[8],&flag[9],&flag[10]);
1260 Warning("ReadTransPar","Error reading file %s\n",fTransParName.Data());
1263 // Check that the module exist
1264 AliModule *mod = GetModule(detName);
1266 // Get the array of media numbers
1267 TArrayI &idtmed = *mod->GetIdtmed();
1268 // Check that the tracking medium code is valid
1269 if(0<=itmed && itmed < 100) {
1270 ktmed=idtmed[itmed];
1272 Warning("ReadTransPar","Invalid tracking medium code %d for %s\n",itmed,mod->GetName());
1275 // Set energy thresholds
1276 for(kz=0;kz<kncuts;kz++) {
1278 printf(" * %-6s set to %10.3E for tracking medium code %4d for %s\n",
1279 kpars[kz],cut[kz],itmed,mod->GetName());
1280 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
1283 // Set transport mechanisms
1284 for(kz=0;kz<knflags;kz++) {
1286 printf(" * %-6s set to %10d for tracking medium code %4d for %s\n",
1287 kpars[kncuts+kz],flag[kz],itmed,mod->GetName());
1288 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
1292 Warning("ReadTransPar","Invalid medium code %d *\n",itmed);
1296 Warning("ReadTransPar","Module %s not present\n",detName);
1302 //_____________________________________________________________________________
1303 TBranch* AliRun::MakeBranchInTree(TTree *tree, const char* name, void* address, Int_t size, char *file)
1306 // Makes branch in given tree and diverts them to a separate file
1309 printf("* MakeBranch * Making Branch %s \n",name);
1311 TBranch *branch = tree->Branch(name,address,size);
1314 char * outFile = new char[strlen(gAlice->GetBaseFile())+strlen(file)+2];
1315 sprintf(outFile,"%s/%s",gAlice->GetBaseFile(),file);
1316 TDirectory *cwd = gDirectory;
1317 branch->SetFile(outFile);
1318 TIter next( branch->GetListOfBranches());
1319 while ((branch=(TBranch*)next())) {
1320 branch->SetFile(outFile);
1323 printf("* MakeBranch * Diverting Branch %s to file %s\n",name,file);
1331 //_____________________________________________________________________________
1332 TBranch* AliRun::MakeBranchInTree(TTree *tree, const char* name, const char *classname, void* address, Int_t size, Int_t splitlevel, char *file)
1335 // Makes branch in given tree and diverts them to a separate file
1337 TDirectory *cwd = gDirectory;
1338 TBranch *branch = tree->Branch(name,classname,address,size,splitlevel);
1341 printf("* MakeBranch * Making Branch %s \n",name);
1343 char * outFile = new char[strlen(gAlice->GetBaseFile())+strlen(file)+2];
1344 sprintf(outFile,"%s/%s",gAlice->GetBaseFile(),file);
1345 branch->SetFile(outFile);
1346 TIter next( branch->GetListOfBranches());
1347 while ((branch=(TBranch*)next())) {
1348 branch->SetFile(outFile);
1351 printf("* MakeBranch * Diverting Branch %s to file %s\n",name,file);
1357 //_____________________________________________________________________________
1358 void AliRun::MakeTree(Option_t *option, char *file)
1361 // Create the ROOT trees
1362 // Loop on all detectors to create the Root branch (if any)
1368 const char *oK = strstr(option,"K");
1369 const char *oH = strstr(option,"H");
1370 const char *oE = strstr(option,"E");
1371 const char *oD = strstr(option,"D");
1372 const char *oR = strstr(option,"R");
1373 const char *oS = strstr(option,"S");
1376 if (oK && !fTreeK) {
1377 sprintf(hname,"TreeK%d",fEvent);
1378 fTreeK = new TTree(hname,"Kinematics");
1379 // Create a branch for particles
1380 MakeBranchInTree(fTreeK,
1381 "Particles", "TParticle", &fParticleBuffer, 4000, 1, file) ;
1384 if (oH && !fTreeH) {
1385 sprintf(hname,"TreeH%d",fEvent);
1386 fTreeH = new TTree(hname,"Hits");
1387 fTreeH->SetAutoSave(1000000000); //no autosave
1390 if (oD && !fTreeD) {
1391 sprintf(hname,"TreeD%d",fEvent);
1392 fTreeD = new TTree(hname,"Digits");
1395 if (oS && !fTreeS) {
1396 sprintf(hname,"TreeS%d",fEvent);
1397 fTreeS = new TTree(hname,"SDigits");
1400 if (oR && !fTreeR) {
1401 sprintf(hname,"TreeR%d",fEvent);
1402 fTreeR = new TTree(hname,"Reconstruction");
1405 if (oE && !fTreeE) {
1406 fTreeE = new TTree("TE","Header");
1408 = MakeBranchInTree(fTreeE,
1409 "Header", "AliHeader", &gAliHeader, 4000, 0, file) ;
1410 branch->SetAutoDelete(kFALSE);
1415 // Create a branch for hits/digits for each detector
1416 // Each branch is a TClonesArray. Each data member of the Hits classes
1417 // will be in turn a subbranch of the detector master branch
1418 TIter next(fModules);
1419 AliModule *detector;
1420 while((detector = (AliModule*)next())) {
1421 if (oH) detector->MakeBranch(option,file);
1425 //_____________________________________________________________________________
1426 Int_t AliRun::PurifyKine(Int_t lastSavedTrack, Int_t nofTracks)
1429 // PurifyKine with external parameters
1431 fHgwmk = lastSavedTrack;
1432 fNtrack = nofTracks;
1437 //_____________________________________________________________________________
1438 TParticle* AliRun::Particle(Int_t i)
1440 if(!(*fParticleMap)[i]) {
1441 Int_t nentries = fParticles->GetEntries();
1443 // algorithmic way of getting entry index
1444 // (primary particles are filled after secondaries)
1446 if (i<fHeader.GetNprimary())
1447 entry = i+fHeader.GetNsecondary();
1449 entry = i-fHeader.GetNprimary();
1451 // only check the algorithmic way and give
1452 // the fatal error if it is wrong
1453 if (entry != fParticleFileMap[i]) {
1455 "!!!! The algorithmic way is WRONG: !!!\n entry: %d map: %d",
1456 entry, fParticleFileMap[i]);
1459 fTreeK->GetEntry(fParticleFileMap[i]);
1460 //fTreeK->GetEntry(entry);
1461 new ((*fParticles)[nentries]) TParticle(*fParticleBuffer);
1462 fParticleMap->AddAt((*fParticles)[nentries],i);
1464 return (TParticle *) (*fParticleMap)[i];
1467 //_____________________________________________________________________________
1468 void AliRun::PurifyKine()
1471 // Compress kinematic tree keeping only flagged particles
1472 // and renaming the particle id's in all the hits
1474 // TClonesArray &particles = *fParticles;
1475 TObjArray &particles = *fParticleMap;
1476 int nkeep=fHgwmk+1, parent, i;
1477 TParticle *part, *father;
1478 TArrayI map(particles.GetLast()+1);
1480 // Save in Header total number of tracks before compression
1481 fHeader.SetNtrack(fHeader.GetNtrack()+fNtrack-fHgwmk);
1483 // If no tracks generated return now
1484 if(fHgwmk+1 == fNtrack) return;
1486 Int_t toshrink = fNtrack-fHgwmk-1;
1488 // First pass, invalid Daughter information
1489 for(i=0; i<fNtrack; i++) {
1490 // Preset map, to be removed later
1491 if(i<=fHgwmk) map[i]=i ;
1494 // particles.UncheckedAt(i)->ResetBit(kDaughtersBit);
1495 if((part=(TParticle*) particles.At(i))) part->ResetBit(kDaughtersBit);
1498 // Invalid daughter information for the parent of the first particle
1499 // generated. This may or may not be the current primary according to
1500 // whether decays have been recorded among the primaries
1501 part = (TParticle *)particles.At(fHgwmk+1);
1502 particles.At(part->GetFirstMother())->ResetBit(kDaughtersBit);
1503 // Second pass, build map between old and new numbering
1504 for(i=fHgwmk+1; i<fNtrack; i++) {
1505 if(particles.At(i)->TestBit(kKeepBit)) {
1507 // This particle has to be kept
1509 // If old and new are different, have to move the pointer
1510 if(i!=nkeep) particles[nkeep]=particles.At(i);
1511 part = (TParticle*) particles.At(nkeep);
1513 // as the parent is always *before*, it must be already
1514 // in place. This is what we are checking anyway!
1515 if((parent=part->GetFirstMother())>fHgwmk)
1516 if(map[parent]==-99) Fatal("PurifyKine","map[%d] = -99!\n",parent);
1517 else part->SetFirstMother(map[parent]);
1523 // Fix daughters information
1524 for (i=fHgwmk+1; i<nkeep; i++) {
1525 part = (TParticle *)particles.At(i);
1526 parent = part->GetFirstMother();
1528 father = (TParticle *)particles.At(parent);
1529 if(father->TestBit(kDaughtersBit)) {
1531 if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
1532 if(i>father->GetLastDaughter()) father->SetLastDaughter(i);
1534 // Initialise daughters info for first pass
1535 father->SetFirstDaughter(i);
1536 father->SetLastDaughter(i);
1537 father->SetBit(kDaughtersBit);
1542 // Now loop on all registered hit lists
1543 TIter next(fHitLists);
1544 TCollection *hitList;
1545 while((hitList = (TCollection*)next())) {
1546 TIter nexthit(hitList);
1548 while((hit = (AliHit*)nexthit())) {
1549 hit->SetTrack(map[hit->GetTrack()]);
1554 // This for detectors which have a special mapping mechanism
1555 // for hits, such as TPC and TRD
1558 TIter nextmod(fModules);
1559 AliModule *detector;
1560 while((detector = (AliModule*)nextmod())) {
1561 detector->RemapTrackHitIDs(map.GetArray());
1564 // Now the output bit, from fHgwmk to nkeep we write everything and we erase
1565 if(nkeep>fParticleFileMap.GetSize()) fParticleFileMap.Set(Int_t (nkeep*1.5));
1568 for (i=fHgwmk+1; i<nkeep; ++i) {
1569 fParticleBuffer = (TParticle*) particles.At(i);
1570 fParticleFileMap[i]=(Int_t) fTreeK->GetEntries();
1575 for (i=nkeep; i<fNtrack; ++i) particles[i]=0;
1577 fLoadPoint-=toshrink;
1578 for(i=fLoadPoint; i<fLoadPoint+toshrink; ++i) fParticles->RemoveAt(i);
1585 //_____________________________________________________________________________
1586 void AliRun::BeginEvent()
1589 // Reset all Detectors & kinematics & trees
1596 fLego->BeginEvent();
1606 // Initialise event header
1607 fHeader.Reset(fRun,fEvent);
1611 sprintf(hname,"TreeK%d",fEvent);
1612 fTreeK->SetName(hname);
1616 sprintf(hname,"TreeH%d",fEvent);
1617 fTreeH->SetName(hname);
1621 sprintf(hname,"TreeD%d",fEvent);
1622 fTreeD->SetName(hname);
1626 sprintf(hname,"TreeS%d",fEvent);
1627 fTreeS->SetName(hname);
1631 sprintf(hname,"TreeR%d",fEvent);
1632 fTreeR->SetName(hname);
1635 //_____________________________________________________________________________
1636 void AliRun::ResetDigits()
1639 // Reset all Detectors digits
1641 TIter next(fModules);
1642 AliModule *detector;
1643 while((detector = (AliModule*)next())) {
1644 detector->ResetDigits();
1648 //_____________________________________________________________________________
1649 void AliRun::ResetSDigits()
1652 // Reset all Detectors digits
1654 TIter next(fModules);
1655 AliModule *detector;
1656 while((detector = (AliModule*)next())) {
1657 detector->ResetSDigits();
1661 //_____________________________________________________________________________
1662 void AliRun::ResetHits()
1665 // Reset all Detectors hits
1667 TIter next(fModules);
1668 AliModule *detector;
1669 while((detector = (AliModule*)next())) {
1670 detector->ResetHits();
1674 //_____________________________________________________________________________
1675 void AliRun::ResetPoints()
1678 // Reset all Detectors points
1680 TIter next(fModules);
1681 AliModule *detector;
1682 while((detector = (AliModule*)next())) {
1683 detector->ResetPoints();
1687 //_____________________________________________________________________________
1688 void AliRun::RunMC(Int_t nevent, const char *setup)
1691 // Main function to be called to process a galice run
1693 // Root > gAlice.Run();
1694 // a positive number of events will cause the finish routine
1698 // check if initialisation has been done
1699 if (!fInitDone) InitMC(setup);
1701 // Create the Root Tree with one branch per detector
1705 if (gSystem->Getenv("CONFIG_SPLIT_FILE")) {
1706 MakeTree("K","Kine.root");
1707 MakeTree("H","Hits.root");
1712 gMC->ProcessRun(nevent);
1714 // End of this run, close files
1715 if(nevent>0) FinishRun();
1718 //_____________________________________________________________________________
1719 void AliRun::RunReco(const char *detector)
1722 // Main function to be called to reconstruct Alice event
1726 Digits2Reco(detector);
1729 //_____________________________________________________________________________
1731 void AliRun::Hits2Digits(const char *selected)
1733 // Convert Hits to sumable digits
1735 Hits2SDigits(selected);
1736 SDigits2Digits(selected);
1740 //_____________________________________________________________________________
1742 void AliRun::Tree2Tree(Option_t *option, const char *selected)
1745 // Function to transform the content of
1747 // - TreeH to TreeS (option "S")
1748 // - TreeS to TreeD (option "D")
1749 // - TreeD to TreeR (option "R")
1751 // If multiple options are specified ("SDR"), transformation will be done in sequence for
1752 // selected detector and for all detectors if none is selected (detector string
1753 // can contain blank separated list of detector names).
1756 const char *oS = strstr(option,"S");
1757 const char *oD = strstr(option,"D");
1758 const char *oR = strstr(option,"R");
1760 gAlice->GetEvent(0);
1762 TObjArray *detectors = gAlice->Detectors();
1764 TIter next(detectors);
1766 AliDetector *detector = 0;
1768 TDirectory *cwd = gDirectory;
1772 while((detector = (AliDetector*)next())) {
1774 if (strcmp(detector->GetName(),selected)) continue;
1775 if (detector->IsActive()){
1777 cout << "Processing " << detector->GetName() << "..." << endl;
1778 if (gSystem->Getenv("CONFIG_SPLIT_FILE")) {
1780 sprintf(outFile,"SDigits.%s.root",detector->GetName());
1781 detector->MakeBranch("S",outFile);
1784 sprintf(outFile,"Digits.%s.root",detector->GetName());
1785 detector->MakeBranch("D",outFile);
1788 sprintf(outFile,"Reco.%s.root",detector->GetName());
1789 detector->MakeBranch("R",outFile);
1792 detector->MakeBranch(option);
1798 detector->Hits2SDigits();
1800 detector->SDigits2Digits();
1802 detector->Digits2Reco();
1811 //_____________________________________________________________________________
1812 void AliRun::RunLego(const char *setup, Int_t nc1, Float_t c1min,
1813 Float_t c1max,Int_t nc2,Float_t c2min,Float_t c2max,
1814 Float_t rmin,Float_t rmax,Float_t zmax, AliLegoGenerator* gener)
1817 // Generates lego plots of:
1818 // - radiation length map phi vs theta
1819 // - radiation length map phi vs eta
1820 // - interaction length map
1821 // - g/cm2 length map
1823 // ntheta bins in theta, eta
1824 // themin minimum angle in theta (degrees)
1825 // themax maximum angle in theta (degrees)
1827 // phimin minimum angle in phi (degrees)
1828 // phimax maximum angle in phi (degrees)
1829 // rmin minimum radius
1830 // rmax maximum radius
1833 // The number of events generated = ntheta*nphi
1834 // run input parameters in macro setup (default="Config.C")
1836 // Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
1839 <img src="picts/AliRunLego1.gif">
1844 <img src="picts/AliRunLego2.gif">
1849 <img src="picts/AliRunLego3.gif">
1854 // check if initialisation has been done
1855 if (!fInitDone) InitMC(setup);
1856 //Save current generator
1857 AliGenerator *gen=Generator();
1859 // Set new generator
1860 if (!gener) gener = new AliLegoGenerator();
1861 ResetGenerator(gener);
1863 // Configure Generator
1864 gener->SetRadiusRange(rmin, rmax);
1865 gener->SetZMax(zmax);
1866 gener->SetCoor1Range(nc1, c1min, c1max);
1867 gener->SetCoor2Range(nc2, c2min, c2max);
1870 //Create Lego object
1871 fLego = new AliLego("lego",gener);
1873 //Prepare MC for Lego Run
1878 gMC->ProcessRun(nc1*nc2+1);
1880 // Create only the Root event Tree
1883 // End of this run, close files
1885 // Restore current generator
1886 ResetGenerator(gen);
1887 // Delete Lego Object
1888 delete fLego; fLego=0;
1891 //_____________________________________________________________________________
1892 void AliRun::SetConfigFunction(const char * config)
1895 // Set the signature of the function contained in Config.C to configure
1898 fConfigFunction=config;
1901 //_____________________________________________________________________________
1902 void AliRun::SetCurrentTrack(Int_t track)
1905 // Set current track number
1910 //_____________________________________________________________________________
1911 void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
1912 Float_t *vpos, Float_t *polar, Float_t tof,
1913 AliMCProcess mech, Int_t &ntr, Float_t weight)
1916 // Load a track on the stack
1918 // done 0 if the track has to be transported
1920 // parent identifier of the parent track. -1 for a primary
1921 // pdg particle code
1922 // pmom momentum GeV/c
1924 // polar polarisation
1925 // tof time of flight in seconds
1926 // mecha production mechanism
1927 // ntr on output the number of the track stored
1929 TClonesArray &particles = *fParticles;
1930 TParticle *particle;
1932 const Int_t kfirstdaughter=-1;
1933 const Int_t klastdaughter=-1;
1935 // const Float_t tlife=0;
1938 // Here we get the static mass
1939 // For MC is ok, but a more sophisticated method could be necessary
1940 // if the calculated mass is required
1941 // also, this method is potentially dangerous if the mass
1942 // used in the MC is not the same of the PDG database
1944 mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
1945 Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
1946 pmom[1]*pmom[1]+pmom[2]*pmom[2]);
1948 //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",
1949 //pname,mass,e,fNtrack,pdg,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],kS,mecha);
1951 particle=new(particles[fLoadPoint++]) TParticle(pdg,kS,parent,-1,kfirstdaughter,
1952 klastdaughter,pmom[0],pmom[1],pmom[2],
1953 e,vpos[0],vpos[1],vpos[2],tof);
1954 particle->SetPolarisation(TVector3(polar[0],polar[1],polar[2]));
1955 particle->SetWeight(weight);
1956 particle->SetUniqueID(mech);
1957 if(!done) particle->SetBit(kDoneBit);
1958 // Declare that the daughter information is valid
1959 particle->SetBit(kDaughtersBit);
1960 // Add the particle to the stack
1961 fParticleMap->AddAtAndExpand(particle,fNtrack);
1964 particle=(TParticle*) fParticleMap->At(parent);
1965 particle->SetLastDaughter(fNtrack);
1966 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
1969 // This is a primary track. Set high water mark for this event
1972 // Set also number if primary tracks
1973 fHeader.SetNprimary(fHgwmk+1);
1974 fHeader.SetNtrack(fHgwmk+1);
1980 // Here we get the static mass
1981 // For MC is ok, but a more sophisticated method could be necessary
1982 // if the calculated mass is required
1983 // also, this method is potentially dangerous if the mass
1984 // used in the MC is not the same of the PDG database
1986 Float_t mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
1987 Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
1988 pmom[1]*pmom[1]+pmom[2]*pmom[2]);
1990 SetTrack(done, parent, pdg, pmom[0], pmom[1], pmom[2], e,
1991 vpos[0], vpos[1], vpos[2], tof, polar[0],polar[1],polar[2],
1996 //_____________________________________________________________________________
1997 void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg,
1998 Double_t px, Double_t py, Double_t pz, Double_t e,
1999 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
2000 Double_t polx, Double_t poly, Double_t polz,
2001 AliMCProcess mech, Int_t &ntr, Float_t weight)
2004 // Load a track on the stack
2006 // done 0 if the track has to be transported
2008 // parent identifier of the parent track. -1 for a primary
2009 // pdg particle code
2010 // kS generation status code
2011 // px, py, pz momentum GeV/c
2012 // vx, vy, vz position
2013 // polar polarisation
2014 // tof time of flight in seconds
2015 // mech production mechanism
2016 // ntr on output the number of the track stored
2018 // New method interface:
2019 // arguments were changed to be in correspondence with TParticle
2021 // Note: the energy is not calculated from the static mass but
2022 // it is passed by argument e.
2024 TClonesArray &particles = *fParticles;
2027 const Int_t kFirstDaughter=-1;
2028 const Int_t kLastDaughter=-1;
2031 = new(particles[fLoadPoint++]) TParticle(pdg, kS, parent, -1,
2032 kFirstDaughter, kLastDaughter,
2033 px, py, pz, e, vx, vy, vz, tof);
2035 particle->SetPolarisation(polx, poly, polz);
2036 particle->SetWeight(weight);
2037 particle->SetUniqueID(mech);
2039 if(!done) particle->SetBit(kDoneBit);
2041 // Declare that the daughter information is valid
2042 particle->SetBit(kDaughtersBit);
2043 // Add the particle to the stack
2044 fParticleMap->AddAtAndExpand(particle,fNtrack);
2047 particle=(TParticle*) fParticleMap->At(parent);
2048 particle->SetLastDaughter(fNtrack);
2049 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
2052 // This is a primary track. Set high water mark for this event
2055 // Set also number if primary tracks
2056 fHeader.SetNprimary(fHgwmk+1);
2057 fHeader.SetNtrack(fHgwmk+1);
2062 //_____________________________________________________________________________
2063 void AliRun::SetHighWaterMark(const Int_t nt)
2066 // Set high water mark for last track in event
2069 // Set also number if primary tracks
2070 fHeader.SetNprimary(fHgwmk+1);
2071 fHeader.SetNtrack(fHgwmk+1);
2074 //_____________________________________________________________________________
2075 void AliRun::KeepTrack(const Int_t track)
2078 // flags a track to be kept
2080 fParticleMap->At(track)->SetBit(kKeepBit);
2083 //_____________________________________________________________________________
2084 void AliRun::StepManager(Int_t id)
2087 // Called at every step during transport
2091 // --- If lego option, do it and leave
2093 fLego->StepManager();
2096 //Update energy deposition tables
2097 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
2099 //Call the appropriate stepping routine;
2100 AliModule *det = (AliModule*)fModules->At(id);
2102 fMCQA->StepManager(id);
2108 //_____________________________________________________________________________
2109 void AliRun::Streamer(TBuffer &R__b)
2111 // Stream an object of class AliRun.
2113 if (R__b.IsReading()) {
2114 if (!gAlice) gAlice = this;
2116 AliRun::Class()->ReadBuffer(R__b, this);
2118 gROOT->GetListOfBrowsables()->Add(this,"Run");
2120 fTreeE = (TTree*)gDirectory->Get("TE");
2121 if (fTreeE) fTreeE->SetBranchAddress("Header", &gAliHeader);
2122 else Error("Streamer","cannot find Header Tree\n");
2123 fTreeE->GetEntry(0);
2127 AliRun::Class()->WriteBuffer(R__b, this);