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