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