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