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