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