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