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.37 2000/06/09 20:05:11 morsch
19 Introduce possibility to chose magnetic field version 3: AliMagFDM + field02.dat
21 Revision 1.36 2000/06/08 14:03:58 hristov
22 Only one initializer for a default argument
24 Revision 1.35 2000/06/07 10:13:14 hristov
25 Delete only existent objects.
27 Revision 1.34 2000/05/18 10:45:38 fca
28 Delete Particle Factory properly
30 Revision 1.33 2000/05/16 13:10:40 fca
31 New method IsNewTrack and fix for a problem in Father-Daughter relations
33 Revision 1.32 2000/04/27 10:38:21 fca
34 Correct termination of Lego Run and introduce Lego getter in AliRun
36 Revision 1.31 2000/04/26 10:17:32 fca
37 Changes in Lego for G4 compatibility
39 Revision 1.30 2000/04/18 19:11:40 fca
40 Introduce variable Config.C function signature
42 Revision 1.29 2000/04/07 11:12:34 fca
43 G4 compatibility changes
45 Revision 1.28 2000/04/05 06:51:06 fca
46 Workaround for an HP compiler problem
48 Revision 1.27 2000/03/22 18:08:07 fca
49 Rationalisation of the virtual MC interfaces
51 Revision 1.26 2000/03/22 13:42:26 fca
52 SetGenerator does not replace an existing generator, ResetGenerator does
54 Revision 1.25 2000/02/23 16:25:22 fca
55 AliVMC and AliGeant3 classes introduced
56 ReadEuclid moved from AliRun to AliModule
58 Revision 1.24 2000/01/19 17:17:20 fca
59 Introducing a list of lists of hits -- more hits allowed for detector now
61 Revision 1.23 1999/12/03 11:14:31 fca
62 Fixing previous wrong checking
64 Revision 1.21 1999/11/25 10:40:08 fca
65 Fixing daughters information also in primary tracks
67 Revision 1.20 1999/10/04 18:08:49 fca
68 Adding protection against inconsistent Euclid files
70 Revision 1.19 1999/09/29 07:50:40 fca
71 Introduction of the Copyright and cvs Log
75 ///////////////////////////////////////////////////////////////////////////////
77 // Control class for Alice C++ //
78 // Only one single instance of this class exists. //
79 // The object is created in main program aliroot //
80 // and is pointed by the global gAlice. //
82 // -Supports the list of all Alice Detectors (fModules). //
83 // -Supports the list of particles (fParticles). //
84 // -Supports the Trees. //
85 // -Supports the geometry. //
86 // -Supports the event display. //
89 <img src="picts/AliRunClass.gif">
94 <img src="picts/alirun.gif">
98 ///////////////////////////////////////////////////////////////////////////////
106 #include <TObjectTable.h>
108 #include "TParticle.h"
110 #include "AliDisplay.h"
120 static AliHeader *header;
124 //_____________________________________________________________________________
128 // Default constructor for AliRun
152 fPDGDB = 0; //Particle factory object!
154 fConfigFunction = "\0";
157 //_____________________________________________________________________________
158 AliRun::AliRun(const char *name, const char *title)
162 // Constructor for the main processor.
163 // Creates the geometry
164 // Creates the list of Detectors.
165 // Creates the list of particles.
181 fConfigFunction = "Config();";
183 gROOT->GetListOfBrowsables()->Add(this,name);
185 // create the support list for the various Detectors
186 fModules = new TObjArray(77);
188 // Create the TNode geometry for the event display
190 BuildSimpleGeometry();
200 // Create the particle stack
201 fParticles = new TClonesArray("TParticle",100);
205 // Create default mag field
210 // Prepare the tracking medium lists
211 fImedia = new TArrayI(1000);
212 for(i=0;i<1000;i++) (*fImedia)[i]=-99;
215 fPDGDB = TDatabasePDG::Instance(); //Particle factory object!
217 // Create HitLists list
218 fHitLists = new TList();
221 //_____________________________________________________________________________
225 // Defaullt AliRun destructor
244 fParticles->Delete();
251 //_____________________________________________________________________________
252 void AliRun::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
255 // Add a hit to detector id
257 TObjArray &dets = *fModules;
258 if(dets[id]) ((AliModule*) dets[id])->AddHit(track,vol,hits);
261 //_____________________________________________________________________________
262 void AliRun::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
265 // Add digit to detector id
267 TObjArray &dets = *fModules;
268 if(dets[id]) ((AliModule*) dets[id])->AddDigit(tracks,digits);
271 //_____________________________________________________________________________
272 void AliRun::Browse(TBrowser *b)
275 // Called when the item "Run" is clicked on the left pane
276 // of the Root browser.
277 // It displays the Root Trees and all detectors.
279 if (fTreeK) b->Add(fTreeK,fTreeK->GetName());
280 if (fTreeH) b->Add(fTreeH,fTreeH->GetName());
281 if (fTreeD) b->Add(fTreeD,fTreeD->GetName());
282 if (fTreeE) b->Add(fTreeE,fTreeE->GetName());
283 if (fTreeR) b->Add(fTreeR,fTreeR->GetName());
285 TIter next(fModules);
287 while((detector = (AliModule*)next())) {
288 b->Add(detector,detector->GetName());
292 //_____________________________________________________________________________
296 // Initialize Alice geometry
301 //_____________________________________________________________________________
302 void AliRun::BuildSimpleGeometry()
305 // Create a simple TNode geometry used by Root display engine
307 // Initialise geometry
309 fGeometry = new TGeometry("AliceGeom","Galice Geometry for Hits");
310 new TMaterial("void","Vacuum",0,0,0); //Everything is void
311 TBRIK *brik = new TBRIK("S_alice","alice volume","void",2000,2000,3000);
312 brik->SetVisibility(0);
313 new TNode("alice","alice","S_alice");
316 //_____________________________________________________________________________
317 void AliRun::CleanDetectors()
320 // Clean Detectors at the end of event
322 TIter next(fModules);
324 while((detector = (AliModule*)next())) {
325 detector->FinishEvent();
329 //_____________________________________________________________________________
330 void AliRun::CleanParents()
333 // Clean Particles stack.
334 // Set parent/daughter relations
336 TClonesArray &particles = *(gAlice->Particles());
339 for(i=0; i<fNtrack; i++) {
340 part = (TParticle *)particles.UncheckedAt(i);
341 if(!part->TestBit(Daughters_Bit)) {
342 part->SetFirstDaughter(-1);
343 part->SetLastDaughter(-1);
348 //_____________________________________________________________________________
349 Int_t AliRun::DistancetoPrimitive(Int_t, Int_t)
352 // Return the distance from the mouse to the AliRun object
358 //_____________________________________________________________________________
359 void AliRun::DumpPart (Int_t i)
362 // Dumps particle i in the stack
364 TClonesArray &particles = *fParticles;
365 ((TParticle*) particles[i])->Print();
368 //_____________________________________________________________________________
369 void AliRun::DumpPStack ()
372 // Dumps the particle stack
374 TClonesArray &particles = *fParticles;
376 "\n\n=======================================================================\n");
377 for (Int_t i=0;i<fNtrack;i++)
379 printf("-> %d ",i); ((TParticle*) particles[i])->Print();
380 printf("--------------------------------------------------------------\n");
383 "\n=======================================================================\n\n");
386 //_____________________________________________________________________________
387 void AliRun::SetField(Int_t type, Int_t version, Float_t scale,
388 Float_t maxField, char* filename)
391 // Set magnetic field parameters
392 // type Magnetic field transport flag 0=no field, 2=helix, 3=Runge Kutta
393 // version Magnetic field map version (only 1 active now)
394 // scale Scale factor for the magnetic field
395 // maxField Maximum value for the magnetic field
398 // --- Sanity check on mag field flags
399 if(type<0 || type > 2) {
401 "Invalid magnetic field flag: %5d; Helix tracking chosen instead\n"
405 if(fField) delete fField;
407 fField = new AliMagFC("Map1"," ",type,version,scale,maxField);
408 } else if(version<=2) {
409 fField = new AliMagFCM("Map2-3",filename,type,version,scale,maxField);
411 } else if(version==3) {
412 fField = new AliMagFDM("Map4",filename,type,version,scale,maxField);
415 Warning("SetField","Invalid map %d\n",version);
419 //_____________________________________________________________________________
420 void AliRun::FillTree()
423 // Fills all AliRun TTrees
425 if (fTreeK) fTreeK->Fill();
426 if (fTreeH) fTreeH->Fill();
427 if (fTreeD) fTreeD->Fill();
428 if (fTreeR) fTreeR->Fill();
431 //_____________________________________________________________________________
432 void AliRun::FinishPrimary()
435 // Called at the end of each primary track
438 // static Int_t count=0;
439 // const Int_t times=10;
440 // This primary is finished, purify stack
443 // Write out hits if any
444 if (gAlice->TreeH()) {
445 gAlice->TreeH()->Fill();
452 // if(++count%times==1) gObjectTable->Print();
455 //_____________________________________________________________________________
456 void AliRun::FinishEvent()
459 // Called at the end of the event.
463 if(fLego) fLego->FinishEvent();
465 //Update the energy deposit tables
467 for(i=0;i<fEventEnergy.GetSize();i++) {
468 fSummEnergy[i]+=fEventEnergy[i];
469 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
471 fEventEnergy.Reset();
473 // Clean detector information
476 // Write out the kinematics
482 // Write out the digits
488 // Write out reconstructed clusters
493 // Write out the event Header information
494 if (fTreeE) fTreeE->Fill();
499 // Write Tree headers
500 // Int_t ievent = fHeader.GetEvent();
502 // sprintf(hname,"TreeK%d",ievent);
503 if (fTreeK) fTreeK->Write(0,TObject::kOverwrite);
504 // sprintf(hname,"TreeH%d",ievent);
505 if (fTreeH) fTreeH->Write(0,TObject::kOverwrite);
506 // sprintf(hname,"TreeD%d",ievent);
507 if (fTreeD) fTreeD->Write(0,TObject::kOverwrite);
508 // sprintf(hname,"TreeR%d",ievent);
509 if (fTreeR) fTreeR->Write(0,TObject::kOverwrite);
514 //_____________________________________________________________________________
515 void AliRun::FinishRun()
518 // Called at the end of the run.
522 if(fLego) fLego->FinishRun();
524 // Clean detector information
525 TIter next(fModules);
527 while((detector = (AliModule*)next())) {
528 detector->FinishRun();
531 //Output energy summary tables
534 // file is retrieved from whatever tree
536 if (fTreeK) File = fTreeK->GetCurrentFile();
537 if ((!File) && (fTreeH)) File = fTreeH->GetCurrentFile();
538 if ((!File) && (fTreeD)) File = fTreeD->GetCurrentFile();
539 if ((!File) && (fTreeE)) File = fTreeE->GetCurrentFile();
541 Error("FinishRun","There isn't root file!");
545 fTreeE->Write(0,TObject::kOverwrite);
547 // Clean tree information
549 delete fTreeK; fTreeK = 0;
552 delete fTreeH; fTreeH = 0;
555 delete fTreeD; fTreeD = 0;
558 delete fTreeR; fTreeR = 0;
561 delete fTreeE; fTreeE = 0;
564 // Write AliRun info and all detectors parameters
571 //_____________________________________________________________________________
572 void AliRun::FlagTrack(Int_t track)
575 // Flags a track and all its family tree to be kept
582 particle=(TParticle*)fParticles->UncheckedAt(curr);
584 // If the particle is flagged the three from here upward is saved already
585 if(particle->TestBit(Keep_Bit)) return;
587 // Save this particle
588 particle->SetBit(Keep_Bit);
590 // Move to father if any
591 if((curr=particle->GetFirstMother())==-1) return;
595 //_____________________________________________________________________________
596 void AliRun::EnergySummary()
599 // Print summary of deposited energy
605 Int_t kn, i, left, j, id;
606 const Float_t zero=0;
607 Int_t ievent=fHeader.GetEvent()+1;
609 // Energy loss information
611 printf("***************** Energy Loss Information per event (GEV) *****************\n");
612 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
615 fEventEnergy[ndep]=kn;
620 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,zero))/ed;
623 fSummEnergy[ndep]=ed;
624 fSum2Energy[ndep]=TMath::Min((Float_t) 99.,TMath::Max(ed2,zero));
629 for(kn=0;kn<(ndep-1)/3+1;kn++) {
631 for(i=0;i<(3<left?3:left);i++) {
633 id=Int_t (fEventEnergy[j]+0.1);
634 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
639 // Relative energy loss in different detectors
640 printf("******************** Relative Energy Loss per event ********************\n");
641 printf("Total energy loss per event %10.3f GeV\n",edtot);
642 for(kn=0;kn<(ndep-1)/5+1;kn++) {
644 for(i=0;i<(5<left?5:left);i++) {
646 id=Int_t (fEventEnergy[j]+0.1);
647 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
651 for(kn=0;kn<75;kn++) printf("*");
655 // Reset the TArray's
656 // fEventEnergy.Set(0);
657 // fSummEnergy.Set(0);
658 // fSum2Energy.Set(0);
661 //_____________________________________________________________________________
662 AliModule *AliRun::GetModule(const char *name)
665 // Return pointer to detector from name
667 return (AliModule*)fModules->FindObject(name);
670 //_____________________________________________________________________________
671 AliDetector *AliRun::GetDetector(const char *name)
674 // Return pointer to detector from name
676 return (AliDetector*)fModules->FindObject(name);
679 //_____________________________________________________________________________
680 Int_t AliRun::GetModuleID(const char *name)
683 // Return galice internal detector identifier from name
686 TObject *mod=fModules->FindObject(name);
687 if(mod) i=fModules->IndexOf(mod);
691 //_____________________________________________________________________________
692 Int_t AliRun::GetEvent(Int_t event)
695 // Connect the Trees Kinematics and Hits for event # event
696 // Set branch addresses
699 // Reset existing structures
704 // Delete Trees already connected
705 if (fTreeK) delete fTreeK;
706 if (fTreeH) delete fTreeH;
707 if (fTreeD) delete fTreeD;
708 if (fTreeR) delete fTreeR;
710 // Get header from file
711 if(fTreeE) fTreeE->GetEntry(event);
712 else Error("GetEvent","Cannot file Header Tree\n");
714 // Get Kine Tree from file
716 sprintf(treeName,"TreeK%d",event);
717 fTreeK = (TTree*)gDirectory->Get(treeName);
718 if (fTreeK) fTreeK->SetBranchAddress("Particles", &fParticles);
719 else Error("GetEvent","cannot find Kine Tree for event:%d\n",event);
721 // Get Hits Tree header from file
722 sprintf(treeName,"TreeH%d",event);
723 fTreeH = (TTree*)gDirectory->Get(treeName);
725 Error("GetEvent","cannot find Hits Tree for event:%d\n",event);
728 // Get Digits Tree header from file
729 sprintf(treeName,"TreeD%d",event);
730 fTreeD = (TTree*)gDirectory->Get(treeName);
732 Warning("GetEvent","cannot find Digits Tree for event:%d\n",event);
736 // Get Reconstruct Tree header from file
737 sprintf(treeName,"TreeR%d",event);
738 fTreeR = (TTree*)gDirectory->Get(treeName);
740 // printf("WARNING: cannot find Reconstructed Tree for event:%d\n",event);
743 // Set Trees branch addresses
744 TIter next(fModules);
746 while((detector = (AliModule*)next())) {
747 detector->SetTreeAddress();
750 if (fTreeK) fTreeK->GetEvent(0);
751 fNtrack = Int_t (fParticles->GetEntries());
755 //_____________________________________________________________________________
756 TGeometry *AliRun::GetGeometry()
759 // Import Alice geometry from current file
760 // Return pointer to geometry object
762 if (!fGeometry) fGeometry = (TGeometry*)gDirectory->Get("AliceGeom");
764 // Unlink and relink nodes in detectors
765 // This is bad and there must be a better way...
768 TIter next(fModules);
770 while((detector = (AliModule*)next())) {
771 detector->SetTreeAddress();
772 TList *dnodes=detector->Nodes();
775 for ( j=0; j<dnodes->GetSize(); j++) {
776 node = (TNode*) dnodes->At(j);
777 node1 = fGeometry->GetNode(node->GetName());
778 dnodes->Remove(node);
779 dnodes->AddAt(node1,j);
785 //_____________________________________________________________________________
786 void AliRun::GetNextTrack(Int_t &mtrack, Int_t &ipart, Float_t *pmom,
787 Float_t &e, Float_t *vpos, Float_t *polar,
791 // Return next track from stack of particles
796 for(Int_t i=fNtrack-1; i>=0; i--) {
797 track=(TParticle*) fParticles->UncheckedAt(i);
798 if(!track->TestBit(Done_Bit)) {
800 // The track has not yet been processed
802 ipart=track->GetPdgCode();
810 track->GetPolarisation(pol);
815 track->SetBit(Done_Bit);
821 // stop and start timer when we start a primary track
822 Int_t nprimaries = fHeader.GetNprimary();
823 if (fCurrent >= nprimaries) return;
824 if (fCurrent < nprimaries-1) {
826 track=(TParticle*) fParticles->UncheckedAt(fCurrent+1);
827 // track->SetProcessTime(fTimer.CpuTime());
832 //_____________________________________________________________________________
833 Int_t AliRun::GetPrimary(Int_t track)
836 // return number of primary that has generated track
844 part = (TParticle *)fParticles->UncheckedAt(current);
845 parent=part->GetFirstMother();
846 if(parent<0) return current;
850 //_____________________________________________________________________________
851 void AliRun::InitMC(const char *setup)
854 // Initialize the Alice setup
858 Warning("Init","Cannot initialise AliRun twice!\n");
862 gROOT->LoadMacro(setup);
863 gInterpreter->ProcessLine(fConfigFunction.Data());
865 gMC->DefineParticles(); //Create standard MC particles
867 TObject *objfirst, *objlast;
869 fNdets = fModules->GetLast()+1;
872 //=================Create Materials and geometry
875 TIter next(fModules);
877 while((detector = (AliModule*)next())) {
878 detector->SetTreeAddress();
879 objlast = gDirectory->GetList()->Last();
881 // Add Detector histograms in Detector list of histograms
882 if (objlast) objfirst = gDirectory->GetList()->After(objlast);
883 else objfirst = gDirectory->GetList()->First();
885 detector->Histograms()->Add(objfirst);
886 objfirst = gDirectory->GetList()->After(objfirst);
889 SetTransPar(); //Read the cuts for all materials
891 MediaTable(); //Build the special IMEDIA table
893 //Initialise geometry deposition table
894 fEventEnergy.Set(gMC->NofVolumes()+1);
895 fSummEnergy.Set(gMC->NofVolumes()+1);
896 fSum2Energy.Set(gMC->NofVolumes()+1);
898 //Compute cross-sections
901 //Write Geometry object to current file.
907 // Save stuff at the beginning of the file to avoid file corruption
911 //_____________________________________________________________________________
912 void AliRun::MediaTable()
915 // Built media table to get from the media number to
918 Int_t kz, nz, idt, lz, i, k, ind;
920 TObjArray &dets = *gAlice->Detectors();
924 for (kz=0;kz<fNdets;kz++) {
925 // If detector is defined
926 if((det=(AliModule*) dets[kz])) {
927 TArrayI &idtmed = *(det->GetIdtmed());
928 for(nz=0;nz<100;nz++) {
929 // Find max and min material number
930 if((idt=idtmed[nz])) {
931 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
932 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
935 if(det->LoMedium() > det->HiMedium()) {
939 if(det->HiMedium() > fImedia->GetSize()) {
940 Error("MediaTable","Increase fImedia from %d to %d",
941 fImedia->GetSize(),det->HiMedium());
944 // Tag all materials in rage as belonging to detector kz
945 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
952 // Print summary table
953 printf(" Traking media ranges:\n");
954 for(i=0;i<(fNdets-1)/6+1;i++) {
955 for(k=0;k< (6<fNdets-i*6?6:fNdets-i*6);k++) {
957 det=(AliModule*)dets[ind];
959 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
962 printf(" %6s: %3d -> %3d;","NULL",0,0);
968 //____________________________________________________________________________
969 void AliRun::SetGenerator(AliGenerator *generator)
972 // Load the event generator
974 if(!fGenerator) fGenerator = generator;
977 //____________________________________________________________________________
978 void AliRun::ResetGenerator(AliGenerator *generator)
981 // Load the event generator
985 Warning("ResetGenerator","Replacing generator %s with %s\n",
986 fGenerator->GetName(),generator->GetName());
988 Warning("ResetGenerator","Replacing generator %s with NULL\n",
989 fGenerator->GetName());
990 fGenerator = generator;
993 //____________________________________________________________________________
994 void AliRun::SetTransPar(char* filename)
997 // Read filename to set the transport parameters
1001 const Int_t ncuts=10;
1002 const Int_t nflags=11;
1003 const Int_t npars=ncuts+nflags;
1004 const char pars[npars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
1005 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
1006 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
1007 "MULS","PAIR","PHOT","RAYL"};
1013 Int_t i, itmed, iret, ktmed, kz;
1016 // See whether the file is there
1017 filtmp=gSystem->ExpandPathName(filename);
1018 lun=fopen(filtmp,"r");
1021 Warning("SetTransPar","File %s does not exist!\n",filename);
1025 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
1026 printf(" *%59s\n","*");
1027 printf(" * Please check carefully what you are doing!%10s\n","*");
1028 printf(" *%59s\n","*");
1031 // Initialise cuts and flags
1032 for(i=0;i<ncuts;i++) cut[i]=-99;
1033 for(i=0;i<nflags;i++) flag[i]=-99;
1035 for(i=0;i<256;i++) line[i]='\0';
1036 // Read up to the end of line excluded
1037 iret=fscanf(lun,"%[^\n]",line);
1041 printf(" *%59s\n","*");
1042 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
1045 // Read the end of line
1048 if(line[0]=='*') continue;
1050 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",
1051 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
1052 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
1053 &flag[8],&flag[9],&flag[10]);
1057 Warning("SetTransPar","Error reading file %s\n",filename);
1060 // Check that the module exist
1061 AliModule *mod = GetModule(detName);
1063 // Get the array of media numbers
1064 TArrayI &idtmed = *mod->GetIdtmed();
1065 // Check that the tracking medium code is valid
1066 if(0<=itmed && itmed < 100) {
1067 ktmed=idtmed[itmed];
1069 Warning("SetTransPar","Invalid tracking medium code %d for %s\n",itmed,mod->GetName());
1072 // Set energy thresholds
1073 for(kz=0;kz<ncuts;kz++) {
1075 printf(" * %-6s set to %10.3E for tracking medium code %4d for %s\n",
1076 pars[kz],cut[kz],itmed,mod->GetName());
1077 gMC->Gstpar(ktmed,pars[kz],cut[kz]);
1080 // Set transport mechanisms
1081 for(kz=0;kz<nflags;kz++) {
1083 printf(" * %-6s set to %10d for tracking medium code %4d for %s\n",
1084 pars[ncuts+kz],flag[kz],itmed,mod->GetName());
1085 gMC->Gstpar(ktmed,pars[ncuts+kz],Float_t(flag[kz]));
1089 Warning("SetTransPar","Invalid medium code %d *\n",itmed);
1093 Warning("SetTransPar","Module %s not present\n",detName);
1099 //_____________________________________________________________________________
1100 void AliRun::MakeTree(Option_t *option)
1103 // Create the ROOT trees
1104 // Loop on all detectors to create the Root branch (if any)
1110 char *K = strstr(option,"K");
1111 char *H = strstr(option,"H");
1112 char *E = strstr(option,"E");
1113 char *D = strstr(option,"D");
1114 char *R = strstr(option,"R");
1117 sprintf(hname,"TreeK%d",fEvent);
1118 fTreeK = new TTree(hname,"Kinematics");
1119 // Create a branch for particles
1120 fTreeK->Branch("Particles",&fParticles,4000);
1124 sprintf(hname,"TreeH%d",fEvent);
1125 fTreeH = new TTree(hname,"Hits");
1126 fTreeH->SetAutoSave(1000000000); //no autosave
1130 sprintf(hname,"TreeD%d",fEvent);
1131 fTreeD = new TTree(hname,"Digits");
1135 sprintf(hname,"TreeR%d",fEvent);
1136 fTreeR = new TTree(hname,"Reconstruction");
1140 fTreeE = new TTree("TE","Header");
1141 // Create a branch for Header
1142 fTreeE->Branch("Header","AliHeader",&header,4000);
1146 // Create a branch for hits/digits for each detector
1147 // Each branch is a TClonesArray. Each data member of the Hits classes
1148 // will be in turn a subbranch of the detector master branch
1149 TIter next(fModules);
1150 AliModule *detector;
1151 while((detector = (AliModule*)next())) {
1152 if (H || D || R) detector->MakeBranch(option);
1156 //_____________________________________________________________________________
1157 Int_t AliRun::PurifyKine(Int_t lastSavedTrack, Int_t nofTracks)
1160 // PurifyKine with external parameters
1162 fHgwmk = lastSavedTrack;
1163 fNtrack = nofTracks;
1168 //_____________________________________________________________________________
1169 void AliRun::PurifyKine()
1172 // Compress kinematic tree keeping only flagged particles
1173 // and renaming the particle id's in all the hits
1175 TClonesArray &particles = *fParticles;
1176 int nkeep=fHgwmk+1, parent, i;
1177 TParticle *part, *partnew, *father;
1178 int *map = new int[particles.GetEntries()];
1180 // Save in Header total number of tracks before compression
1181 fHeader.SetNtrack(fHeader.GetNtrack()+fNtrack-fHgwmk);
1183 // First pass, invalid Daughter information
1184 for(i=0; i<fNtrack; i++) {
1185 // Preset map, to be removed later
1186 if(i<=fHgwmk) map[i]=i ; else map[i] = -99;
1187 ((TParticle *)particles.UncheckedAt(i))->ResetBit(Daughters_Bit);
1189 // Second pass, build map between old and new numbering
1190 for(i=fHgwmk+1; i<fNtrack; i++) {
1191 part = (TParticle *)particles.UncheckedAt(i);
1192 if(part->TestBit(Keep_Bit)) {
1194 // This particle has to be kept
1198 // Old and new are different, have to copy
1199 partnew = (TParticle *)particles.UncheckedAt(nkeep);
1200 // Change due to a bug in the HP compiler
1201 // *partnew = *part;
1202 memcpy(partnew,part,sizeof(TParticle));
1203 } else partnew = part;
1205 // as the parent is always *before*, it must be already
1206 // in place. This is what we are checking anyway!
1207 if((parent=partnew->GetFirstMother())>fHgwmk) {
1208 if(map[parent]==-99) printf("map[%d] = -99!\n",parent);
1209 partnew->SetFirstMother(map[parent]);
1216 // Fix daughters information
1217 for (i=0; i<fNtrack; i++) {
1218 part = (TParticle *)particles.UncheckedAt(i);
1219 parent = part->GetFirstMother();
1221 father = (TParticle *)particles.UncheckedAt(parent);
1222 if(father->TestBit(Daughters_Bit)) {
1224 if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
1225 if(i>father->GetLastDaughter()) father->SetLastDaughter(i);
1227 // Iitialise daughters info for first pass
1228 father->SetFirstDaughter(i);
1229 father->SetLastDaughter(i);
1230 father->SetBit(Daughters_Bit);
1236 // Now loop on all detectors and reset the hits
1238 TIter next(fModules);
1239 AliModule *detector;
1240 while((detector = (AliModule*)next())) {
1241 if (!detector->Hits()) continue;
1242 TClonesArray &vHits=*(detector->Hits());
1243 if(vHits.GetEntries() != detector->GetNhits())
1244 printf("vHits.GetEntries()!=detector->GetNhits(): %d != %d\n",
1245 vHits.GetEntries(),detector->GetNhits());
1246 for (i=0; i<detector->GetNhits(); i++) {
1247 OneHit = (AliHit *)vHits.UncheckedAt(i);
1248 OneHit->SetTrack(map[OneHit->GetTrack()]);
1253 // Now loop on all registered hit lists
1254 TIter next(fHitLists);
1255 TCollection *hitList;
1256 while((hitList = (TCollection*)next())) {
1257 TIter nexthit(hitList);
1259 while((hit = (AliHit*)nexthit())) {
1260 hit->SetTrack(map[hit->GetTrack()]);
1266 particles.SetLast(fHgwmk);
1270 //_____________________________________________________________________________
1271 void AliRun::BeginEvent()
1274 // Reset all Detectors & kinematics & trees
1281 fLego->BeginEvent();
1290 // Initialise event header
1291 fHeader.Reset(fRun,fEvent);
1295 sprintf(hname,"TreeK%d",fEvent);
1296 fTreeK->SetName(hname);
1300 sprintf(hname,"TreeH%d",fEvent);
1301 fTreeH->SetName(hname);
1305 sprintf(hname,"TreeD%d",fEvent);
1306 fTreeD->SetName(hname);
1310 sprintf(hname,"TreeR%d",fEvent);
1311 fTreeR->SetName(hname);
1315 //_____________________________________________________________________________
1316 void AliRun::ResetDigits()
1319 // Reset all Detectors digits
1321 TIter next(fModules);
1322 AliModule *detector;
1323 while((detector = (AliModule*)next())) {
1324 detector->ResetDigits();
1328 //_____________________________________________________________________________
1329 void AliRun::ResetHits()
1332 // Reset all Detectors hits
1334 TIter next(fModules);
1335 AliModule *detector;
1336 while((detector = (AliModule*)next())) {
1337 detector->ResetHits();
1341 //_____________________________________________________________________________
1342 void AliRun::ResetPoints()
1345 // Reset all Detectors points
1347 TIter next(fModules);
1348 AliModule *detector;
1349 while((detector = (AliModule*)next())) {
1350 detector->ResetPoints();
1354 //_____________________________________________________________________________
1355 void AliRun::RunMC(Int_t nevent, const char *setup)
1358 // Main function to be called to process a galice run
1360 // Root > gAlice.Run();
1361 // a positive number of events will cause the finish routine
1365 // check if initialisation has been done
1366 if (!fInitDone) InitMC(setup);
1368 // Create the Root Tree with one branch per detector
1371 gMC->ProcessRun(nevent);
1373 // End of this run, close files
1374 if(nevent>0) FinishRun();
1377 //_____________________________________________________________________________
1378 void AliRun::RunLego(const char *setup,Int_t ntheta,Float_t themin,
1379 Float_t themax,Int_t nphi,Float_t phimin,Float_t phimax,
1380 Float_t rmin,Float_t rmax,Float_t zmax)
1383 // Generates lego plots of:
1384 // - radiation length map phi vs theta
1385 // - radiation length map phi vs eta
1386 // - interaction length map
1387 // - g/cm2 length map
1389 // ntheta bins in theta, eta
1390 // themin minimum angle in theta (degrees)
1391 // themax maximum angle in theta (degrees)
1393 // phimin minimum angle in phi (degrees)
1394 // phimax maximum angle in phi (degrees)
1395 // rmin minimum radius
1396 // rmax maximum radius
1399 // The number of events generated = ntheta*nphi
1400 // run input parameters in macro setup (default="Config.C")
1402 // Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
1405 <img src="picts/AliRunLego1.gif">
1410 <img src="picts/AliRunLego2.gif">
1415 <img src="picts/AliRunLego3.gif">
1420 // check if initialisation has been done
1421 if (!fInitDone) InitMC(setup);
1423 //Save current generator
1424 AliGenerator *gen=Generator();
1426 //Create Lego object
1427 fLego = new AliLego("lego",ntheta,themin,themax,nphi,phimin,phimax,rmin,rmax,zmax);
1429 //Prepare MC for Lego Run
1433 gMC->ProcessRun(ntheta*nphi+1);
1435 // Create only the Root event Tree
1438 // End of this run, close files
1441 // Delete Lego Object
1442 delete fLego; fLego=0;
1444 // Restore current generator
1448 //_____________________________________________________________________________
1449 void AliRun::SetCurrentTrack(Int_t track)
1452 // Set current track number
1457 //_____________________________________________________________________________
1458 void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
1459 Float_t *vpos, Float_t *polar, Float_t tof,
1460 const char *mecha, Int_t &ntr, Float_t weight)
1463 // Load a track on the stack
1465 // done 0 if the track has to be transported
1467 // parent identifier of the parent track. -1 for a primary
1468 // pdg particle code
1469 // pmom momentum GeV/c
1471 // polar polarisation
1472 // tof time of flight in seconds
1473 // mecha production mechanism
1474 // ntr on output the number of the track stored
1476 TClonesArray &particles = *fParticles;
1477 TParticle *particle;
1479 const Int_t firstdaughter=-1;
1480 const Int_t lastdaughter=-1;
1482 // const Float_t tlife=0;
1485 // Here we get the static mass
1486 // For MC is ok, but a more sophisticated method could be necessary
1487 // if the calculated mass is required
1488 // also, this method is potentially dangerous if the mass
1489 // used in the MC is not the same of the PDG database
1491 mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
1492 Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
1493 pmom[1]*pmom[1]+pmom[2]*pmom[2]);
1495 //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",
1496 //pname,mass,e,fNtrack,pdg,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],KS,mecha);
1498 particle=new(particles[fNtrack]) TParticle(pdg,KS,parent,-1,firstdaughter,
1499 lastdaughter,pmom[0],pmom[1],pmom[2],
1500 e,vpos[0],vpos[1],vpos[2],tof);
1501 // polar[0],polar[1],polar[2],tof,
1503 ((TParticle*)particles[fNtrack])->SetPolarisation(TVector3(polar[0],polar[1],polar[2]));
1504 ((TParticle*)particles[fNtrack])->SetWeight(weight);
1505 if(!done) particle->SetBit(Done_Bit);
1506 //Declare that the daughter information is valid
1507 ((TParticle*)particles[fNtrack])->SetBit(Daughters_Bit);
1510 particle=(TParticle*) fParticles->UncheckedAt(parent);
1511 particle->SetLastDaughter(fNtrack);
1512 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
1515 // This is a primary track. Set high water mark for this event
1518 // Set also number if primary tracks
1519 fHeader.SetNprimary(fHgwmk+1);
1520 fHeader.SetNtrack(fHgwmk+1);
1525 //_____________________________________________________________________________
1526 void AliRun::KeepTrack(const Int_t track)
1529 // flags a track to be kept
1531 TClonesArray &particles = *fParticles;
1532 ((TParticle*)particles[track])->SetBit(Keep_Bit);
1535 //_____________________________________________________________________________
1536 void AliRun::StepManager(Int_t id)
1539 // Called at every step during transport
1543 // --- If lego option, do it and leave
1545 fLego->StepManager();
1548 //Update energy deposition tables
1549 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
1551 //Call the appropriate stepping routine;
1552 AliModule *det = (AliModule*)fModules->At(id);
1553 if(det) det->StepManager();
1557 //_____________________________________________________________________________
1558 void AliRun::Streamer(TBuffer &R__b)
1561 // Stream an object of class AliRun.
1563 if (R__b.IsReading()) {
1564 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
1565 TNamed::Streamer(R__b);
1566 if (!gAlice) gAlice = this;
1567 gROOT->GetListOfBrowsables()->Add(this,"Run");
1568 fTreeE = (TTree*)gDirectory->Get("TE");
1569 if (fTreeE) fTreeE->SetBranchAddress("Header", &header);
1570 else Error("Streamer","cannot find Header Tree\n");
1574 fHeader.Streamer(R__b);
1584 R__b >> fPDGDB; //Particle factory object!
1585 fTreeE->GetEntry(0);
1587 fHeader.SetEvent(0);
1588 fPDGDB = TDatabasePDG::Instance(); //Particle factory object!
1591 fConfigFunction.Streamer(R__b);
1593 fConfigFunction="Config();";
1596 R__b.WriteVersion(AliRun::IsA());
1597 TNamed::Streamer(R__b);
1601 fHeader.Streamer(R__b);
1610 R__b << fPDGDB; //Particle factory object!
1611 fConfigFunction.Streamer(R__b);