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