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