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