]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliRun.cxx
Correct problem when getting geometry
[u/mrichter/AliRoot.git] / STEER / AliRun.cxx
1 ///////////////////////////////////////////////////////////////////////////////
2 //                                                                           //
3 //  Control class for Alice C++                                              //
4 //  Only one single instance of this class exists.                           //
5 //  The object is created in main program aliroot                            //
6 //  and is pointed by the global gAlice.                                     //
7 //                                                                           //
8 //   -Supports the list of all Alice Detectors (fModules).                 //
9 //   -Supports the list of particles (fParticles).                           //
10 //   -Supports the Trees.                                                    //
11 //   -Supports the geometry.                                                 //
12 //   -Supports the event display.                                            //
13 //Begin_Html
14 /*
15 <img src="picts/AliRunClass.gif">
16 */
17 //End_Html
18 //Begin_Html
19 /*
20 <img src="picts/alirun.gif">
21 */
22 //End_Html
23 //                                                                           //
24 ///////////////////////////////////////////////////////////////////////////////
25
26 #include <TFile.h>
27 #include <TRandom.h>
28 #include <TBRIK.h> 
29 #include <TNode.h> 
30 #include <TCint.h> 
31 #include <TSystem.h>
32 #include <TObjectTable.h>
33
34 #include "TParticle.h"
35 #include "AliRun.h"
36 #include "AliDisplay.h"
37
38 #include "AliCallf77.h" 
39  
40 #include <stdlib.h>
41 #include <stdio.h>
42 #include <string.h>
43  
44 AliRun *gAlice;
45
46 static AliHeader *header;
47
48 #ifndef WIN32
49
50 # define rxgtrak rxgtrak_
51 # define rxstrak rxstrak_
52 # define rxkeep  rxkeep_ 
53 # define rxouth  rxouth_
54 #else
55
56 # define rxgtrak RXGTRAK 
57 # define rxstrak RXSTRAK 
58 # define rxkeep  RXKEEP  
59 # define rxouth  RXOUTH
60 #endif
61
62 static TArrayF sEventEnergy;
63 static TArrayF sSummEnergy;
64 static TArrayF sSum2Energy;
65
66 ClassImp(AliRun)
67
68 //_____________________________________________________________________________
69 AliRun::AliRun()
70 {
71   //
72   // Default constructor for AliRun
73   //
74   header=&fHeader;
75   fRun       = 0;
76   fEvent     = 0;
77   fCurrent   = -1;
78   fModules = 0;
79   fGenerator = 0;
80   fTreeD     = 0;
81   fTreeK     = 0;
82   fTreeH     = 0;
83   fTreeE     = 0;
84   fTreeR     = 0;
85   fParticles = 0;
86   fGeometry  = 0;
87   fDisplay   = 0;
88   fField     = 0;
89   fMC        = 0;
90   fNdets     = 0;
91   fImedia    = 0;
92   fTrRmax    = 1.e10;
93   fTrZmax    = 1.e10;
94   fInitDone  = kFALSE;
95   fLego      = 0;
96   fPDGDB     = 0;        //Particle factory object!
97 }
98
99 //_____________________________________________________________________________
100 AliRun::AliRun(const char *name, const char *title)
101   : TNamed(name,title)
102 {
103   //
104   //  Constructor for the main processor.
105   //  Creates the geometry
106   //  Creates the list of Detectors.
107   //  Creates the list of particles.
108   //
109   Int_t i;
110   
111   gAlice     = this;
112   fTreeD     = 0;
113   fTreeK     = 0;
114   fTreeH     = 0;
115   fTreeE     = 0;
116   fTreeR     = 0;
117   fTrRmax    = 1.e10;
118   fTrZmax    = 1.e10;
119   fGenerator = 0;
120   fInitDone  = kFALSE;
121   fLego      = 0;
122   fField     = 0;
123   
124   gROOT->GetListOfBrowsables()->Add(this,name);
125   //
126   // create the support list for the various Detectors
127   fModules = new TObjArray(77);
128   //
129   // Create the TNode geometry for the event display
130   
131   BuildSimpleGeometry();
132   
133
134   fNtrack=0;
135   fHgwmk=0;
136   fCurrent=-1;
137   header=&fHeader;
138   fRun       = 0;
139   fEvent     = 0;
140   //
141   // Create the particle stack
142   fParticles = new TClonesArray("TParticle",100);
143   
144   fDisplay = 0;
145   //
146   // Create default mag field
147   SetField();
148   //
149   fMC      = gMC;
150   //
151   // Prepare the tracking medium lists
152   fImedia = new TArrayI(1000);
153   for(i=0;i<1000;i++) (*fImedia)[i]=-99;
154   //
155   // Make particles
156   fPDGDB     = TDatabasePDG::Instance();        //Particle factory object!
157 }
158
159 //_____________________________________________________________________________
160 AliRun::~AliRun()
161 {
162   //
163   // Defaullt AliRun destructor
164   //
165   delete fImedia;
166   delete fField;
167   delete fMC;
168   delete fGeometry;
169   delete fDisplay;
170   delete fGenerator;
171   delete fLego;
172   delete fTreeD;
173   delete fTreeK;
174   delete fTreeH;
175   delete fTreeE;
176   delete fTreeR;
177   if (fModules) {
178     fModules->Delete();
179     delete fModules;
180   }
181   if (fParticles) {
182     fParticles->Delete();
183     delete fParticles;
184   }
185 }
186
187 //_____________________________________________________________________________
188 void AliRun::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
189 {
190   //
191   //  Add a hit to detector id
192   //
193   TObjArray &dets = *fModules;
194   if(dets[id]) ((AliModule*) dets[id])->AddHit(track,vol,hits);
195 }
196
197 //_____________________________________________________________________________
198 void AliRun::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
199 {
200   //
201   // Add digit to detector id
202   //
203   TObjArray &dets = *fModules;
204   if(dets[id]) ((AliModule*) dets[id])->AddDigit(tracks,digits);
205 }
206
207 //_____________________________________________________________________________
208 void AliRun::Browse(TBrowser *b)
209 {
210   //
211   // Called when the item "Run" is clicked on the left pane
212   // of the Root browser.
213   // It displays the Root Trees and all detectors.
214   //
215   if (fTreeK) b->Add(fTreeK,fTreeK->GetName());
216   if (fTreeH) b->Add(fTreeH,fTreeH->GetName());
217   if (fTreeD) b->Add(fTreeD,fTreeD->GetName());
218   if (fTreeE) b->Add(fTreeE,fTreeE->GetName());
219   if (fTreeR) b->Add(fTreeR,fTreeR->GetName());
220   
221   TIter next(fModules);
222   AliModule *detector;
223   while((detector = (AliModule*)next())) {
224     b->Add(detector,detector->GetName());
225   }
226 }
227
228 //_____________________________________________________________________________
229 void AliRun::Build()
230 {
231   //
232   // Initialize Alice geometry
233   // Dummy routine
234   //
235 }
236  
237 //_____________________________________________________________________________
238 void AliRun::BuildSimpleGeometry()
239 {
240   //
241   // Create a simple TNode geometry used by Root display engine
242   //
243   // Initialise geometry
244   //
245   fGeometry = new TGeometry("AliceGeom","Galice Geometry for Hits");
246   new TMaterial("void","Vacuum",0,0,0);  //Everything is void
247   TBRIK *brik = new TBRIK("S_alice","alice volume","void",2000,2000,3000);
248   brik->SetVisibility(0);
249   new TNode("alice","alice","S_alice");
250 }
251
252 //_____________________________________________________________________________
253 void AliRun::CleanDetectors()
254 {
255   //
256   // Clean Detectors at the end of event
257   //
258   TIter next(fModules);
259   AliModule *detector;
260   while((detector = (AliModule*)next())) {
261     detector->FinishEvent();
262   }
263 }
264
265 //_____________________________________________________________________________
266 void AliRun::CleanParents()
267 {
268   //
269   // Clean Particles stack.
270   // Set parent/daughter relations
271   //
272   TClonesArray &particles = *(gAlice->Particles());
273   TParticle *part;
274   int i;
275   for(i=0; i<fNtrack; i++) {
276     part = (TParticle *)particles.UncheckedAt(i);
277     if(!part->TestBit(Daughters_Bit)) {
278       part->SetFirstDaughter(-1);
279       part->SetLastDaughter(-1);
280     }
281   }
282 }
283
284 //_____________________________________________________________________________
285 Int_t AliRun::DistancetoPrimitive(Int_t, Int_t)
286 {
287   //
288   // Return the distance from the mouse to the AliRun object
289   // Dummy routine
290   //
291   return 9999;
292 }
293
294 //_____________________________________________________________________________
295 void AliRun::DumpPart (Int_t i)
296 {
297   //
298   // Dumps particle i in the stack
299   //
300   TClonesArray &particles = *fParticles;
301   ((TParticle*) particles[i])->Print();
302 }
303
304 //_____________________________________________________________________________
305 void AliRun::DumpPStack ()
306 {
307   //
308   // Dumps the particle stack
309   //
310   TClonesArray &particles = *fParticles;
311   printf(
312          "\n\n=======================================================================\n");
313   for (Int_t i=0;i<fNtrack;i++) 
314     {
315       printf("-> %d ",i); ((TParticle*) particles[i])->Print();
316       printf("--------------------------------------------------------------\n");
317     }
318   printf(
319          "\n=======================================================================\n\n");
320 }
321
322 //_____________________________________________________________________________
323 void AliRun::SetField(Int_t type, Int_t version, Float_t scale,
324                       Float_t maxField, char* filename)
325 {
326   //
327   //  Set magnetic field parameters
328   //  type      Magnetic field transport flag 0=no field, 2=helix, 3=Runge Kutta
329   //  version   Magnetic field map version (only 1 active now)
330   //  scale     Scale factor for the magnetic field
331   //  maxField  Maximum value for the magnetic field
332
333   //
334   // --- Sanity check on mag field flags
335   if(type<0 || type > 2) {
336     Warning("SetField",
337             "Invalid magnetic field flag: %5d; Helix tracking chosen instead\n"
338            ,type);
339     type=2;
340   }
341   if(fField) delete fField;
342   if(version==1) {
343     fField = new AliMagFC("Map1"," ",type,version,scale,maxField);
344   } else if(version<=3) {
345     fField = new AliMagFCM("Map2-3",filename,type,version,scale,maxField);
346     fField->ReadField();
347   } else {
348     Warning("SetField","Invalid map %d\n",version);
349   }
350 }
351
352 //_____________________________________________________________________________
353 void AliRun::FillTree()
354 {
355   //
356   // Fills all AliRun TTrees
357   //
358   if (fTreeK) fTreeK->Fill();
359   if (fTreeH) fTreeH->Fill();
360   if (fTreeD) fTreeD->Fill();
361   if (fTreeR) fTreeR->Fill();
362 }
363  
364 //_____________________________________________________________________________
365 void AliRun::FinishPrimary()
366 {
367   //
368   // Called  at the end of each primary track
369   //
370   
371   //  static Int_t count=0;
372   //  const Int_t times=10;
373   // This primary is finished, purify stack
374   gAlice->PurifyKine();
375
376   // Write out hits if any
377   if (gAlice->TreeH()) {
378     gAlice->TreeH()->Fill();
379   }
380   
381   // Reset Hits info
382   gAlice->ResetHits();
383
384   //
385   //  if(++count%times==1) gObjectTable->Print();
386 }
387
388 //_____________________________________________________________________________
389 void AliRun::FinishEvent()
390 {
391   //
392   // Called at the end of the event.
393   //
394   
395   //Update the energy deposit tables
396   Int_t i;
397   for(i=0;i<sEventEnergy.GetSize();i++) {
398     sSummEnergy[i]+=sEventEnergy[i];
399     sSum2Energy[i]+=sEventEnergy[i]*sEventEnergy[i];
400   }
401   sEventEnergy.Reset();
402   
403   // Clean detector information
404   CleanDetectors();
405   
406   // Write out the kinematics
407   if (fTreeK) {
408     CleanParents();
409     fTreeK->Fill();
410   }
411   
412   // Write out the digits
413   if (fTreeD) {
414     fTreeD->Fill();
415     ResetDigits();
416   }
417   
418   // Write out reconstructed clusters  
419   if (fTreeR) {
420     fTreeR->Fill();
421   }
422
423   // Write out the event Header information
424   if (fTreeE) fTreeE->Fill();
425   
426   // Reset stack info
427   ResetStack();
428   
429   // Write Tree headers
430   //  Int_t ievent = fHeader.GetEvent();
431   //  char hname[30];
432   //  sprintf(hname,"TreeK%d",ievent);
433   if (fTreeK) fTreeK->Write();
434   //  sprintf(hname,"TreeH%d",ievent);
435   if (fTreeH) fTreeH->Write();
436   //  sprintf(hname,"TreeD%d",ievent);
437   if (fTreeD) fTreeD->Write();
438   //  sprintf(hname,"TreeR%d",ievent);
439   if (fTreeR) fTreeR->Write();
440 }
441
442 //_____________________________________________________________________________
443 void AliRun::FinishRun()
444 {
445   //
446   // Called at the end of the run.
447   //
448
449   // Clean detector information
450   TIter next(fModules);
451   AliModule *detector;
452   while((detector = (AliModule*)next())) {
453     detector->FinishRun();
454   }
455   
456   //Output energy summary tables
457   EnergySummary();
458   
459   // file is retrieved from whatever tree
460   TFile *File = 0;
461   if (fTreeK) File = fTreeK->GetCurrentFile();
462   if ((!File) && (fTreeH)) File = fTreeH->GetCurrentFile();
463   if ((!File) && (fTreeD)) File = fTreeD->GetCurrentFile();
464   if ((!File) && (fTreeE)) File = fTreeE->GetCurrentFile();
465   if( NULL==File ) {
466     Error("FinishRun","There isn't root file!");
467     exit(1);
468   }
469   File->cd();
470   fTreeE->Write();
471   
472   // Clean tree information
473   delete fTreeK; fTreeK = 0;
474   delete fTreeH; fTreeH = 0;
475   delete fTreeD; fTreeD = 0;
476   delete fTreeR; fTreeR = 0;
477   delete fTreeE; fTreeE = 0;
478   
479   // Write AliRun info and all detectors parameters
480   Write();
481   
482   // Close output file
483   File->Write();
484   File->Close();
485 }
486
487 //_____________________________________________________________________________
488 void AliRun::FlagTrack(Int_t track)
489 {
490   //
491   // Flags a track and all its family tree to be kept
492   //
493   int curr;
494   TParticle *particle;
495
496   curr=track;
497   while(1) {
498     particle=(TParticle*)fParticles->UncheckedAt(curr);
499     
500     // If the particle is flagged the three from here upward is saved already
501     if(particle->TestBit(Keep_Bit)) return;
502     
503     // Save this particle
504     particle->SetBit(Keep_Bit);
505     
506     // Move to father if any
507     if((curr=particle->GetFirstMother())==-1) return;
508   }
509 }
510  
511 //_____________________________________________________________________________
512 void AliRun::EnergySummary()
513 {
514   //
515   // Print summary of deposited energy
516   //
517
518   Int_t ndep=0;
519   Float_t edtot=0;
520   Float_t ed, ed2;
521   Int_t kn, i, left, j, id;
522   const Float_t zero=0;
523   Int_t ievent=fHeader.GetEvent()+1;
524   //
525   // Energy loss information
526   if(ievent) {
527     printf("***************** Energy Loss Information per event (GEV) *****************\n");
528     for(kn=1;kn<sEventEnergy.GetSize();kn++) {
529       ed=sSummEnergy[kn];
530       if(ed>0) {
531         sEventEnergy[ndep]=kn;
532         if(ievent>1) {
533           ed=ed/ievent;
534           ed2=sSum2Energy[kn];
535           ed2=ed2/ievent;
536           ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,zero))/ed;
537         } else 
538           ed2=99;
539         sSummEnergy[ndep]=ed;
540         sSum2Energy[ndep]=TMath::Min((Float_t) 99.,TMath::Max(ed2,zero));
541         edtot+=ed;
542         ndep++;
543       }
544     }
545     for(kn=0;kn<(ndep-1)/3+1;kn++) {
546       left=ndep-kn*3;
547       for(i=0;i<(3<left?3:left);i++) {
548         j=kn*3+i;
549         id=Int_t (sEventEnergy[j]+0.1);
550         printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),sSummEnergy[j],sSum2Energy[j]);
551       }
552       printf("\n");
553     }
554     //
555     // Relative energy loss in different detectors
556     printf("******************** Relative Energy Loss per event ********************\n");
557     printf("Total energy loss per event %10.3f GeV\n",edtot);
558     for(kn=0;kn<(ndep-1)/5+1;kn++) {
559       left=ndep-kn*5;
560       for(i=0;i<(5<left?5:left);i++) {
561         j=kn*5+i;
562         id=Int_t (sEventEnergy[j]+0.1);
563         printf(" %s %10.3f%%;",gMC->VolName(id),100*sSummEnergy[j]/edtot);
564       }
565       printf("\n");
566     }
567     for(kn=0;kn<75;kn++) printf("*"); 
568     printf("\n");
569   }
570   //
571   // Reset the TArray's
572   sEventEnergy.Set(0);
573   sSummEnergy.Set(0);
574   sSum2Energy.Set(0);
575 }
576
577 //_____________________________________________________________________________
578 AliModule *AliRun::GetModule(const char *name)
579 {
580   //
581   // Return pointer to detector from name
582   //
583   return (AliModule*)fModules->FindObject(name);
584 }
585  
586 //_____________________________________________________________________________
587 AliDetector *AliRun::GetDetector(const char *name)
588 {
589   //
590   // Return pointer to detector from name
591   //
592   return (AliDetector*)fModules->FindObject(name);
593 }
594  
595 //_____________________________________________________________________________
596 Int_t AliRun::GetModuleID(const char *name)
597 {
598   //
599   // Return galice internal detector identifier from name
600   //
601   Int_t i=-1;
602   TObject *mod=fModules->FindObject(name);
603   if(mod) i=fModules->IndexOf(mod);
604   return i;
605 }
606  
607 //_____________________________________________________________________________
608 Int_t AliRun::GetEvent(Int_t event)
609 {
610   //
611   // Connect the Trees Kinematics and Hits for event # event
612   // Set branch addresses
613   //
614
615   // Reset existing structures
616   ResetStack();
617   ResetHits();
618   ResetDigits();
619   
620   // Delete Trees already connected
621   if (fTreeK) delete fTreeK;
622   if (fTreeH) delete fTreeH;
623   if (fTreeD) delete fTreeD;
624   if (fTreeR) delete fTreeR;
625
626   // Get header from file
627   if(fTreeE) fTreeE->GetEntry(event);
628   else Error("GetEvent","Cannot file Header Tree\n");
629   
630   // Get Kine Tree from file
631   char treeName[20];
632   sprintf(treeName,"TreeK%d",event);
633   fTreeK = (TTree*)gDirectory->Get(treeName);
634   if (fTreeK) fTreeK->SetBranchAddress("Particles", &fParticles);
635   else    Error("GetEvent","cannot find Kine Tree for event:%d\n",event);
636   
637   // Get Hits Tree header from file
638   sprintf(treeName,"TreeH%d",event);
639   fTreeH = (TTree*)gDirectory->Get(treeName);
640   if (!fTreeH) {
641     Error("GetEvent","cannot find Hits Tree for event:%d\n",event);
642   }
643   
644   // Get Digits Tree header from file
645   sprintf(treeName,"TreeD%d",event);
646   fTreeD = (TTree*)gDirectory->Get(treeName);
647   if (!fTreeD) {
648     Warning("GetEvent","cannot find Digits Tree for event:%d\n",event);
649   }
650   
651   
652   // Get Reconstruct Tree header from file
653   sprintf(treeName,"TreeR%d",event);
654   fTreeR = (TTree*)gDirectory->Get(treeName);
655   if (!fTreeR) {
656     //    printf("WARNING: cannot find Reconstructed Tree for event:%d\n",event);
657   }
658  
659   // Set Trees branch addresses
660   TIter next(fModules);
661   AliModule *detector;
662   while((detector = (AliModule*)next())) {
663     detector->SetTreeAddress();
664   }
665   
666   if (fTreeK) fTreeK->GetEvent(0);
667   fNtrack = Int_t (fParticles->GetEntries());
668   return fNtrack;
669 }
670
671 //_____________________________________________________________________________
672 TGeometry *AliRun::GetGeometry()
673 {
674   //
675   // Import Alice geometry from current file
676   // Return pointer to geometry object
677   //
678   if (!fGeometry) fGeometry = (TGeometry*)gDirectory->Get("AliceGeom");
679   //
680   // Unlink and relink nodes in detectors
681   // This is bad and there must be a better way...
682   //
683   
684   TIter next(fModules);
685   AliModule *detector;
686   while((detector = (AliModule*)next())) {
687     detector->SetTreeAddress();
688     TList *dnodes=detector->Nodes();
689     Int_t j;
690     TNode *node, *node1;
691     for ( j=0; j<dnodes->GetSize(); j++) {
692       node = (TNode*) dnodes->At(j);
693       node1 = fGeometry->GetNode(node->GetName());
694       dnodes->Remove(node);
695       dnodes->AddAt(node1,j);
696     }
697   }
698   return fGeometry;
699 }
700
701 //_____________________________________________________________________________
702 void AliRun::GetNextTrack(Int_t &mtrack, Int_t &ipart, Float_t *pmom,
703                           Float_t &e, Float_t *vpos, Float_t *polar,
704                           Float_t &tof)
705 {
706   //
707   // Return next track from stack of particles
708   //
709   TVector3 pol;
710   fCurrent=-1;
711   TParticle *track;
712   for(Int_t i=fNtrack-1; i>=0; i--) {
713     track=(TParticle*) fParticles->UncheckedAt(i);
714     if(!track->TestBit(Done_Bit)) {
715       //
716       // The track has not yet been processed
717       fCurrent=i;
718       ipart=track->GetPdgCode();
719       pmom[0]=track->Px();
720       pmom[1]=track->Py(); 
721       pmom[2]=track->Pz();
722       e     =track->Energy();
723       vpos[0]=track->Vx();
724       vpos[1]=track->Vy();
725       vpos[2]=track->Vz();
726       track->GetPolarisation(pol);
727       polar[0]=pol.X();
728       polar[1]=pol.Y();
729       polar[2]=pol.Z();
730       tof=track->T();
731       track->SetBit(Done_Bit);
732       break;
733     }
734   }
735   mtrack=fCurrent;
736   //
737   // stop and start timer when we start a primary track
738   Int_t nprimaries = fHeader.GetNprimary();
739   if (fCurrent >= nprimaries) return;
740   if (fCurrent < nprimaries-1) {
741     fTimer.Stop();
742     track=(TParticle*) fParticles->UncheckedAt(fCurrent+1);
743     //    track->SetProcessTime(fTimer.CpuTime());
744   }
745   fTimer.Start();
746 }
747
748 //_____________________________________________________________________________
749 Int_t AliRun::GetPrimary(Int_t track)
750 {
751   //
752   // return number of primary that has generated track
753   //
754   int current, parent;
755   TParticle *part;
756   //
757   parent=track;
758   while (1) {
759     current=parent;
760     part = (TParticle *)fParticles->UncheckedAt(current);
761     parent=part->GetFirstMother();
762     if(parent<0) return current;
763   }
764 }
765  
766 //_____________________________________________________________________________
767 void AliRun::Init(const char *setup)
768 {
769   //
770   // Initialize the Alice setup
771   //
772
773   gROOT->LoadMacro(setup);
774   gInterpreter->ProcessLine("Config();");
775
776   gMC->DefineParticles();  //Create standard MC particles
777
778   TObject *objfirst, *objlast;
779
780   fNdets = fModules->GetLast()+1;
781
782   //
783   //=================Create Materials, geometry, histograms, etc
784    TIter next(fModules);
785    AliModule *detector;
786    while((detector = (AliModule*)next())) {
787       detector->SetTreeAddress();
788       objlast = gDirectory->GetList()->Last();
789       
790       // Initialise detector materials, geometry, histograms,etc
791       detector->CreateMaterials();
792       detector->CreateGeometry();
793       detector->BuildGeometry();
794       detector->Init();
795       
796       // Add Detector histograms in Detector list of histograms
797       if (objlast) objfirst = gDirectory->GetList()->After(objlast);
798       else         objfirst = gDirectory->GetList()->First();
799       while (objfirst) {
800         detector->Histograms()->Add(objfirst);
801         objfirst = gDirectory->GetList()->After(objfirst);
802       }
803    }
804    SetTransPar(); //Read the cuts for all materials
805    
806    MediaTable(); //Build the special IMEDIA table
807    
808    //Close the geometry structure
809    gMC->Ggclos();
810    
811    //Initialise geometry deposition table
812    sEventEnergy.Set(gMC->NofVolumes()+1);
813    sSummEnergy.Set(gMC->NofVolumes()+1);
814    sSum2Energy.Set(gMC->NofVolumes()+1);
815    
816    //Create the color table
817    gMC->SetColors();
818    
819    //Compute cross-sections
820    gMC->Gphysi();
821    
822    //Write Geometry object to current file.
823    fGeometry->Write();
824    
825    fInitDone = kTRUE;
826 }
827
828 //_____________________________________________________________________________
829 void AliRun::MediaTable()
830 {
831   //
832   // Built media table to get from the media number to
833   // the detector id
834   //
835   Int_t kz, nz, idt, lz, i, k, ind;
836   //  Int_t ibeg;
837   TObjArray &dets = *gAlice->Detectors();
838   AliModule *det;
839   //
840   // For all detectors
841   for (kz=0;kz<fNdets;kz++) {
842     // If detector is defined
843     if((det=(AliModule*) dets[kz])) {
844         TArrayI &idtmed = *(det->GetIdtmed()); 
845         for(nz=0;nz<100;nz++) {
846         // Find max and min material number
847         if((idt=idtmed[nz])) {
848           det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
849           det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
850         }
851       }
852       if(det->LoMedium() > det->HiMedium()) {
853         det->LoMedium() = 0;
854         det->HiMedium() = 0;
855       } else {
856         if(det->HiMedium() > fImedia->GetSize()) {
857           Error("MediaTable","Increase fImedia from %d to %d",
858                 fImedia->GetSize(),det->HiMedium());
859           return;
860         }
861         // Tag all materials in rage as belonging to detector kz
862         for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
863           (*fImedia)[lz]=kz;
864         }
865       }
866     }
867   }
868   //
869   // Print summary table
870   printf(" Traking media ranges:\n");
871   for(i=0;i<(fNdets-1)/6+1;i++) {
872     for(k=0;k< (6<fNdets-i*6?6:fNdets-i*6);k++) {
873       ind=i*6+k;
874       det=(AliModule*)dets[ind];
875       if(det)
876         printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
877                det->HiMedium());
878       else
879         printf(" %6s: %3d -> %3d;","NULL",0,0);
880     }
881     printf("\n");
882   }
883 }
884
885 //____________________________________________________________________________
886 void AliRun::SetGenerator(AliGenerator *generator)
887 {
888   //
889   // Load the event generator
890   //
891   if(!fGenerator) fGenerator = generator;
892 }
893
894 //____________________________________________________________________________
895 void AliRun::SetTransPar(char* filename)
896 {
897   //
898   // Read filename to set the transport parameters
899   //
900
901
902   const Int_t ncuts=10;
903   const Int_t nflags=11;
904   const Int_t npars=ncuts+nflags;
905   const char pars[npars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
906                                "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
907                                "BREM","COMP","DCAY","DRAY","HADR","LOSS",
908                                "MULS","PAIR","PHOT","RAYL"};
909   char line[256];
910   char detName[7];
911   char* filtmp;
912   Float_t cut[ncuts];
913   Int_t flag[nflags];
914   Int_t i, itmed, iret, ktmed, kz;
915   FILE *lun;
916   //
917   // See whether the file is there
918   filtmp=gSystem->ExpandPathName(filename);
919   lun=fopen(filtmp,"r");
920   delete [] filtmp;
921   if(!lun) {
922     Warning("SetTransPar","File %s does not exist!\n",filename);
923     return;
924   }
925   //
926   printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
927   printf(" *%59s\n","*");
928   printf(" *       Please check carefully what you are doing!%10s\n","*");
929   printf(" *%59s\n","*");
930   //
931   while(1) {
932     // Initialise cuts and flags
933     for(i=0;i<ncuts;i++) cut[i]=-99;
934     for(i=0;i<nflags;i++) flag[i]=-99;
935     itmed=0;
936     for(i=0;i<256;i++) line[i]='\0';
937     // Read up to the end of line excluded
938     iret=fscanf(lun,"%[^\n]",line);
939     if(iret<0) {
940       //End of file
941       fclose(lun);
942       printf(" *%59s\n","*");
943       printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
944       return;
945     }
946     // Read the end of line
947     fscanf(lun,"%*c");
948     if(!iret) continue;
949     if(line[0]=='*') continue;
950     // Read the numbers
951     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",
952                 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
953                 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
954                 &flag[8],&flag[9],&flag[10]);
955     if(!iret) continue;
956     if(iret<0) {
957       //reading error
958       Warning("SetTransPar","Error reading file %s\n",filename);
959       continue;
960     }
961     // Check that the module exist
962     AliModule *mod = GetModule(detName);
963     if(mod) {
964       // Get the array of media numbers
965       TArrayI &idtmed = *mod->GetIdtmed();
966       // Check that the tracking medium code is valid
967       if(0<=itmed && itmed < 100) {
968         ktmed=idtmed[itmed];
969         if(!ktmed) {
970           Warning("SetTransPar","Invalid tracking medium code %d for %s\n",itmed,mod->GetName());
971           continue;
972         }
973         // Set energy thresholds
974         for(kz=0;kz<ncuts;kz++) {
975           if(cut[kz]>=0) {
976             printf(" *  %-6s set to %10.3E for tracking medium code %4d for %s\n",
977                    pars[kz],cut[kz],itmed,mod->GetName());
978             gMC->Gstpar(ktmed,pars[kz],cut[kz]);
979           }
980         }
981         // Set transport mechanisms
982         for(kz=0;kz<nflags;kz++) {
983           if(flag[kz]>=0) {
984             printf(" *  %-6s set to %10d for tracking medium code %4d for %s\n",
985                    pars[ncuts+kz],flag[kz],itmed,mod->GetName());
986             gMC->Gstpar(ktmed,pars[ncuts+kz],Float_t(flag[kz]));
987           }
988         }
989       } else {
990         Warning("SetTransPar","Invalid medium code %d *\n",itmed);
991         continue;
992       }
993     } else {
994       Warning("SetTransPar","Module %s not present\n",detName);
995       continue;
996     }
997   }
998 }
999
1000 //_____________________________________________________________________________
1001 void AliRun::MakeTree(Option_t *option)
1002 {
1003   //
1004   //  Create the ROOT trees
1005   //  Loop on all detectors to create the Root branch (if any)
1006   //
1007
1008   //
1009   // Analyse options
1010   char *K = strstr(option,"K");
1011   char *H = strstr(option,"H");
1012   char *E = strstr(option,"E");
1013   char *D = strstr(option,"D");
1014   char *R = strstr(option,"R");
1015   //
1016   if (K && !fTreeK) fTreeK = new TTree("TreeK0","Kinematics");
1017   if (H && !fTreeH) fTreeH = new TTree("TreeH0","Hits");
1018   if (D && !fTreeD) fTreeD = new TTree("TreeD0","Digits");
1019   if (E && !fTreeE) fTreeE = new TTree("TE","Header");
1020   if (R && !fTreeR) fTreeR = new TTree("TreeR0","Reconstruction");
1021   if (fTreeH) fTreeH->SetAutoSave(1000000000); //no autosave
1022   //
1023   // Create a branch for hits/digits for each detector
1024   // Each branch is a TClonesArray. Each data member of the Hits classes
1025   // will be in turn a subbranch of the detector master branch
1026   TIter next(fModules);
1027   AliModule *detector;
1028   while((detector = (AliModule*)next())) {
1029      if (H || D || R) detector->MakeBranch(option);
1030   }
1031   //  Create a branch for particles
1032   if (fTreeK && K) fTreeK->Branch("Particles",&fParticles,4000);
1033   
1034   //  Create a branch for Header
1035   if (fTreeE && E) fTreeE->Branch("Header","AliHeader",&header,4000);
1036 }
1037
1038 //_____________________________________________________________________________
1039 Int_t AliRun::PurifyKine(Int_t lastSavedTrack, Int_t nofTracks)
1040 {
1041   //
1042   // PurifyKine with external parameters
1043   //
1044   fHgwmk = lastSavedTrack;
1045   fNtrack = nofTracks;
1046   PurifyKine();
1047   return fHgwmk;
1048 }
1049
1050 //_____________________________________________________________________________
1051 void AliRun::PurifyKine()
1052 {
1053   //
1054   // Compress kinematic tree keeping only flagged particles
1055   // and renaming the particle id's in all the hits
1056   //
1057   TClonesArray &particles = *fParticles;
1058   int nkeep=fHgwmk+1, parent, i;
1059   TParticle *part, *partnew, *father;
1060   AliHit *OneHit;
1061   int *map = new int[particles.GetEntries()];
1062
1063   // Save in Header total number of tracks before compression
1064   fHeader.SetNtrack(fHeader.GetNtrack()+fNtrack-fHgwmk);
1065
1066   // Preset map, to be removed later
1067   for(i=0; i<fNtrack; i++) {
1068     if(i<=fHgwmk) map[i]=i ; else map[i] = -99 ;}
1069   // Second pass, build map between old and new numbering
1070   for(i=fHgwmk+1; i<fNtrack; i++) {
1071     part = (TParticle *)particles.UncheckedAt(i);
1072     if(part->TestBit(Keep_Bit)) {
1073       
1074       // This particle has to be kept
1075       map[i]=nkeep;
1076       if(i!=nkeep) {
1077         
1078         // Old and new are different, have to copy
1079         partnew = (TParticle *)particles.UncheckedAt(nkeep);
1080         *partnew = *part;
1081       } else partnew = part;
1082       
1083       // as the parent is always *before*, it must be already
1084       // in place. This is what we are checking anyway!
1085       if((parent=partnew->GetFirstMother())>fHgwmk) {
1086         if(map[parent]==-99) printf("map[%d] = -99!\n",parent);
1087         partnew->SetFirstMother(map[parent]);
1088       }
1089       nkeep++;
1090     }
1091   }
1092   fNtrack=nkeep;
1093   
1094   // Fix daughters information
1095   for (i=fHgwmk+1; i<fNtrack; i++) {
1096     part = (TParticle *)particles.UncheckedAt(i);
1097     parent = part->GetFirstMother();
1098     father = (TParticle *)particles.UncheckedAt(parent);
1099     if(father->TestBit(Daughters_Bit)) {
1100       
1101       if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
1102       if(i>father->GetLastDaughter())  father->SetLastDaughter(i);
1103     } else {
1104       // Iitialise daughters info for first pass
1105       father->SetFirstDaughter(i);
1106       father->SetLastDaughter(i);
1107       father->SetBit(Daughters_Bit);
1108     }
1109   }
1110   
1111   // Now loop on all detectors and reset the hits
1112   TIter next(fModules);
1113   AliModule *detector;
1114   while((detector = (AliModule*)next())) {
1115     if (!detector->Hits()) continue;
1116     TClonesArray &vHits=*(detector->Hits());
1117     if(vHits.GetEntries() != detector->GetNhits())
1118       printf("vHits.GetEntries()!=detector->GetNhits(): %d != %d\n",
1119              vHits.GetEntries(),detector->GetNhits());
1120     for (i=0; i<detector->GetNhits(); i++) {
1121       OneHit = (AliHit *)vHits.UncheckedAt(i);
1122       OneHit->SetTrack(map[OneHit->GetTrack()]);
1123     }
1124   }
1125
1126   fHgwmk=nkeep-1;
1127   particles.SetLast(fHgwmk);
1128   delete [] map;
1129 }
1130
1131 //_____________________________________________________________________________
1132 void AliRun::Reset(Int_t run, Int_t idevent)
1133 {
1134   //
1135   //  Reset all Detectors & kinematics & trees
1136   //
1137   char hname[30];
1138   //
1139   ResetStack();
1140   ResetHits();
1141   ResetDigits();
1142
1143   // Initialise event header
1144   fHeader.Reset(run,idevent);
1145
1146   if(fTreeK) {
1147     fTreeK->Reset();
1148     sprintf(hname,"TreeK%d",idevent);
1149     fTreeK->SetName(hname);
1150   }
1151   if(fTreeH) {
1152     fTreeH->Reset();
1153     sprintf(hname,"TreeH%d",idevent);
1154     fTreeH->SetName(hname);
1155   }
1156   if(fTreeD) {
1157     fTreeD->Reset();
1158     sprintf(hname,"TreeD%d",idevent);
1159     fTreeD->SetName(hname);
1160   }
1161   if(fTreeR) {
1162     fTreeR->Reset();
1163     sprintf(hname,"TreeR%d",idevent);
1164     fTreeR->SetName(hname);
1165   }
1166 }
1167
1168 //_____________________________________________________________________________
1169 void AliRun::ResetDigits()
1170 {
1171   //
1172   //  Reset all Detectors digits
1173   //
1174   TIter next(fModules);
1175   AliModule *detector;
1176   while((detector = (AliModule*)next())) {
1177      detector->ResetDigits();
1178   }
1179 }
1180
1181 //_____________________________________________________________________________
1182 void AliRun::ResetHits()
1183 {
1184   //
1185   //  Reset all Detectors hits
1186   //
1187   TIter next(fModules);
1188   AliModule *detector;
1189   while((detector = (AliModule*)next())) {
1190      detector->ResetHits();
1191   }
1192 }
1193
1194 //_____________________________________________________________________________
1195 void AliRun::ResetPoints()
1196 {
1197   //
1198   // Reset all Detectors points
1199   //
1200   TIter next(fModules);
1201   AliModule *detector;
1202   while((detector = (AliModule*)next())) {
1203      detector->ResetPoints();
1204   }
1205 }
1206
1207 //_____________________________________________________________________________
1208 void AliRun::Run(Int_t nevent, const char *setup)
1209 {
1210   //
1211   // Main function to be called to process a galice run
1212   // example
1213   //    Root > gAlice.Run(); 
1214   // a positive number of events will cause the finish routine
1215   // to be called
1216   //
1217
1218   Int_t i, todo;
1219   // check if initialisation has been done
1220   if (!fInitDone) Init(setup);
1221   
1222   // Create the Root Tree with one branch per detector
1223   if(!fEvent) {
1224     gAlice->MakeTree("KHDER");
1225   }
1226
1227   todo = TMath::Abs(nevent);
1228   for (i=0; i<todo; i++) {
1229   // Process one run (one run = one event)
1230      gAlice->Reset(fRun, fEvent);
1231      gMC->Gtrigi();
1232      gMC->Gtrigc();
1233      gMC->Gtrig();
1234      gAlice->FinishEvent();
1235      fEvent++;
1236   }
1237   
1238   // End of this run, close files
1239   if(nevent>0) gAlice->FinishRun();
1240 }
1241
1242 //_____________________________________________________________________________
1243 void AliRun::RunLego(const char *setup,Int_t ntheta,Float_t themin,
1244                      Float_t themax,Int_t nphi,Float_t phimin,Float_t phimax,
1245                      Float_t rmin,Float_t rmax,Float_t zmax)
1246 {
1247   //
1248   // Generates lego plots of:
1249   //    - radiation length map phi vs theta
1250   //    - radiation length map phi vs eta
1251   //    - interaction length map
1252   //    - g/cm2 length map
1253   //
1254   //  ntheta    bins in theta, eta
1255   //  themin    minimum angle in theta (degrees)
1256   //  themax    maximum angle in theta (degrees)
1257   //  nphi      bins in phi
1258   //  phimin    minimum angle in phi (degrees)
1259   //  phimax    maximum angle in phi (degrees)
1260   //  rmin      minimum radius
1261   //  rmax      maximum radius
1262   //  
1263   //
1264   //  The number of events generated = ntheta*nphi
1265   //  run input parameters in macro setup (default="Config.C")
1266   //
1267   //  Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
1268   //Begin_Html
1269   /*
1270     <img src="picts/AliRunLego1.gif">
1271   */
1272   //End_Html
1273   //Begin_Html
1274   /*
1275     <img src="picts/AliRunLego2.gif">
1276   */
1277   //End_Html
1278   //Begin_Html
1279   /*
1280     <img src="picts/AliRunLego3.gif">
1281   */
1282   //End_Html
1283   //
1284
1285   // check if initialisation has been done
1286   if (!fInitDone) Init(setup);
1287   
1288   fLego = new AliLego("lego","lego");
1289   fLego->Init(ntheta,themin,themax,nphi,phimin,phimax,rmin,rmax,zmax);
1290   fLego->Run();
1291   
1292   // Create only the Root event Tree
1293   gAlice->MakeTree("E");
1294   
1295   // End of this run, close files
1296   gAlice->FinishRun();
1297 }
1298
1299 //_____________________________________________________________________________
1300 void AliRun::SetCurrentTrack(Int_t track)
1301
1302   //
1303   // Set current track number
1304   //
1305   fCurrent = track; 
1306 }
1307  
1308 //_____________________________________________________________________________
1309 void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
1310                       Float_t *vpos, Float_t *polar, Float_t tof,
1311                       const char *mecha, Int_t &ntr, Float_t weight)
1312
1313   //
1314   // Load a track on the stack
1315   //
1316   // done     0 if the track has to be transported
1317   //          1 if not
1318   // parent   identifier of the parent track. -1 for a primary
1319   // pdg    particle code
1320   // pmom     momentum GeV/c
1321   // vpos     position 
1322   // polar    polarisation 
1323   // tof      time of flight in seconds
1324   // mecha    production mechanism
1325   // ntr      on output the number of the track stored
1326   //
1327   TClonesArray &particles = *fParticles;
1328   TParticle *particle;
1329   Float_t mass;
1330   const Int_t firstdaughter=-1;
1331   const Int_t lastdaughter=-1;
1332   const Int_t KS=0;
1333   //  const Float_t tlife=0;
1334   
1335   //
1336   // Here we get the static mass
1337   // For MC is ok, but a more sophisticated method could be necessary
1338   // if the calculated mass is required
1339   // also, this method is potentially dangerous if the mass
1340   // used in the MC is not the same of the PDG database
1341   //
1342   mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
1343   Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
1344                         pmom[1]*pmom[1]+pmom[2]*pmom[2]);
1345   
1346   //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",
1347   //pname,mass,e,fNtrack,pdg,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],KS,mecha);
1348   
1349   particle=new(particles[fNtrack]) TParticle(pdg,KS,parent,-1,firstdaughter,
1350                                              lastdaughter,pmom[0],pmom[1],pmom[2],
1351                                              e,vpos[0],vpos[1],vpos[2],tof);
1352   //                                         polar[0],polar[1],polar[2],tof,
1353   //                                         mecha,weight);
1354   ((TParticle*)particles[fNtrack])->SetPolarisation(TVector3(polar[0],polar[1],polar[2]));
1355   ((TParticle*)particles[fNtrack])->SetWeight(weight);
1356   if(!done) particle->SetBit(Done_Bit);
1357   
1358   if(parent>=0) {
1359     particle=(TParticle*) fParticles->UncheckedAt(parent);
1360     particle->SetLastDaughter(fNtrack);
1361     if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
1362   } else { 
1363     //
1364     // This is a primary track. Set high water mark for this event
1365     fHgwmk=fNtrack;
1366     //
1367     // Set also number if primary tracks
1368     fHeader.SetNprimary(fHgwmk+1);
1369     fHeader.SetNtrack(fHgwmk+1);
1370   }
1371   ntr = fNtrack++;
1372 }
1373
1374 //_____________________________________________________________________________
1375 void AliRun::KeepTrack(const Int_t track)
1376
1377   //
1378   // flags a track to be kept
1379   //
1380   TClonesArray &particles = *fParticles;
1381   ((TParticle*)particles[track])->SetBit(Keep_Bit);
1382 }
1383  
1384 //_____________________________________________________________________________
1385 void AliRun::StepManager(Int_t id) const
1386 {
1387   //
1388   // Called at every step during transport
1389   //
1390
1391   Int_t copy;
1392   //
1393   // --- If lego option, do it and leave 
1394   if (fLego) {
1395     fLego->StepManager();
1396     return;
1397   }
1398   //Update energy deposition tables
1399   sEventEnergy[gMC->CurrentVolID(copy)]+=gMC->Edep();
1400   
1401   //Call the appropriate stepping routine;
1402   AliModule *det = (AliModule*)fModules->At(id);
1403   if(det) det->StepManager();
1404 }
1405
1406 //_____________________________________________________________________________
1407 void AliRun::ReadEuclid(const char* filnam, const AliModule *det, char* topvol)
1408 {
1409   //                                                                     
1410   //       read in the geometry of the detector in euclid file format    
1411   //                                                                     
1412   //        id_det : the detector identification (2=its,...)            
1413   //        topvol : return parameter describing the name of the top    
1414   //        volume of geometry.                                          
1415   //                                                                     
1416   //            author : m. maire                                        
1417   //                                                                     
1418   //     28.07.98
1419   //     several changes have been made by miroslav helbich
1420   //     subroutine is rewrited to follow the new established way of memory
1421   //     booking for tracking medias and rotation matrices.
1422   //     all used tracking media have to be defined first, for this you can use
1423   //     subroutine  greutmed.
1424   //     top volume is searched as only volume not positioned into another 
1425   //
1426
1427   Int_t i, nvol, iret, itmed, irot, numed, npar, ndiv, iaxe;
1428   Int_t ndvmx, nr, flag;
1429   char key[5], card[77], natmed[21];
1430   char name[5], mother[5], shape[5], konly[5], volst[7000][5];
1431   char *filtmp;
1432   Float_t par[50];
1433   Float_t teta1, phi1, teta2, phi2, teta3, phi3, orig, step;
1434   Float_t xo, yo, zo;
1435   Int_t idrot[5000],istop[7000];
1436   FILE *lun;
1437   //
1438   // *** The input filnam name will be with extension '.euc'
1439   filtmp=gSystem->ExpandPathName(filnam);
1440   lun=fopen(filtmp,"r");
1441   delete [] filtmp;
1442   if(!lun) {
1443     Error("ReadEuclid","Could not open file %s\n",filnam);
1444     return;
1445   }
1446   //* --- definition of rotation matrix 0 ---  
1447   TArrayI &idtmed = *(det->GetIdtmed());
1448   idrot[0]=0;
1449   nvol=0;
1450  L10:
1451   for(i=0;i<77;i++) card[i]=0;
1452   iret=fscanf(lun,"%77[^\n]",card);
1453   if(iret<=0) goto L20;
1454   fscanf(lun,"%*c");
1455   //*
1456   strncpy(key,card,4);
1457   key[4]='\0';
1458   if (!strcmp(key,"TMED")) {
1459     sscanf(&card[5],"%d '%[^']'",&itmed,natmed);
1460     //Pad the string with blanks
1461     i=-1;
1462     while(natmed[++i]);
1463     while(i<20) natmed[i++]=' ';
1464     natmed[i]='\0';
1465     //
1466     gMC->Gckmat(idtmed[itmed],natmed);
1467     //*
1468   } else if (!strcmp(key,"ROTM")) {
1469     sscanf(&card[4],"%d %f %f %f %f %f %f",&irot,&teta1,&phi1,&teta2,&phi2,&teta3,&phi3);
1470     det->AliMatrix(idrot[irot],teta1,phi1,teta2,phi2,teta3,phi3);
1471     //*
1472   } else if (!strcmp(key,"VOLU")) {
1473     sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, shape, &numed, &npar);
1474     if (npar>0) {
1475       for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
1476       fscanf(lun,"%*c");
1477     }
1478     gMC->Gsvolu( name, shape, idtmed[numed], par, npar);
1479     //*     save the defined volumes
1480     strcpy(volst[++nvol],name);
1481     istop[nvol]=1;
1482     //*
1483   } else if (!strcmp(key,"DIVN")) {
1484     sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, mother, &ndiv, &iaxe);
1485     gMC->Gsdvn  ( name, mother, ndiv, iaxe );
1486     //*
1487   } else if (!strcmp(key,"DVN2")) {
1488     sscanf(&card[5],"'%[^']' '%[^']' %d %d %f %d",name, mother, &ndiv, &iaxe, &orig, &numed);
1489     gMC->Gsdvn2( name, mother, ndiv, iaxe, orig,idtmed[numed]);
1490     //*
1491   } else if (!strcmp(key,"DIVT")) {
1492     sscanf(&card[5],"'%[^']' '%[^']' %f %d %d %d", name, mother, &step, &iaxe, &numed, &ndvmx);
1493     gMC->Gsdvt ( name, mother, step, iaxe, idtmed[numed], ndvmx);
1494     //*
1495   } else if (!strcmp(key,"DVT2")) {
1496     sscanf(&card[5],"'%[^']' '%[^']' %f %d %f %d %d", name, mother, &step, &iaxe, &orig, &numed, &ndvmx);
1497     gMC->Gsdvt2 ( name, mother, step, iaxe, orig, idtmed[numed], ndvmx );
1498     //*
1499   } else if (!strcmp(key,"POSI")) {
1500     sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']'", name, &nr, mother, &xo, &yo, &zo, &irot, konly);
1501     //*** volume name cannot be the top volume
1502     for(i=1;i<=nvol;i++) {
1503       if (!strcmp(volst[i],name)) istop[i]=0;
1504     }
1505     //*
1506     gMC->Gspos  ( name, nr, mother, xo, yo, zo, idrot[irot], konly );
1507     //*
1508   } else if (!strcmp(key,"POSP")) {
1509     sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']' %d", name, &nr, mother, &xo, &yo, &zo, &irot, konly, &npar);
1510     if (npar > 0) {
1511       for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
1512       fscanf(lun,"%*c");
1513     }
1514     //*** volume name cannot be the top volume
1515     for(i=1;i<=nvol;i++) {
1516       if (!strcmp(volst[i],name)) istop[i]=0;
1517     }
1518     //*
1519     gMC->Gsposp ( name, nr, mother, xo,yo,zo, idrot[irot], konly, par, npar);
1520   }
1521   //*
1522   if (strcmp(key,"END")) goto L10;
1523   //* find top volume in the geometry
1524   flag=0;
1525   for(i=1;i<=nvol;i++) {
1526     if (istop[i] && flag) {
1527       Warning("ReadEuclid"," %s is another possible top volume\n",volst[i]);
1528     }
1529     if (istop[i] && !flag) {
1530       strcpy(topvol,volst[i]);
1531       printf(" *** GREUCL *** volume %s taken as a top volume\n",topvol);
1532       flag=1;
1533     }
1534   }
1535   if (!flag) {
1536     Warning("ReadEuclid","top volume not found\n");
1537   }
1538   fclose (lun);
1539   //*
1540   //*     commented out only for the not cernlib version
1541   printf(" *** GREUCL *** file: %s is now read in\n",filnam);
1542   //
1543   return;
1544   //*
1545   L20:
1546   Error("ReadEuclid","reading error or premature end of file\n");
1547 }
1548
1549 //_____________________________________________________________________________
1550 void AliRun::ReadEuclidMedia(const char* filnam, const AliModule *det)
1551 {
1552   //                                                                     
1553   //       read in the materials and tracking media for the detector     
1554   //                   in euclid file format                             
1555   //                                                                     
1556   //       filnam: name of the input file                                
1557   //       id_det: id_det is the detector identification (2=its,...)     
1558   //                                                                     
1559   //            author : miroslav helbich                                
1560   //
1561   Float_t sxmgmx = gAlice->Field()->Max();
1562   Int_t   isxfld = gAlice->Field()->Integ();
1563   Int_t end, i, iret, itmed;
1564   char key[5], card[130], natmed[21], namate[21];
1565   Float_t ubuf[50];
1566   char* filtmp;
1567   FILE *lun;
1568   Int_t imate;
1569   Int_t nwbuf, isvol, ifield, nmat;
1570   Float_t a, z, dens, radl, absl, fieldm, tmaxfd, stemax, deemax, epsil, stmin;
1571   //
1572   end=strlen(filnam);
1573   for(i=0;i<end;i++) if(filnam[i]=='.') {
1574     end=i;
1575     break;
1576   }
1577   //
1578   // *** The input filnam name will be with extension '.euc'
1579   printf("The file name is %s\n",filnam); //Debug
1580   filtmp=gSystem->ExpandPathName(filnam);
1581   lun=fopen(filtmp,"r");
1582   delete [] filtmp;
1583   if(!lun) {
1584     Warning("ReadEuclidMedia","Could not open file %s\n",filnam);
1585     return;
1586   }
1587   //
1588   // Retrieve Mag Field parameters
1589   Int_t ISXFLD=gAlice->Field()->Integ();
1590   Float_t SXMGMX=gAlice->Field()->Max();
1591   //  TArrayI &idtmed = *(det->GetIdtmed());
1592   //
1593  L10:
1594   for(i=0;i<130;i++) card[i]=0;
1595   iret=fscanf(lun,"%4s %[^\n]",key,card);
1596   if(iret<=0) goto L20;
1597   fscanf(lun,"%*c");
1598   //*
1599   //* read material
1600   if (!strcmp(key,"MATE")) {
1601     sscanf(card,"%d '%[^']' %f %f %f %f %f %d",&imate,namate,&a,&z,&dens,&radl,&absl,&nwbuf);
1602     if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
1603     //Pad the string with blanks
1604     i=-1;
1605     while(namate[++i]);
1606     while(i<20) namate[i++]=' ';
1607     namate[i]='\0';
1608     //
1609     det->AliMaterial(imate,namate,a,z,dens,radl,absl,ubuf,nwbuf);
1610     //* read tracking medium
1611   } else if (!strcmp(key,"TMED")) {
1612     sscanf(card,"%d '%[^']' %d %d %d %f %f %f %f %f %f %d",
1613            &itmed,natmed,&nmat,&isvol,&ifield,&fieldm,&tmaxfd,
1614            &stemax,&deemax,&epsil,&stmin,&nwbuf);
1615     if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
1616     if (ifield<0) ifield=isxfld;
1617     if (fieldm<0) fieldm=sxmgmx;
1618     //Pad the string with blanks
1619     i=-1;
1620     while(natmed[++i]);
1621     while(i<20) natmed[i++]=' ';
1622     natmed[i]='\0';
1623     //
1624     det->AliMedium(itmed,natmed,nmat,isvol,ISXFLD,SXMGMX,tmaxfd,
1625                    stemax,deemax,epsil,stmin,ubuf,nwbuf);
1626     //    (*fImedia)[idtmed[itmed]-1]=id_det;
1627     //*
1628   }
1629   //*
1630   if (strcmp(key,"END")) goto L10;
1631   fclose (lun);
1632   //*
1633   //*     commented out only for the not cernlib version
1634   Warning("ReadEuclidMedia","file: %s is now read in\n",filnam);
1635   //*
1636   return;
1637   //*
1638  L20:
1639   Warning("ReadEuclidMedia","reading error or premature end of file\n");
1640
1641  
1642 //_____________________________________________________________________________
1643 void AliRun::Streamer(TBuffer &R__b)
1644 {
1645   //
1646   // Stream an object of class AliRun.
1647   //
1648   if (R__b.IsReading()) {
1649     Version_t R__v = R__b.ReadVersion(); if (R__v) { }
1650     TNamed::Streamer(R__b);
1651     if (!gAlice) gAlice = this;
1652     gROOT->GetListOfBrowsables()->Add(this,"Run");
1653     fTreeE = (TTree*)gDirectory->Get("TE");
1654     if (fTreeE) fTreeE->SetBranchAddress("Header", &header);
1655     else    Error("Streamer","cannot find Header Tree\n");
1656     R__b >> fNtrack;
1657     R__b >> fHgwmk;
1658     R__b >> fDebug;
1659     fHeader.Streamer(R__b);
1660     R__b >> fModules;
1661     R__b >> fParticles;
1662     R__b >> fField; 
1663     //    R__b >> fMC;
1664     R__b >> fNdets;
1665     R__b >> fTrRmax;
1666     R__b >> fTrZmax;
1667     R__b >> fGenerator;
1668     if(R__v>1) {
1669       R__b >> fPDGDB;        //Particle factory object!
1670       fTreeE->GetEntry(0);
1671     } else {
1672       fHeader.SetEvent(0);
1673       fPDGDB     = TDatabasePDG::Instance();        //Particle factory object!
1674     }
1675   } else {
1676     R__b.WriteVersion(AliRun::IsA());
1677     TNamed::Streamer(R__b);
1678     R__b << fNtrack;
1679     R__b << fHgwmk;
1680     R__b << fDebug;
1681     fHeader.Streamer(R__b);
1682     R__b << fModules;
1683     R__b << fParticles;
1684     R__b << fField;
1685     //    R__b << fMC;
1686     R__b << fNdets;
1687     R__b << fTrRmax;
1688     R__b << fTrZmax;
1689     R__b << fGenerator;
1690     R__b << fPDGDB;        //Particle factory object!
1691   }
1692
1693
1694
1695 //_____________________________________________________________________________
1696 //
1697 //                 Interfaces to Fortran
1698 //
1699 //_____________________________________________________________________________
1700
1701 extern "C" void type_of_call  rxgtrak (Int_t &mtrack, Int_t &ipart, Float_t *pmom, 
1702                                        Float_t &e, Float_t *vpos, Float_t *polar,
1703                                        Float_t &tof)
1704 {
1705   //
1706   //     Fetches next track from the ROOT stack for transport. Called by the
1707   //     modified version of GTREVE.
1708   //
1709   //              Track number in the ROOT stack. If MTRACK=0 no
1710   //      mtrack  more tracks are left in the stack to be
1711   //              transported.
1712   //      ipart   Particle code in the GEANT conventions.
1713   //      pmom[3] Particle momentum in GeV/c
1714   //      e       Particle energy in GeV
1715   //      vpos[3] Particle position
1716   //      tof     Particle time of flight in seconds
1717   //
1718   Int_t pdg;
1719   gAlice->GetNextTrack(mtrack, pdg, pmom, e, vpos, polar, tof);
1720   ipart = gMC->IdFromPDG(pdg);
1721   mtrack++;
1722 }
1723
1724 //_____________________________________________________________________________
1725 extern "C" void type_of_call 
1726 #ifndef WIN32
1727 rxstrak (Int_t &keep, Int_t &parent, Int_t &ipart, Float_t *pmom, 
1728                Float_t *vpos, Float_t &tof, const char* cmech, Int_t &ntr, const int cmlen)
1729 #else
1730 rxstrak (Int_t &keep, Int_t &parent, Int_t &ipart, Float_t *pmom,
1731          Float_t *vpos, Float_t &tof, const char* cmech, const int cmlen,
1732          Int_t &ntr)
1733 #endif
1734 {
1735   //
1736   //     Fetches next track from the ROOT stack for transport. Called by GUKINE
1737   //     and GUSTEP.
1738   //
1739   //              Status of the track. If keep=0 the track is put
1740   //      keep    on the ROOT stack but it is not fetched for
1741   //              transport.
1742   //      parent  Parent track. If parent=0 the track is a primary.
1743   //              In GUSTEP the routine is normally called to store
1744   //              secondaries generated by the current track whose
1745   //              ROOT stack number is MTRACK (common SCKINE.
1746   //      ipart   Particle code in the GEANT conventions.
1747   //      pmom[3] Particle momentum in GeV/c
1748   //      vpos[3] Particle position
1749   //      tof     Particle time of flight in seconds
1750   //
1751   //      cmech   (CHARACTER*10) Particle origin. This field is user
1752   //              defined and it is not used inside the GALICE code.
1753   //      ntr     Number assigned to the particle in the ROOT stack.
1754   //
1755   char mecha[11];
1756   Float_t polar[3]={0.,0.,0.};
1757   for(int i=0; i<10 && i<cmlen; i++) mecha[i]=cmech[i];
1758   mecha[10]=0;
1759   Int_t pdg=gMC->PDGFromId(ipart);
1760   gAlice->SetTrack(keep, parent-1, pdg, pmom, vpos, polar, tof, mecha, ntr);
1761   ntr++;
1762 }
1763
1764 //_____________________________________________________________________________
1765 extern "C" void type_of_call  rxkeep(const Int_t &n)
1766 {
1767   if( NULL==gAlice ) exit(1);
1768   
1769   if( n<=0 || n>gAlice->Particles()->GetEntries() )
1770     {
1771       printf("  Bad index n=%d must be 0<n<=%d\n",
1772              n,gAlice->Particles()->GetEntries());
1773       exit(1);
1774     }
1775   
1776   ((TParticle*)(gAlice->Particles()->UncheckedAt(n-1)))->SetBit(Keep_Bit);
1777 }
1778
1779 //_____________________________________________________________________________
1780 extern "C" void type_of_call  rxouth ()
1781 {
1782   //
1783   // Called by Gtreve at the end of each primary track
1784   //
1785   gAlice->FinishPrimary();
1786 }
1787