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.20 1999/10/04 18:08:49 fca
19 Adding protection against inconsistent Euclid files
21 Revision 1.19 1999/09/29 07:50:40 fca
22 Introduction of the Copyright and cvs Log
26 ///////////////////////////////////////////////////////////////////////////////
28 // Control class for Alice C++ //
29 // Only one single instance of this class exists. //
30 // The object is created in main program aliroot //
31 // and is pointed by the global gAlice. //
33 // -Supports the list of all Alice Detectors (fModules). //
34 // -Supports the list of particles (fParticles). //
35 // -Supports the Trees. //
36 // -Supports the geometry. //
37 // -Supports the event display. //
40 <img src="picts/AliRunClass.gif">
45 <img src="picts/alirun.gif">
49 ///////////////////////////////////////////////////////////////////////////////
57 #include <TObjectTable.h>
59 #include "TParticle.h"
61 #include "AliDisplay.h"
63 #include "AliCallf77.h"
71 static AliHeader *header;
75 # define rxgtrak rxgtrak_
76 # define rxstrak rxstrak_
77 # define rxkeep rxkeep_
78 # define rxouth rxouth_
81 # define rxgtrak RXGTRAK
82 # define rxstrak RXSTRAK
83 # define rxkeep RXKEEP
84 # define rxouth RXOUTH
87 static TArrayF sEventEnergy;
88 static TArrayF sSummEnergy;
89 static TArrayF sSum2Energy;
93 //_____________________________________________________________________________
97 // Default constructor for AliRun
121 fPDGDB = 0; //Particle factory object!
124 //_____________________________________________________________________________
125 AliRun::AliRun(const char *name, const char *title)
129 // Constructor for the main processor.
130 // Creates the geometry
131 // Creates the list of Detectors.
132 // Creates the list of particles.
149 gROOT->GetListOfBrowsables()->Add(this,name);
151 // create the support list for the various Detectors
152 fModules = new TObjArray(77);
154 // Create the TNode geometry for the event display
156 BuildSimpleGeometry();
166 // Create the particle stack
167 fParticles = new TClonesArray("TParticle",100);
171 // Create default mag field
176 // Prepare the tracking medium lists
177 fImedia = new TArrayI(1000);
178 for(i=0;i<1000;i++) (*fImedia)[i]=-99;
181 fPDGDB = TDatabasePDG::Instance(); //Particle factory object!
184 //_____________________________________________________________________________
188 // Defaullt AliRun destructor
207 fParticles->Delete();
212 //_____________________________________________________________________________
213 void AliRun::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
216 // Add a hit to detector id
218 TObjArray &dets = *fModules;
219 if(dets[id]) ((AliModule*) dets[id])->AddHit(track,vol,hits);
222 //_____________________________________________________________________________
223 void AliRun::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
226 // Add digit to detector id
228 TObjArray &dets = *fModules;
229 if(dets[id]) ((AliModule*) dets[id])->AddDigit(tracks,digits);
232 //_____________________________________________________________________________
233 void AliRun::Browse(TBrowser *b)
236 // Called when the item "Run" is clicked on the left pane
237 // of the Root browser.
238 // It displays the Root Trees and all detectors.
240 if (fTreeK) b->Add(fTreeK,fTreeK->GetName());
241 if (fTreeH) b->Add(fTreeH,fTreeH->GetName());
242 if (fTreeD) b->Add(fTreeD,fTreeD->GetName());
243 if (fTreeE) b->Add(fTreeE,fTreeE->GetName());
244 if (fTreeR) b->Add(fTreeR,fTreeR->GetName());
246 TIter next(fModules);
248 while((detector = (AliModule*)next())) {
249 b->Add(detector,detector->GetName());
253 //_____________________________________________________________________________
257 // Initialize Alice geometry
262 //_____________________________________________________________________________
263 void AliRun::BuildSimpleGeometry()
266 // Create a simple TNode geometry used by Root display engine
268 // Initialise geometry
270 fGeometry = new TGeometry("AliceGeom","Galice Geometry for Hits");
271 new TMaterial("void","Vacuum",0,0,0); //Everything is void
272 TBRIK *brik = new TBRIK("S_alice","alice volume","void",2000,2000,3000);
273 brik->SetVisibility(0);
274 new TNode("alice","alice","S_alice");
277 //_____________________________________________________________________________
278 void AliRun::CleanDetectors()
281 // Clean Detectors at the end of event
283 TIter next(fModules);
285 while((detector = (AliModule*)next())) {
286 detector->FinishEvent();
290 //_____________________________________________________________________________
291 void AliRun::CleanParents()
294 // Clean Particles stack.
295 // Set parent/daughter relations
297 TClonesArray &particles = *(gAlice->Particles());
300 for(i=0; i<fNtrack; i++) {
301 part = (TParticle *)particles.UncheckedAt(i);
302 if(!part->TestBit(Daughters_Bit)) {
303 part->SetFirstDaughter(-1);
304 part->SetLastDaughter(-1);
309 //_____________________________________________________________________________
310 Int_t AliRun::DistancetoPrimitive(Int_t, Int_t)
313 // Return the distance from the mouse to the AliRun object
319 //_____________________________________________________________________________
320 void AliRun::DumpPart (Int_t i)
323 // Dumps particle i in the stack
325 TClonesArray &particles = *fParticles;
326 ((TParticle*) particles[i])->Print();
329 //_____________________________________________________________________________
330 void AliRun::DumpPStack ()
333 // Dumps the particle stack
335 TClonesArray &particles = *fParticles;
337 "\n\n=======================================================================\n");
338 for (Int_t i=0;i<fNtrack;i++)
340 printf("-> %d ",i); ((TParticle*) particles[i])->Print();
341 printf("--------------------------------------------------------------\n");
344 "\n=======================================================================\n\n");
347 //_____________________________________________________________________________
348 void AliRun::SetField(Int_t type, Int_t version, Float_t scale,
349 Float_t maxField, char* filename)
352 // Set magnetic field parameters
353 // type Magnetic field transport flag 0=no field, 2=helix, 3=Runge Kutta
354 // version Magnetic field map version (only 1 active now)
355 // scale Scale factor for the magnetic field
356 // maxField Maximum value for the magnetic field
359 // --- Sanity check on mag field flags
360 if(type<0 || type > 2) {
362 "Invalid magnetic field flag: %5d; Helix tracking chosen instead\n"
366 if(fField) delete fField;
368 fField = new AliMagFC("Map1"," ",type,version,scale,maxField);
369 } else if(version<=3) {
370 fField = new AliMagFCM("Map2-3",filename,type,version,scale,maxField);
373 Warning("SetField","Invalid map %d\n",version);
377 //_____________________________________________________________________________
378 void AliRun::FillTree()
381 // Fills all AliRun TTrees
383 if (fTreeK) fTreeK->Fill();
384 if (fTreeH) fTreeH->Fill();
385 if (fTreeD) fTreeD->Fill();
386 if (fTreeR) fTreeR->Fill();
389 //_____________________________________________________________________________
390 void AliRun::FinishPrimary()
393 // Called at the end of each primary track
396 // static Int_t count=0;
397 // const Int_t times=10;
398 // This primary is finished, purify stack
399 gAlice->PurifyKine();
401 // Write out hits if any
402 if (gAlice->TreeH()) {
403 gAlice->TreeH()->Fill();
410 // if(++count%times==1) gObjectTable->Print();
413 //_____________________________________________________________________________
414 void AliRun::FinishEvent()
417 // Called at the end of the event.
420 //Update the energy deposit tables
422 for(i=0;i<sEventEnergy.GetSize();i++) {
423 sSummEnergy[i]+=sEventEnergy[i];
424 sSum2Energy[i]+=sEventEnergy[i]*sEventEnergy[i];
426 sEventEnergy.Reset();
428 // Clean detector information
431 // Write out the kinematics
437 // Write out the digits
443 // Write out reconstructed clusters
448 // Write out the event Header information
449 if (fTreeE) fTreeE->Fill();
454 // Write Tree headers
455 // Int_t ievent = fHeader.GetEvent();
457 // sprintf(hname,"TreeK%d",ievent);
458 if (fTreeK) fTreeK->Write();
459 // sprintf(hname,"TreeH%d",ievent);
460 if (fTreeH) fTreeH->Write();
461 // sprintf(hname,"TreeD%d",ievent);
462 if (fTreeD) fTreeD->Write();
463 // sprintf(hname,"TreeR%d",ievent);
464 if (fTreeR) fTreeR->Write();
467 //_____________________________________________________________________________
468 void AliRun::FinishRun()
471 // Called at the end of the run.
474 // Clean detector information
475 TIter next(fModules);
477 while((detector = (AliModule*)next())) {
478 detector->FinishRun();
481 //Output energy summary tables
484 // file is retrieved from whatever tree
486 if (fTreeK) File = fTreeK->GetCurrentFile();
487 if ((!File) && (fTreeH)) File = fTreeH->GetCurrentFile();
488 if ((!File) && (fTreeD)) File = fTreeD->GetCurrentFile();
489 if ((!File) && (fTreeE)) File = fTreeE->GetCurrentFile();
491 Error("FinishRun","There isn't root file!");
497 // Clean tree information
498 delete fTreeK; fTreeK = 0;
499 delete fTreeH; fTreeH = 0;
500 delete fTreeD; fTreeD = 0;
501 delete fTreeR; fTreeR = 0;
502 delete fTreeE; fTreeE = 0;
504 // Write AliRun info and all detectors parameters
512 //_____________________________________________________________________________
513 void AliRun::FlagTrack(Int_t track)
516 // Flags a track and all its family tree to be kept
523 particle=(TParticle*)fParticles->UncheckedAt(curr);
525 // If the particle is flagged the three from here upward is saved already
526 if(particle->TestBit(Keep_Bit)) return;
528 // Save this particle
529 particle->SetBit(Keep_Bit);
531 // Move to father if any
532 if((curr=particle->GetFirstMother())==-1) return;
536 //_____________________________________________________________________________
537 void AliRun::EnergySummary()
540 // Print summary of deposited energy
546 Int_t kn, i, left, j, id;
547 const Float_t zero=0;
548 Int_t ievent=fHeader.GetEvent()+1;
550 // Energy loss information
552 printf("***************** Energy Loss Information per event (GEV) *****************\n");
553 for(kn=1;kn<sEventEnergy.GetSize();kn++) {
556 sEventEnergy[ndep]=kn;
561 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,zero))/ed;
564 sSummEnergy[ndep]=ed;
565 sSum2Energy[ndep]=TMath::Min((Float_t) 99.,TMath::Max(ed2,zero));
570 for(kn=0;kn<(ndep-1)/3+1;kn++) {
572 for(i=0;i<(3<left?3:left);i++) {
574 id=Int_t (sEventEnergy[j]+0.1);
575 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),sSummEnergy[j],sSum2Energy[j]);
580 // Relative energy loss in different detectors
581 printf("******************** Relative Energy Loss per event ********************\n");
582 printf("Total energy loss per event %10.3f GeV\n",edtot);
583 for(kn=0;kn<(ndep-1)/5+1;kn++) {
585 for(i=0;i<(5<left?5:left);i++) {
587 id=Int_t (sEventEnergy[j]+0.1);
588 printf(" %s %10.3f%%;",gMC->VolName(id),100*sSummEnergy[j]/edtot);
592 for(kn=0;kn<75;kn++) printf("*");
596 // Reset the TArray's
602 //_____________________________________________________________________________
603 AliModule *AliRun::GetModule(const char *name)
606 // Return pointer to detector from name
608 return (AliModule*)fModules->FindObject(name);
611 //_____________________________________________________________________________
612 AliDetector *AliRun::GetDetector(const char *name)
615 // Return pointer to detector from name
617 return (AliDetector*)fModules->FindObject(name);
620 //_____________________________________________________________________________
621 Int_t AliRun::GetModuleID(const char *name)
624 // Return galice internal detector identifier from name
627 TObject *mod=fModules->FindObject(name);
628 if(mod) i=fModules->IndexOf(mod);
632 //_____________________________________________________________________________
633 Int_t AliRun::GetEvent(Int_t event)
636 // Connect the Trees Kinematics and Hits for event # event
637 // Set branch addresses
640 // Reset existing structures
645 // Delete Trees already connected
646 if (fTreeK) delete fTreeK;
647 if (fTreeH) delete fTreeH;
648 if (fTreeD) delete fTreeD;
649 if (fTreeR) delete fTreeR;
651 // Get header from file
652 if(fTreeE) fTreeE->GetEntry(event);
653 else Error("GetEvent","Cannot file Header Tree\n");
655 // Get Kine Tree from file
657 sprintf(treeName,"TreeK%d",event);
658 fTreeK = (TTree*)gDirectory->Get(treeName);
659 if (fTreeK) fTreeK->SetBranchAddress("Particles", &fParticles);
660 else Error("GetEvent","cannot find Kine Tree for event:%d\n",event);
662 // Get Hits Tree header from file
663 sprintf(treeName,"TreeH%d",event);
664 fTreeH = (TTree*)gDirectory->Get(treeName);
666 Error("GetEvent","cannot find Hits Tree for event:%d\n",event);
669 // Get Digits Tree header from file
670 sprintf(treeName,"TreeD%d",event);
671 fTreeD = (TTree*)gDirectory->Get(treeName);
673 Warning("GetEvent","cannot find Digits Tree for event:%d\n",event);
677 // Get Reconstruct Tree header from file
678 sprintf(treeName,"TreeR%d",event);
679 fTreeR = (TTree*)gDirectory->Get(treeName);
681 // printf("WARNING: cannot find Reconstructed Tree for event:%d\n",event);
684 // Set Trees branch addresses
685 TIter next(fModules);
687 while((detector = (AliModule*)next())) {
688 detector->SetTreeAddress();
691 if (fTreeK) fTreeK->GetEvent(0);
692 fNtrack = Int_t (fParticles->GetEntries());
696 //_____________________________________________________________________________
697 TGeometry *AliRun::GetGeometry()
700 // Import Alice geometry from current file
701 // Return pointer to geometry object
703 if (!fGeometry) fGeometry = (TGeometry*)gDirectory->Get("AliceGeom");
705 // Unlink and relink nodes in detectors
706 // This is bad and there must be a better way...
709 TIter next(fModules);
711 while((detector = (AliModule*)next())) {
712 detector->SetTreeAddress();
713 TList *dnodes=detector->Nodes();
716 for ( j=0; j<dnodes->GetSize(); j++) {
717 node = (TNode*) dnodes->At(j);
718 node1 = fGeometry->GetNode(node->GetName());
719 dnodes->Remove(node);
720 dnodes->AddAt(node1,j);
726 //_____________________________________________________________________________
727 void AliRun::GetNextTrack(Int_t &mtrack, Int_t &ipart, Float_t *pmom,
728 Float_t &e, Float_t *vpos, Float_t *polar,
732 // Return next track from stack of particles
737 for(Int_t i=fNtrack-1; i>=0; i--) {
738 track=(TParticle*) fParticles->UncheckedAt(i);
739 if(!track->TestBit(Done_Bit)) {
741 // The track has not yet been processed
743 ipart=track->GetPdgCode();
751 track->GetPolarisation(pol);
756 track->SetBit(Done_Bit);
762 // stop and start timer when we start a primary track
763 Int_t nprimaries = fHeader.GetNprimary();
764 if (fCurrent >= nprimaries) return;
765 if (fCurrent < nprimaries-1) {
767 track=(TParticle*) fParticles->UncheckedAt(fCurrent+1);
768 // track->SetProcessTime(fTimer.CpuTime());
773 //_____________________________________________________________________________
774 Int_t AliRun::GetPrimary(Int_t track)
777 // return number of primary that has generated track
785 part = (TParticle *)fParticles->UncheckedAt(current);
786 parent=part->GetFirstMother();
787 if(parent<0) return current;
791 //_____________________________________________________________________________
792 void AliRun::Init(const char *setup)
795 // Initialize the Alice setup
798 gROOT->LoadMacro(setup);
799 gInterpreter->ProcessLine("Config();");
801 gMC->DefineParticles(); //Create standard MC particles
803 TObject *objfirst, *objlast;
805 fNdets = fModules->GetLast()+1;
808 //=================Create Materials, geometry, histograms, etc
809 TIter next(fModules);
811 while((detector = (AliModule*)next())) {
812 detector->SetTreeAddress();
813 objlast = gDirectory->GetList()->Last();
815 // Initialise detector materials, geometry, histograms,etc
816 detector->CreateMaterials();
817 detector->CreateGeometry();
818 detector->BuildGeometry();
821 // Add Detector histograms in Detector list of histograms
822 if (objlast) objfirst = gDirectory->GetList()->After(objlast);
823 else objfirst = gDirectory->GetList()->First();
825 detector->Histograms()->Add(objfirst);
826 objfirst = gDirectory->GetList()->After(objfirst);
829 SetTransPar(); //Read the cuts for all materials
831 MediaTable(); //Build the special IMEDIA table
833 //Close the geometry structure
836 //Initialise geometry deposition table
837 sEventEnergy.Set(gMC->NofVolumes()+1);
838 sSummEnergy.Set(gMC->NofVolumes()+1);
839 sSum2Energy.Set(gMC->NofVolumes()+1);
841 //Create the color table
844 //Compute cross-sections
847 //Write Geometry object to current file.
853 //_____________________________________________________________________________
854 void AliRun::MediaTable()
857 // Built media table to get from the media number to
860 Int_t kz, nz, idt, lz, i, k, ind;
862 TObjArray &dets = *gAlice->Detectors();
866 for (kz=0;kz<fNdets;kz++) {
867 // If detector is defined
868 if((det=(AliModule*) dets[kz])) {
869 TArrayI &idtmed = *(det->GetIdtmed());
870 for(nz=0;nz<100;nz++) {
871 // Find max and min material number
872 if((idt=idtmed[nz])) {
873 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
874 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
877 if(det->LoMedium() > det->HiMedium()) {
881 if(det->HiMedium() > fImedia->GetSize()) {
882 Error("MediaTable","Increase fImedia from %d to %d",
883 fImedia->GetSize(),det->HiMedium());
886 // Tag all materials in rage as belonging to detector kz
887 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
894 // Print summary table
895 printf(" Traking media ranges:\n");
896 for(i=0;i<(fNdets-1)/6+1;i++) {
897 for(k=0;k< (6<fNdets-i*6?6:fNdets-i*6);k++) {
899 det=(AliModule*)dets[ind];
901 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
904 printf(" %6s: %3d -> %3d;","NULL",0,0);
910 //____________________________________________________________________________
911 void AliRun::SetGenerator(AliGenerator *generator)
914 // Load the event generator
916 if(!fGenerator) fGenerator = generator;
919 //____________________________________________________________________________
920 void AliRun::SetTransPar(char* filename)
923 // Read filename to set the transport parameters
927 const Int_t ncuts=10;
928 const Int_t nflags=11;
929 const Int_t npars=ncuts+nflags;
930 const char pars[npars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
931 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
932 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
933 "MULS","PAIR","PHOT","RAYL"};
939 Int_t i, itmed, iret, ktmed, kz;
942 // See whether the file is there
943 filtmp=gSystem->ExpandPathName(filename);
944 lun=fopen(filtmp,"r");
947 Warning("SetTransPar","File %s does not exist!\n",filename);
951 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
952 printf(" *%59s\n","*");
953 printf(" * Please check carefully what you are doing!%10s\n","*");
954 printf(" *%59s\n","*");
957 // Initialise cuts and flags
958 for(i=0;i<ncuts;i++) cut[i]=-99;
959 for(i=0;i<nflags;i++) flag[i]=-99;
961 for(i=0;i<256;i++) line[i]='\0';
962 // Read up to the end of line excluded
963 iret=fscanf(lun,"%[^\n]",line);
967 printf(" *%59s\n","*");
968 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
971 // Read the end of line
974 if(line[0]=='*') continue;
976 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",
977 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
978 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
979 &flag[8],&flag[9],&flag[10]);
983 Warning("SetTransPar","Error reading file %s\n",filename);
986 // Check that the module exist
987 AliModule *mod = GetModule(detName);
989 // Get the array of media numbers
990 TArrayI &idtmed = *mod->GetIdtmed();
991 // Check that the tracking medium code is valid
992 if(0<=itmed && itmed < 100) {
995 Warning("SetTransPar","Invalid tracking medium code %d for %s\n",itmed,mod->GetName());
998 // Set energy thresholds
999 for(kz=0;kz<ncuts;kz++) {
1001 printf(" * %-6s set to %10.3E for tracking medium code %4d for %s\n",
1002 pars[kz],cut[kz],itmed,mod->GetName());
1003 gMC->Gstpar(ktmed,pars[kz],cut[kz]);
1006 // Set transport mechanisms
1007 for(kz=0;kz<nflags;kz++) {
1009 printf(" * %-6s set to %10d for tracking medium code %4d for %s\n",
1010 pars[ncuts+kz],flag[kz],itmed,mod->GetName());
1011 gMC->Gstpar(ktmed,pars[ncuts+kz],Float_t(flag[kz]));
1015 Warning("SetTransPar","Invalid medium code %d *\n",itmed);
1019 Warning("SetTransPar","Module %s not present\n",detName);
1025 //_____________________________________________________________________________
1026 void AliRun::MakeTree(Option_t *option)
1029 // Create the ROOT trees
1030 // Loop on all detectors to create the Root branch (if any)
1035 char *K = strstr(option,"K");
1036 char *H = strstr(option,"H");
1037 char *E = strstr(option,"E");
1038 char *D = strstr(option,"D");
1039 char *R = strstr(option,"R");
1041 if (K && !fTreeK) fTreeK = new TTree("TreeK0","Kinematics");
1042 if (H && !fTreeH) fTreeH = new TTree("TreeH0","Hits");
1043 if (D && !fTreeD) fTreeD = new TTree("TreeD0","Digits");
1044 if (E && !fTreeE) fTreeE = new TTree("TE","Header");
1045 if (R && !fTreeR) fTreeR = new TTree("TreeR0","Reconstruction");
1046 if (fTreeH) fTreeH->SetAutoSave(1000000000); //no autosave
1048 // Create a branch for hits/digits for each detector
1049 // Each branch is a TClonesArray. Each data member of the Hits classes
1050 // will be in turn a subbranch of the detector master branch
1051 TIter next(fModules);
1052 AliModule *detector;
1053 while((detector = (AliModule*)next())) {
1054 if (H || D || R) detector->MakeBranch(option);
1056 // Create a branch for particles
1057 if (fTreeK && K) fTreeK->Branch("Particles",&fParticles,4000);
1059 // Create a branch for Header
1060 if (fTreeE && E) fTreeE->Branch("Header","AliHeader",&header,4000);
1063 //_____________________________________________________________________________
1064 Int_t AliRun::PurifyKine(Int_t lastSavedTrack, Int_t nofTracks)
1067 // PurifyKine with external parameters
1069 fHgwmk = lastSavedTrack;
1070 fNtrack = nofTracks;
1075 //_____________________________________________________________________________
1076 void AliRun::PurifyKine()
1079 // Compress kinematic tree keeping only flagged particles
1080 // and renaming the particle id's in all the hits
1082 TClonesArray &particles = *fParticles;
1083 int nkeep=fHgwmk+1, parent, i;
1084 TParticle *part, *partnew, *father;
1086 int *map = new int[particles.GetEntries()];
1088 // Save in Header total number of tracks before compression
1089 fHeader.SetNtrack(fHeader.GetNtrack()+fNtrack-fHgwmk);
1091 // Preset map, to be removed later
1092 for(i=0; i<fNtrack; i++) {
1093 if(i<=fHgwmk) map[i]=i ; else map[i] = -99 ;}
1094 // Second pass, build map between old and new numbering
1095 for(i=fHgwmk+1; i<fNtrack; i++) {
1096 part = (TParticle *)particles.UncheckedAt(i);
1097 if(part->TestBit(Keep_Bit)) {
1099 // This particle has to be kept
1103 // Old and new are different, have to copy
1104 partnew = (TParticle *)particles.UncheckedAt(nkeep);
1106 } else partnew = part;
1108 // as the parent is always *before*, it must be already
1109 // in place. This is what we are checking anyway!
1110 if((parent=partnew->GetFirstMother())>fHgwmk) {
1111 if(map[parent]==-99) printf("map[%d] = -99!\n",parent);
1112 partnew->SetFirstMother(map[parent]);
1119 // Fix daughters information
1120 for (i=0; i<fNtrack; i++) {
1121 part = (TParticle *)particles.UncheckedAt(i);
1122 parent = part->GetFirstMother();
1124 father = (TParticle *)particles.UncheckedAt(parent);
1125 if(father->TestBit(Daughters_Bit)) {
1127 if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
1128 if(i>father->GetLastDaughter()) father->SetLastDaughter(i);
1130 // Iitialise daughters info for first pass
1131 father->SetFirstDaughter(i);
1132 father->SetLastDaughter(i);
1133 father->SetBit(Daughters_Bit);
1138 // Now loop on all detectors and reset the hits
1139 TIter next(fModules);
1140 AliModule *detector;
1141 while((detector = (AliModule*)next())) {
1142 if (!detector->Hits()) continue;
1143 TClonesArray &vHits=*(detector->Hits());
1144 if(vHits.GetEntries() != detector->GetNhits())
1145 printf("vHits.GetEntries()!=detector->GetNhits(): %d != %d\n",
1146 vHits.GetEntries(),detector->GetNhits());
1147 for (i=0; i<detector->GetNhits(); i++) {
1148 OneHit = (AliHit *)vHits.UncheckedAt(i);
1149 OneHit->SetTrack(map[OneHit->GetTrack()]);
1154 particles.SetLast(fHgwmk);
1158 //_____________________________________________________________________________
1159 void AliRun::Reset(Int_t run, Int_t idevent)
1162 // Reset all Detectors & kinematics & trees
1170 // Initialise event header
1171 fHeader.Reset(run,idevent);
1175 sprintf(hname,"TreeK%d",idevent);
1176 fTreeK->SetName(hname);
1180 sprintf(hname,"TreeH%d",idevent);
1181 fTreeH->SetName(hname);
1185 sprintf(hname,"TreeD%d",idevent);
1186 fTreeD->SetName(hname);
1190 sprintf(hname,"TreeR%d",idevent);
1191 fTreeR->SetName(hname);
1195 //_____________________________________________________________________________
1196 void AliRun::ResetDigits()
1199 // Reset all Detectors digits
1201 TIter next(fModules);
1202 AliModule *detector;
1203 while((detector = (AliModule*)next())) {
1204 detector->ResetDigits();
1208 //_____________________________________________________________________________
1209 void AliRun::ResetHits()
1212 // Reset all Detectors hits
1214 TIter next(fModules);
1215 AliModule *detector;
1216 while((detector = (AliModule*)next())) {
1217 detector->ResetHits();
1221 //_____________________________________________________________________________
1222 void AliRun::ResetPoints()
1225 // Reset all Detectors points
1227 TIter next(fModules);
1228 AliModule *detector;
1229 while((detector = (AliModule*)next())) {
1230 detector->ResetPoints();
1234 //_____________________________________________________________________________
1235 void AliRun::Run(Int_t nevent, const char *setup)
1238 // Main function to be called to process a galice run
1240 // Root > gAlice.Run();
1241 // a positive number of events will cause the finish routine
1246 // check if initialisation has been done
1247 if (!fInitDone) Init(setup);
1249 // Create the Root Tree with one branch per detector
1251 gAlice->MakeTree("KHDER");
1254 todo = TMath::Abs(nevent);
1255 for (i=0; i<todo; i++) {
1256 // Process one run (one run = one event)
1257 gAlice->Reset(fRun, fEvent);
1261 gAlice->FinishEvent();
1265 // End of this run, close files
1266 if(nevent>0) gAlice->FinishRun();
1269 //_____________________________________________________________________________
1270 void AliRun::RunLego(const char *setup,Int_t ntheta,Float_t themin,
1271 Float_t themax,Int_t nphi,Float_t phimin,Float_t phimax,
1272 Float_t rmin,Float_t rmax,Float_t zmax)
1275 // Generates lego plots of:
1276 // - radiation length map phi vs theta
1277 // - radiation length map phi vs eta
1278 // - interaction length map
1279 // - g/cm2 length map
1281 // ntheta bins in theta, eta
1282 // themin minimum angle in theta (degrees)
1283 // themax maximum angle in theta (degrees)
1285 // phimin minimum angle in phi (degrees)
1286 // phimax maximum angle in phi (degrees)
1287 // rmin minimum radius
1288 // rmax maximum radius
1291 // The number of events generated = ntheta*nphi
1292 // run input parameters in macro setup (default="Config.C")
1294 // Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
1297 <img src="picts/AliRunLego1.gif">
1302 <img src="picts/AliRunLego2.gif">
1307 <img src="picts/AliRunLego3.gif">
1312 // check if initialisation has been done
1313 if (!fInitDone) Init(setup);
1315 fLego = new AliLego("lego","lego");
1316 fLego->Init(ntheta,themin,themax,nphi,phimin,phimax,rmin,rmax,zmax);
1319 // Create only the Root event Tree
1320 gAlice->MakeTree("E");
1322 // End of this run, close files
1323 gAlice->FinishRun();
1326 //_____________________________________________________________________________
1327 void AliRun::SetCurrentTrack(Int_t track)
1330 // Set current track number
1335 //_____________________________________________________________________________
1336 void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
1337 Float_t *vpos, Float_t *polar, Float_t tof,
1338 const char *mecha, Int_t &ntr, Float_t weight)
1341 // Load a track on the stack
1343 // done 0 if the track has to be transported
1345 // parent identifier of the parent track. -1 for a primary
1346 // pdg particle code
1347 // pmom momentum GeV/c
1349 // polar polarisation
1350 // tof time of flight in seconds
1351 // mecha production mechanism
1352 // ntr on output the number of the track stored
1354 TClonesArray &particles = *fParticles;
1355 TParticle *particle;
1357 const Int_t firstdaughter=-1;
1358 const Int_t lastdaughter=-1;
1360 // const Float_t tlife=0;
1363 // Here we get the static mass
1364 // For MC is ok, but a more sophisticated method could be necessary
1365 // if the calculated mass is required
1366 // also, this method is potentially dangerous if the mass
1367 // used in the MC is not the same of the PDG database
1369 mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
1370 Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
1371 pmom[1]*pmom[1]+pmom[2]*pmom[2]);
1373 //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",
1374 //pname,mass,e,fNtrack,pdg,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],KS,mecha);
1376 particle=new(particles[fNtrack]) TParticle(pdg,KS,parent,-1,firstdaughter,
1377 lastdaughter,pmom[0],pmom[1],pmom[2],
1378 e,vpos[0],vpos[1],vpos[2],tof);
1379 // polar[0],polar[1],polar[2],tof,
1381 ((TParticle*)particles[fNtrack])->SetPolarisation(TVector3(polar[0],polar[1],polar[2]));
1382 ((TParticle*)particles[fNtrack])->SetWeight(weight);
1383 if(!done) particle->SetBit(Done_Bit);
1386 particle=(TParticle*) fParticles->UncheckedAt(parent);
1387 particle->SetLastDaughter(fNtrack);
1388 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
1391 // This is a primary track. Set high water mark for this event
1394 // Set also number if primary tracks
1395 fHeader.SetNprimary(fHgwmk+1);
1396 fHeader.SetNtrack(fHgwmk+1);
1401 //_____________________________________________________________________________
1402 void AliRun::KeepTrack(const Int_t track)
1405 // flags a track to be kept
1407 TClonesArray &particles = *fParticles;
1408 ((TParticle*)particles[track])->SetBit(Keep_Bit);
1411 //_____________________________________________________________________________
1412 void AliRun::StepManager(Int_t id) const
1415 // Called at every step during transport
1420 // --- If lego option, do it and leave
1422 fLego->StepManager();
1425 //Update energy deposition tables
1426 sEventEnergy[gMC->CurrentVolID(copy)]+=gMC->Edep();
1428 //Call the appropriate stepping routine;
1429 AliModule *det = (AliModule*)fModules->At(id);
1430 if(det) det->StepManager();
1433 //_____________________________________________________________________________
1434 void AliRun::ReadEuclid(const char* filnam, const AliModule *det, char* topvol)
1437 // read in the geometry of the detector in euclid file format
1439 // id_det : the detector identification (2=its,...)
1440 // topvol : return parameter describing the name of the top
1441 // volume of geometry.
1443 // author : m. maire
1446 // several changes have been made by miroslav helbich
1447 // subroutine is rewrited to follow the new established way of memory
1448 // booking for tracking medias and rotation matrices.
1449 // all used tracking media have to be defined first, for this you can use
1450 // subroutine greutmed.
1451 // top volume is searched as only volume not positioned into another
1454 Int_t i, nvol, iret, itmed, irot, numed, npar, ndiv, iaxe;
1455 Int_t ndvmx, nr, flag;
1456 char key[5], card[77], natmed[21];
1457 char name[5], mother[5], shape[5], konly[5], volst[7000][5];
1460 Float_t teta1, phi1, teta2, phi2, teta3, phi3, orig, step;
1462 const Int_t maxrot=5000;
1463 Int_t idrot[maxrot],istop[7000];
1466 // *** The input filnam name will be with extension '.euc'
1467 filtmp=gSystem->ExpandPathName(filnam);
1468 lun=fopen(filtmp,"r");
1471 Error("ReadEuclid","Could not open file %s\n",filnam);
1474 //* --- definition of rotation matrix 0 ---
1475 TArrayI &idtmed = *(det->GetIdtmed());
1476 for(i=1; i<maxrot; ++i) idrot[i]=-99;
1480 for(i=0;i<77;i++) card[i]=0;
1481 iret=fscanf(lun,"%77[^\n]",card);
1482 if(iret<=0) goto L20;
1485 strncpy(key,card,4);
1487 if (!strcmp(key,"TMED")) {
1488 sscanf(&card[5],"%d '%[^']'",&itmed,natmed);
1489 if( itmed<0 || itmed>=100 ) {
1490 Error("ReadEuclid","TMED illegal medium number %d for %s\n",itmed,natmed);
1493 //Pad the string with blanks
1496 while(i<20) natmed[i++]=' ';
1499 if( idtmed[itmed]<=0 ) {
1500 Error("ReadEuclid","TMED undefined medium number %d for %s\n",itmed,natmed);
1503 gMC->Gckmat(idtmed[itmed],natmed);
1505 } else if (!strcmp(key,"ROTM")) {
1506 sscanf(&card[4],"%d %f %f %f %f %f %f",&irot,&teta1,&phi1,&teta2,&phi2,&teta3,&phi3);
1507 if( irot<=0 || irot>=maxrot ) {
1508 Error("ReadEuclid","ROTM rotation matrix number %d illegal\n",irot);
1511 det->AliMatrix(idrot[irot],teta1,phi1,teta2,phi2,teta3,phi3);
1513 } else if (!strcmp(key,"VOLU")) {
1514 sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, shape, &numed, &npar);
1516 for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
1519 gMC->Gsvolu( name, shape, idtmed[numed], par, npar);
1520 //* save the defined volumes
1521 strcpy(volst[++nvol],name);
1524 } else if (!strcmp(key,"DIVN")) {
1525 sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, mother, &ndiv, &iaxe);
1526 gMC->Gsdvn ( name, mother, ndiv, iaxe );
1528 } else if (!strcmp(key,"DVN2")) {
1529 sscanf(&card[5],"'%[^']' '%[^']' %d %d %f %d",name, mother, &ndiv, &iaxe, &orig, &numed);
1530 gMC->Gsdvn2( name, mother, ndiv, iaxe, orig,idtmed[numed]);
1532 } else if (!strcmp(key,"DIVT")) {
1533 sscanf(&card[5],"'%[^']' '%[^']' %f %d %d %d", name, mother, &step, &iaxe, &numed, &ndvmx);
1534 gMC->Gsdvt ( name, mother, step, iaxe, idtmed[numed], ndvmx);
1536 } else if (!strcmp(key,"DVT2")) {
1537 sscanf(&card[5],"'%[^']' '%[^']' %f %d %f %d %d", name, mother, &step, &iaxe, &orig, &numed, &ndvmx);
1538 gMC->Gsdvt2 ( name, mother, step, iaxe, orig, idtmed[numed], ndvmx );
1540 } else if (!strcmp(key,"POSI")) {
1541 sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']'", name, &nr, mother, &xo, &yo, &zo, &irot, konly);
1542 if( irot<0 || irot>=maxrot ) {
1543 Error("ReadEuclid","POSI %s#%d rotation matrix number %d illegal\n",name,nr,irot);
1546 if( idrot[irot] == -99) {
1547 Error("ReadEuclid","POSI %s#%d undefined matrix number %d\n",name,nr,irot);
1550 //*** volume name cannot be the top volume
1551 for(i=1;i<=nvol;i++) {
1552 if (!strcmp(volst[i],name)) istop[i]=0;
1555 gMC->Gspos ( name, nr, mother, xo, yo, zo, idrot[irot], konly );
1557 } else if (!strcmp(key,"POSP")) {
1558 sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']' %d", name, &nr, mother, &xo, &yo, &zo, &irot, konly, &npar);
1559 if( irot<0 || irot>=maxrot ) {
1560 Error("ReadEuclid","POSP %s#%d rotation matrix number %d illegal\n",name,nr,irot);
1563 if( idrot[irot] == -99) {
1564 Error("ReadEuclid","POSP %s#%d undefined matrix number %d\n",name,nr,irot);
1568 for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
1571 //*** volume name cannot be the top volume
1572 for(i=1;i<=nvol;i++) {
1573 if (!strcmp(volst[i],name)) istop[i]=0;
1576 gMC->Gsposp ( name, nr, mother, xo,yo,zo, idrot[irot], konly, par, npar);
1579 if (strcmp(key,"END")) goto L10;
1580 //* find top volume in the geometry
1582 for(i=1;i<=nvol;i++) {
1583 if (istop[i] && flag) {
1584 Warning("ReadEuclid"," %s is another possible top volume\n",volst[i]);
1586 if (istop[i] && !flag) {
1587 strcpy(topvol,volst[i]);
1588 printf(" *** GREUCL *** volume %s taken as a top volume\n",topvol);
1593 Warning("ReadEuclid","top volume not found\n");
1597 //* commented out only for the not cernlib version
1598 printf(" *** GREUCL *** file: %s is now read in\n",filnam);
1603 Error("ReadEuclid","reading error or premature end of file\n");
1606 //_____________________________________________________________________________
1607 void AliRun::ReadEuclidMedia(const char* filnam, const AliModule *det)
1610 // read in the materials and tracking media for the detector
1611 // in euclid file format
1613 // filnam: name of the input file
1614 // id_det: id_det is the detector identification (2=its,...)
1616 // author : miroslav helbich
1618 Float_t sxmgmx = gAlice->Field()->Max();
1619 Int_t isxfld = gAlice->Field()->Integ();
1620 Int_t end, i, iret, itmed;
1621 char key[5], card[130], natmed[21], namate[21];
1626 Int_t nwbuf, isvol, ifield, nmat;
1627 Float_t a, z, dens, radl, absl, fieldm, tmaxfd, stemax, deemax, epsil, stmin;
1630 for(i=0;i<end;i++) if(filnam[i]=='.') {
1635 // *** The input filnam name will be with extension '.euc'
1636 printf("The file name is %s\n",filnam); //Debug
1637 filtmp=gSystem->ExpandPathName(filnam);
1638 lun=fopen(filtmp,"r");
1641 Warning("ReadEuclidMedia","Could not open file %s\n",filnam);
1645 // Retrieve Mag Field parameters
1646 Int_t ISXFLD=gAlice->Field()->Integ();
1647 Float_t SXMGMX=gAlice->Field()->Max();
1648 // TArrayI &idtmed = *(det->GetIdtmed());
1651 for(i=0;i<130;i++) card[i]=0;
1652 iret=fscanf(lun,"%4s %[^\n]",key,card);
1653 if(iret<=0) goto L20;
1657 if (!strcmp(key,"MATE")) {
1658 sscanf(card,"%d '%[^']' %f %f %f %f %f %d",&imate,namate,&a,&z,&dens,&radl,&absl,&nwbuf);
1659 if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
1660 //Pad the string with blanks
1663 while(i<20) namate[i++]=' ';
1666 det->AliMaterial(imate,namate,a,z,dens,radl,absl,ubuf,nwbuf);
1667 //* read tracking medium
1668 } else if (!strcmp(key,"TMED")) {
1669 sscanf(card,"%d '%[^']' %d %d %d %f %f %f %f %f %f %d",
1670 &itmed,natmed,&nmat,&isvol,&ifield,&fieldm,&tmaxfd,
1671 &stemax,&deemax,&epsil,&stmin,&nwbuf);
1672 if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
1673 if (ifield<0) ifield=isxfld;
1674 if (fieldm<0) fieldm=sxmgmx;
1675 //Pad the string with blanks
1678 while(i<20) natmed[i++]=' ';
1681 det->AliMedium(itmed,natmed,nmat,isvol,ISXFLD,SXMGMX,tmaxfd,
1682 stemax,deemax,epsil,stmin,ubuf,nwbuf);
1683 // (*fImedia)[idtmed[itmed]-1]=id_det;
1687 if (strcmp(key,"END")) goto L10;
1690 //* commented out only for the not cernlib version
1691 Warning("ReadEuclidMedia","file: %s is now read in\n",filnam);
1696 Warning("ReadEuclidMedia","reading error or premature end of file\n");
1699 //_____________________________________________________________________________
1700 void AliRun::Streamer(TBuffer &R__b)
1703 // Stream an object of class AliRun.
1705 if (R__b.IsReading()) {
1706 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
1707 TNamed::Streamer(R__b);
1708 if (!gAlice) gAlice = this;
1709 gROOT->GetListOfBrowsables()->Add(this,"Run");
1710 fTreeE = (TTree*)gDirectory->Get("TE");
1711 if (fTreeE) fTreeE->SetBranchAddress("Header", &header);
1712 else Error("Streamer","cannot find Header Tree\n");
1716 fHeader.Streamer(R__b);
1726 R__b >> fPDGDB; //Particle factory object!
1727 fTreeE->GetEntry(0);
1729 fHeader.SetEvent(0);
1730 fPDGDB = TDatabasePDG::Instance(); //Particle factory object!
1733 R__b.WriteVersion(AliRun::IsA());
1734 TNamed::Streamer(R__b);
1738 fHeader.Streamer(R__b);
1747 R__b << fPDGDB; //Particle factory object!
1752 //_____________________________________________________________________________
1754 // Interfaces to Fortran
1756 //_____________________________________________________________________________
1758 extern "C" void type_of_call rxgtrak (Int_t &mtrack, Int_t &ipart, Float_t *pmom,
1759 Float_t &e, Float_t *vpos, Float_t *polar,
1763 // Fetches next track from the ROOT stack for transport. Called by the
1764 // modified version of GTREVE.
1766 // Track number in the ROOT stack. If MTRACK=0 no
1767 // mtrack more tracks are left in the stack to be
1769 // ipart Particle code in the GEANT conventions.
1770 // pmom[3] Particle momentum in GeV/c
1771 // e Particle energy in GeV
1772 // vpos[3] Particle position
1773 // tof Particle time of flight in seconds
1776 gAlice->GetNextTrack(mtrack, pdg, pmom, e, vpos, polar, tof);
1777 ipart = gMC->IdFromPDG(pdg);
1781 //_____________________________________________________________________________
1782 extern "C" void type_of_call
1784 rxstrak (Int_t &keep, Int_t &parent, Int_t &ipart, Float_t *pmom,
1785 Float_t *vpos, Float_t &tof, const char* cmech, Int_t &ntr, const int cmlen)
1787 rxstrak (Int_t &keep, Int_t &parent, Int_t &ipart, Float_t *pmom,
1788 Float_t *vpos, Float_t &tof, const char* cmech, const int cmlen,
1793 // Fetches next track from the ROOT stack for transport. Called by GUKINE
1796 // Status of the track. If keep=0 the track is put
1797 // keep on the ROOT stack but it is not fetched for
1799 // parent Parent track. If parent=0 the track is a primary.
1800 // In GUSTEP the routine is normally called to store
1801 // secondaries generated by the current track whose
1802 // ROOT stack number is MTRACK (common SCKINE.
1803 // ipart Particle code in the GEANT conventions.
1804 // pmom[3] Particle momentum in GeV/c
1805 // vpos[3] Particle position
1806 // tof Particle time of flight in seconds
1808 // cmech (CHARACTER*10) Particle origin. This field is user
1809 // defined and it is not used inside the GALICE code.
1810 // ntr Number assigned to the particle in the ROOT stack.
1813 Float_t polar[3]={0.,0.,0.};
1814 for(int i=0; i<10 && i<cmlen; i++) mecha[i]=cmech[i];
1816 Int_t pdg=gMC->PDGFromId(ipart);
1817 gAlice->SetTrack(keep, parent-1, pdg, pmom, vpos, polar, tof, mecha, ntr);
1821 //_____________________________________________________________________________
1822 extern "C" void type_of_call rxkeep(const Int_t &n)
1824 if( NULL==gAlice ) exit(1);
1826 if( n<=0 || n>gAlice->Particles()->GetEntries() )
1828 printf(" Bad index n=%d must be 0<n<=%d\n",
1829 n,gAlice->Particles()->GetEntries());
1833 ((TParticle*)(gAlice->Particles()->UncheckedAt(n-1)))->SetBit(Keep_Bit);
1836 //_____________________________________________________________________________
1837 extern "C" void type_of_call rxouth ()
1840 // Called by Gtreve at the end of each primary track
1842 gAlice->FinishPrimary();