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