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