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