Adding GEANT3 common blocks to known structures for CINT
[u/mrichter/AliRoot.git] / STEER / AliRun.cxx
CommitLineData
99d554c8 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16/*
17$Log$
b13db077 18Revision 1.24 2000/01/19 17:17:20 fca
19Introducing a list of lists of hits -- more hits allowed for detector now
20
1cedd08a 21Revision 1.23 1999/12/03 11:14:31 fca
22Fixing previous wrong checking
23
00719c1b 24Revision 1.21 1999/11/25 10:40:08 fca
25Fixing daughters information also in primary tracks
26
ae23d366 27Revision 1.20 1999/10/04 18:08:49 fca
28Adding protection against inconsistent Euclid files
29
3fcc96a1 30Revision 1.19 1999/09/29 07:50:40 fca
31Introduction of the Copyright and cvs Log
32
99d554c8 33*/
34
fe4da5cc 35///////////////////////////////////////////////////////////////////////////////
36// //
37// Control class for Alice C++ //
38// Only one single instance of this class exists. //
39// The object is created in main program aliroot //
40// and is pointed by the global gAlice. //
41// //
8494b010 42// -Supports the list of all Alice Detectors (fModules). //
fe4da5cc 43// -Supports the list of particles (fParticles). //
44// -Supports the Trees. //
45// -Supports the geometry. //
46// -Supports the event display. //
47//Begin_Html
48/*
1439f98e 49<img src="picts/AliRunClass.gif">
fe4da5cc 50*/
51//End_Html
52//Begin_Html
53/*
1439f98e 54<img src="picts/alirun.gif">
fe4da5cc 55*/
56//End_Html
57// //
58///////////////////////////////////////////////////////////////////////////////
59
60#include <TFile.h>
61#include <TRandom.h>
62#include <TBRIK.h>
63#include <TNode.h>
fe4da5cc 64#include <TCint.h>
65#include <TSystem.h>
a8f1fb7c 66#include <TObjectTable.h>
fe4da5cc 67
1578254f 68#include "TParticle.h"
fe4da5cc 69#include "AliRun.h"
fe4da5cc 70#include "AliDisplay.h"
b13db077 71#include "AliVMC.h"
fe4da5cc 72
fe4da5cc 73#include <stdlib.h>
74#include <stdio.h>
75#include <string.h>
76
77AliRun *gAlice;
78
79static AliHeader *header;
80
fe4da5cc 81static TArrayF sEventEnergy;
82static TArrayF sSummEnergy;
83static TArrayF sSum2Energy;
84
fe4da5cc 85ClassImp(AliRun)
86
87//_____________________________________________________________________________
88AliRun::AliRun()
89{
90 //
91 // Default constructor for AliRun
92 //
93 header=&fHeader;
94 fRun = 0;
95 fEvent = 0;
96 fCurrent = -1;
8494b010 97 fModules = 0;
fe4da5cc 98 fGenerator = 0;
99 fTreeD = 0;
100 fTreeK = 0;
101 fTreeH = 0;
102 fTreeE = 0;
103 fTreeR = 0;
104 fParticles = 0;
105 fGeometry = 0;
106 fDisplay = 0;
107 fField = 0;
108 fMC = 0;
109 fNdets = 0;
110 fImedia = 0;
111 fTrRmax = 1.e10;
112 fTrZmax = 1.e10;
fe4da5cc 113 fInitDone = kFALSE;
114 fLego = 0;
1578254f 115 fPDGDB = 0; //Particle factory object!
1cedd08a 116 fHitLists = 0;
fe4da5cc 117}
118
119//_____________________________________________________________________________
120AliRun::AliRun(const char *name, const char *title)
121 : TNamed(name,title)
122{
123 //
124 // Constructor for the main processor.
125 // Creates the geometry
126 // Creates the list of Detectors.
127 // Creates the list of particles.
128 //
129 Int_t i;
130
131 gAlice = this;
132 fTreeD = 0;
133 fTreeK = 0;
134 fTreeH = 0;
135 fTreeE = 0;
136 fTreeR = 0;
137 fTrRmax = 1.e10;
138 fTrZmax = 1.e10;
1141f8e4 139 fGenerator = 0;
fe4da5cc 140 fInitDone = kFALSE;
141 fLego = 0;
142 fField = 0;
143
144 gROOT->GetListOfBrowsables()->Add(this,name);
145 //
146 // create the support list for the various Detectors
8494b010 147 fModules = new TObjArray(77);
fe4da5cc 148 //
149 // Create the TNode geometry for the event display
150
151 BuildSimpleGeometry();
152
153
154 fNtrack=0;
155 fHgwmk=0;
156 fCurrent=-1;
157 header=&fHeader;
158 fRun = 0;
159 fEvent = 0;
160 //
161 // Create the particle stack
1578254f 162 fParticles = new TClonesArray("TParticle",100);
fe4da5cc 163
164 fDisplay = 0;
165 //
166 // Create default mag field
167 SetField();
168 //
cfce8870 169 fMC = gMC;
fe4da5cc 170 //
171 // Prepare the tracking medium lists
172 fImedia = new TArrayI(1000);
173 for(i=0;i<1000;i++) (*fImedia)[i]=-99;
1578254f 174 //
175 // Make particles
176 fPDGDB = TDatabasePDG::Instance(); //Particle factory object!
1cedd08a 177 //
178 // Create HitLists list
179 fHitLists = new TList();
fe4da5cc 180}
181
182//_____________________________________________________________________________
183AliRun::~AliRun()
184{
185 //
186 // Defaullt AliRun destructor
187 //
fe4da5cc 188 delete fImedia;
189 delete fField;
190 delete fMC;
191 delete fGeometry;
192 delete fDisplay;
193 delete fGenerator;
194 delete fLego;
195 delete fTreeD;
196 delete fTreeK;
197 delete fTreeH;
198 delete fTreeE;
199 delete fTreeR;
8494b010 200 if (fModules) {
201 fModules->Delete();
202 delete fModules;
fe4da5cc 203 }
204 if (fParticles) {
205 fParticles->Delete();
206 delete fParticles;
207 }
1cedd08a 208 delete fHitLists;
fe4da5cc 209}
210
211//_____________________________________________________________________________
212void AliRun::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
213{
214 //
215 // Add a hit to detector id
216 //
8494b010 217 TObjArray &dets = *fModules;
218 if(dets[id]) ((AliModule*) dets[id])->AddHit(track,vol,hits);
fe4da5cc 219}
220
221//_____________________________________________________________________________
222void AliRun::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
223{
224 //
225 // Add digit to detector id
226 //
8494b010 227 TObjArray &dets = *fModules;
228 if(dets[id]) ((AliModule*) dets[id])->AddDigit(tracks,digits);
fe4da5cc 229}
230
231//_____________________________________________________________________________
232void AliRun::Browse(TBrowser *b)
233{
234 //
235 // Called when the item "Run" is clicked on the left pane
236 // of the Root browser.
237 // It displays the Root Trees and all detectors.
238 //
239 if (fTreeK) b->Add(fTreeK,fTreeK->GetName());
240 if (fTreeH) b->Add(fTreeH,fTreeH->GetName());
241 if (fTreeD) b->Add(fTreeD,fTreeD->GetName());
242 if (fTreeE) b->Add(fTreeE,fTreeE->GetName());
243 if (fTreeR) b->Add(fTreeR,fTreeR->GetName());
244
8494b010 245 TIter next(fModules);
246 AliModule *detector;
247 while((detector = (AliModule*)next())) {
fe4da5cc 248 b->Add(detector,detector->GetName());
249 }
250}
251
252//_____________________________________________________________________________
253void AliRun::Build()
254{
255 //
256 // Initialize Alice geometry
257 // Dummy routine
258 //
259}
260
261//_____________________________________________________________________________
262void AliRun::BuildSimpleGeometry()
263{
264 //
265 // Create a simple TNode geometry used by Root display engine
266 //
267 // Initialise geometry
268 //
269 fGeometry = new TGeometry("AliceGeom","Galice Geometry for Hits");
270 new TMaterial("void","Vacuum",0,0,0); //Everything is void
271 TBRIK *brik = new TBRIK("S_alice","alice volume","void",2000,2000,3000);
272 brik->SetVisibility(0);
273 new TNode("alice","alice","S_alice");
274}
275
276//_____________________________________________________________________________
277void AliRun::CleanDetectors()
278{
279 //
280 // Clean Detectors at the end of event
281 //
8494b010 282 TIter next(fModules);
283 AliModule *detector;
284 while((detector = (AliModule*)next())) {
fe4da5cc 285 detector->FinishEvent();
286 }
287}
288
289//_____________________________________________________________________________
290void AliRun::CleanParents()
291{
292 //
293 // Clean Particles stack.
1578254f 294 // Set parent/daughter relations
fe4da5cc 295 //
296 TClonesArray &particles = *(gAlice->Particles());
1578254f 297 TParticle *part;
fe4da5cc 298 int i;
299 for(i=0; i<fNtrack; i++) {
1578254f 300 part = (TParticle *)particles.UncheckedAt(i);
301 if(!part->TestBit(Daughters_Bit)) {
302 part->SetFirstDaughter(-1);
303 part->SetLastDaughter(-1);
fe4da5cc 304 }
305 }
306}
307
308//_____________________________________________________________________________
309Int_t AliRun::DistancetoPrimitive(Int_t, Int_t)
310{
311 //
312 // Return the distance from the mouse to the AliRun object
313 // Dummy routine
314 //
315 return 9999;
316}
317
318//_____________________________________________________________________________
319void AliRun::DumpPart (Int_t i)
320{
321 //
322 // Dumps particle i in the stack
323 //
324 TClonesArray &particles = *fParticles;
1578254f 325 ((TParticle*) particles[i])->Print();
fe4da5cc 326}
327
328//_____________________________________________________________________________
329void AliRun::DumpPStack ()
330{
331 //
332 // Dumps the particle stack
333 //
334 TClonesArray &particles = *fParticles;
335 printf(
336 "\n\n=======================================================================\n");
337 for (Int_t i=0;i<fNtrack;i++)
338 {
1578254f 339 printf("-> %d ",i); ((TParticle*) particles[i])->Print();
fe4da5cc 340 printf("--------------------------------------------------------------\n");
341 }
342 printf(
343 "\n=======================================================================\n\n");
344}
345
346//_____________________________________________________________________________
347void AliRun::SetField(Int_t type, Int_t version, Float_t scale,
348 Float_t maxField, char* filename)
349{
350 //
351 // Set magnetic field parameters
352 // type Magnetic field transport flag 0=no field, 2=helix, 3=Runge Kutta
353 // version Magnetic field map version (only 1 active now)
354 // scale Scale factor for the magnetic field
355 // maxField Maximum value for the magnetic field
356
357 //
358 // --- Sanity check on mag field flags
359 if(type<0 || type > 2) {
23370b7a 360 Warning("SetField",
361 "Invalid magnetic field flag: %5d; Helix tracking chosen instead\n"
fe4da5cc 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 {
23370b7a 372 Warning("SetField","Invalid map %d\n",version);
fe4da5cc 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
6c9704e6 395 // static Int_t count=0;
396 // const Int_t times=10;
fe4da5cc 397 // This primary is finished, purify stack
398 gAlice->PurifyKine();
399
400 // Write out hits if any
401 if (gAlice->TreeH()) {
402 gAlice->TreeH()->Fill();
403 }
404
405 // Reset Hits info
406 gAlice->ResetHits();
a8f1fb7c 407
408 //
409 // if(++count%times==1) gObjectTable->Print();
fe4da5cc 410}
411
412//_____________________________________________________________________________
413void AliRun::FinishEvent()
414{
415 //
416 // Called at the end of the event.
417 //
7fb01480 418
fe4da5cc 419 //Update the energy deposit tables
420 Int_t i;
421 for(i=0;i<sEventEnergy.GetSize();i++) {
422 sSummEnergy[i]+=sEventEnergy[i];
423 sSum2Energy[i]+=sEventEnergy[i]*sEventEnergy[i];
424 }
425 sEventEnergy.Reset();
426
427 // Clean detector information
428 CleanDetectors();
429
430 // Write out the kinematics
431 if (fTreeK) {
432 CleanParents();
433 fTreeK->Fill();
434 }
435
436 // Write out the digits
437 if (fTreeD) {
438 fTreeD->Fill();
439 ResetDigits();
440 }
441
442 // Write out reconstructed clusters
443 if (fTreeR) {
444 fTreeR->Fill();
445 }
446
447 // Write out the event Header information
448 if (fTreeE) fTreeE->Fill();
449
450 // Reset stack info
451 ResetStack();
452
453 // Write Tree headers
59fe9bd2 454 // Int_t ievent = fHeader.GetEvent();
455 // char hname[30];
456 // sprintf(hname,"TreeK%d",ievent);
457 if (fTreeK) fTreeK->Write();
458 // sprintf(hname,"TreeH%d",ievent);
459 if (fTreeH) fTreeH->Write();
460 // sprintf(hname,"TreeD%d",ievent);
461 if (fTreeD) fTreeD->Write();
462 // sprintf(hname,"TreeR%d",ievent);
463 if (fTreeR) fTreeR->Write();
fe4da5cc 464}
465
466//_____________________________________________________________________________
467void AliRun::FinishRun()
468{
469 //
470 // Called at the end of the run.
471 //
472
473 // Clean detector information
8494b010 474 TIter next(fModules);
475 AliModule *detector;
476 while((detector = (AliModule*)next())) {
fe4da5cc 477 detector->FinishRun();
478 }
479
480 //Output energy summary tables
481 EnergySummary();
482
483 // file is retrieved from whatever tree
484 TFile *File = 0;
485 if (fTreeK) File = fTreeK->GetCurrentFile();
486 if ((!File) && (fTreeH)) File = fTreeH->GetCurrentFile();
487 if ((!File) && (fTreeD)) File = fTreeD->GetCurrentFile();
488 if ((!File) && (fTreeE)) File = fTreeE->GetCurrentFile();
489 if( NULL==File ) {
490 Error("FinishRun","There isn't root file!");
491 exit(1);
492 }
493 File->cd();
494 fTreeE->Write();
495
496 // Clean tree information
497 delete fTreeK; fTreeK = 0;
498 delete fTreeH; fTreeH = 0;
499 delete fTreeD; fTreeD = 0;
500 delete fTreeR; fTreeR = 0;
501 delete fTreeE; fTreeE = 0;
502
503 // Write AliRun info and all detectors parameters
504 Write();
505
506 // Close output file
507 File->Write();
fe4da5cc 508}
509
510//_____________________________________________________________________________
511void AliRun::FlagTrack(Int_t track)
512{
513 //
514 // Flags a track and all its family tree to be kept
515 //
516 int curr;
1578254f 517 TParticle *particle;
fe4da5cc 518
519 curr=track;
520 while(1) {
1578254f 521 particle=(TParticle*)fParticles->UncheckedAt(curr);
fe4da5cc 522
523 // If the particle is flagged the three from here upward is saved already
524 if(particle->TestBit(Keep_Bit)) return;
525
526 // Save this particle
527 particle->SetBit(Keep_Bit);
528
529 // Move to father if any
1578254f 530 if((curr=particle->GetFirstMother())==-1) return;
fe4da5cc 531 }
532}
533
534//_____________________________________________________________________________
535void AliRun::EnergySummary()
536{
537 //
538 // Print summary of deposited energy
539 //
540
fe4da5cc 541 Int_t ndep=0;
542 Float_t edtot=0;
543 Float_t ed, ed2;
544 Int_t kn, i, left, j, id;
545 const Float_t zero=0;
546 Int_t ievent=fHeader.GetEvent()+1;
547 //
548 // Energy loss information
549 if(ievent) {
550 printf("***************** Energy Loss Information per event (GEV) *****************\n");
551 for(kn=1;kn<sEventEnergy.GetSize();kn++) {
552 ed=sSummEnergy[kn];
553 if(ed>0) {
554 sEventEnergy[ndep]=kn;
555 if(ievent>1) {
556 ed=ed/ievent;
557 ed2=sSum2Energy[kn];
558 ed2=ed2/ievent;
559 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,zero))/ed;
560 } else
561 ed2=99;
562 sSummEnergy[ndep]=ed;
563 sSum2Energy[ndep]=TMath::Min((Float_t) 99.,TMath::Max(ed2,zero));
564 edtot+=ed;
565 ndep++;
566 }
567 }
568 for(kn=0;kn<(ndep-1)/3+1;kn++) {
569 left=ndep-kn*3;
570 for(i=0;i<(3<left?3:left);i++) {
571 j=kn*3+i;
572 id=Int_t (sEventEnergy[j]+0.1);
cfce8870 573 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),sSummEnergy[j],sSum2Energy[j]);
fe4da5cc 574 }
575 printf("\n");
576 }
577 //
578 // Relative energy loss in different detectors
579 printf("******************** Relative Energy Loss per event ********************\n");
580 printf("Total energy loss per event %10.3f GeV\n",edtot);
581 for(kn=0;kn<(ndep-1)/5+1;kn++) {
582 left=ndep-kn*5;
583 for(i=0;i<(5<left?5:left);i++) {
584 j=kn*5+i;
585 id=Int_t (sEventEnergy[j]+0.1);
cfce8870 586 printf(" %s %10.3f%%;",gMC->VolName(id),100*sSummEnergy[j]/edtot);
fe4da5cc 587 }
588 printf("\n");
589 }
590 for(kn=0;kn<75;kn++) printf("*");
591 printf("\n");
592 }
593 //
594 // Reset the TArray's
b13db077 595 // sEventEnergy.Set(0);
596 // sSummEnergy.Set(0);
597 // sSum2Energy.Set(0);
fe4da5cc 598}
599
600//_____________________________________________________________________________
8494b010 601AliModule *AliRun::GetModule(const char *name)
fe4da5cc 602{
603 //
604 // Return pointer to detector from name
605 //
8494b010 606 return (AliModule*)fModules->FindObject(name);
fe4da5cc 607}
608
a68348e9 609//_____________________________________________________________________________
610AliDetector *AliRun::GetDetector(const char *name)
611{
612 //
613 // Return pointer to detector from name
614 //
615 return (AliDetector*)fModules->FindObject(name);
616}
617
fe4da5cc 618//_____________________________________________________________________________
8494b010 619Int_t AliRun::GetModuleID(const char *name)
fe4da5cc 620{
621 //
622 // Return galice internal detector identifier from name
623 //
23370b7a 624 Int_t i=-1;
625 TObject *mod=fModules->FindObject(name);
626 if(mod) i=fModules->IndexOf(mod);
627 return i;
fe4da5cc 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 //
fe4da5cc 637
638 // Reset existing structures
639 ResetStack();
640 ResetHits();
641 ResetDigits();
642
643 // Delete Trees already connected
644 if (fTreeK) delete fTreeK;
645 if (fTreeH) delete fTreeH;
646 if (fTreeD) delete fTreeD;
647 if (fTreeR) delete fTreeR;
59fe9bd2 648
649 // Get header from file
650 if(fTreeE) fTreeE->GetEntry(event);
651 else Error("GetEvent","Cannot file Header Tree\n");
fe4da5cc 652
653 // Get Kine Tree from file
654 char treeName[20];
655 sprintf(treeName,"TreeK%d",event);
656 fTreeK = (TTree*)gDirectory->Get(treeName);
657 if (fTreeK) fTreeK->SetBranchAddress("Particles", &fParticles);
23370b7a 658 else Error("GetEvent","cannot find Kine Tree for event:%d\n",event);
fe4da5cc 659
660 // Get Hits Tree header from file
661 sprintf(treeName,"TreeH%d",event);
662 fTreeH = (TTree*)gDirectory->Get(treeName);
663 if (!fTreeH) {
23370b7a 664 Error("GetEvent","cannot find Hits Tree for event:%d\n",event);
fe4da5cc 665 }
666
667 // Get Digits Tree header from file
668 sprintf(treeName,"TreeD%d",event);
669 fTreeD = (TTree*)gDirectory->Get(treeName);
670 if (!fTreeD) {
07a68c1d 671 Warning("GetEvent","cannot find Digits Tree for event:%d\n",event);
fe4da5cc 672 }
673
674
675 // Get Reconstruct Tree header from file
676 sprintf(treeName,"TreeR%d",event);
677 fTreeR = (TTree*)gDirectory->Get(treeName);
678 if (!fTreeR) {
679 // printf("WARNING: cannot find Reconstructed Tree for event:%d\n",event);
680 }
681
682 // Set Trees branch addresses
8494b010 683 TIter next(fModules);
684 AliModule *detector;
685 while((detector = (AliModule*)next())) {
fe4da5cc 686 detector->SetTreeAddress();
687 }
688
689 if (fTreeK) fTreeK->GetEvent(0);
690 fNtrack = Int_t (fParticles->GetEntries());
691 return fNtrack;
692}
693
694//_____________________________________________________________________________
695TGeometry *AliRun::GetGeometry()
696{
697 //
698 // Import Alice geometry from current file
699 // Return pointer to geometry object
700 //
701 if (!fGeometry) fGeometry = (TGeometry*)gDirectory->Get("AliceGeom");
702 //
703 // Unlink and relink nodes in detectors
704 // This is bad and there must be a better way...
705 //
fe4da5cc 706
8494b010 707 TIter next(fModules);
708 AliModule *detector;
709 while((detector = (AliModule*)next())) {
fe4da5cc 710 detector->SetTreeAddress();
711 TList *dnodes=detector->Nodes();
712 Int_t j;
713 TNode *node, *node1;
714 for ( j=0; j<dnodes->GetSize(); j++) {
715 node = (TNode*) dnodes->At(j);
52d0ab00 716 node1 = fGeometry->GetNode(node->GetName());
fe4da5cc 717 dnodes->Remove(node);
718 dnodes->AddAt(node1,j);
719 }
720 }
721 return fGeometry;
722}
723
724//_____________________________________________________________________________
725void AliRun::GetNextTrack(Int_t &mtrack, Int_t &ipart, Float_t *pmom,
726 Float_t &e, Float_t *vpos, Float_t *polar,
727 Float_t &tof)
728{
729 //
730 // Return next track from stack of particles
731 //
a8f1fb7c 732 TVector3 pol;
fe4da5cc 733 fCurrent=-1;
1578254f 734 TParticle *track;
fe4da5cc 735 for(Int_t i=fNtrack-1; i>=0; i--) {
1578254f 736 track=(TParticle*) fParticles->UncheckedAt(i);
fe4da5cc 737 if(!track->TestBit(Done_Bit)) {
738 //
739 // The track has not yet been processed
740 fCurrent=i;
1578254f 741 ipart=track->GetPdgCode();
742 pmom[0]=track->Px();
743 pmom[1]=track->Py();
744 pmom[2]=track->Pz();
745 e =track->Energy();
746 vpos[0]=track->Vx();
747 vpos[1]=track->Vy();
748 vpos[2]=track->Vz();
a8f1fb7c 749 track->GetPolarisation(pol);
750 polar[0]=pol.X();
751 polar[1]=pol.Y();
752 polar[2]=pol.Z();
1578254f 753 tof=track->T();
fe4da5cc 754 track->SetBit(Done_Bit);
755 break;
756 }
757 }
758 mtrack=fCurrent;
759 //
760 // stop and start timer when we start a primary track
761 Int_t nprimaries = fHeader.GetNprimary();
762 if (fCurrent >= nprimaries) return;
763 if (fCurrent < nprimaries-1) {
764 fTimer.Stop();
1578254f 765 track=(TParticle*) fParticles->UncheckedAt(fCurrent+1);
766 // track->SetProcessTime(fTimer.CpuTime());
fe4da5cc 767 }
768 fTimer.Start();
769}
770
771//_____________________________________________________________________________
772Int_t AliRun::GetPrimary(Int_t track)
773{
774 //
775 // return number of primary that has generated track
776 //
777 int current, parent;
1578254f 778 TParticle *part;
fe4da5cc 779 //
780 parent=track;
781 while (1) {
782 current=parent;
1578254f 783 part = (TParticle *)fParticles->UncheckedAt(current);
784 parent=part->GetFirstMother();
fe4da5cc 785 if(parent<0) return current;
786 }
787}
788
789//_____________________________________________________________________________
790void AliRun::Init(const char *setup)
791{
792 //
793 // Initialize the Alice setup
794 //
795
796 gROOT->LoadMacro(setup);
797 gInterpreter->ProcessLine("Config();");
798
cfce8870 799 gMC->DefineParticles(); //Create standard MC particles
fe4da5cc 800
801 TObject *objfirst, *objlast;
802
23370b7a 803 fNdets = fModules->GetLast()+1;
804
fe4da5cc 805 //
806 //=================Create Materials, geometry, histograms, etc
8494b010 807 TIter next(fModules);
808 AliModule *detector;
809 while((detector = (AliModule*)next())) {
fe4da5cc 810 detector->SetTreeAddress();
811 objlast = gDirectory->GetList()->Last();
812
813 // Initialise detector materials, geometry, histograms,etc
814 detector->CreateMaterials();
815 detector->CreateGeometry();
816 detector->BuildGeometry();
817 detector->Init();
818
819 // Add Detector histograms in Detector list of histograms
820 if (objlast) objfirst = gDirectory->GetList()->After(objlast);
821 else objfirst = gDirectory->GetList()->First();
822 while (objfirst) {
823 detector->Histograms()->Add(objfirst);
824 objfirst = gDirectory->GetList()->After(objfirst);
825 }
826 }
827 SetTransPar(); //Read the cuts for all materials
828
829 MediaTable(); //Build the special IMEDIA table
830
b13db077 831 //Terminate building of geometry
832 printf("%p\n",gVMC);
833 gVMC->FinishGeometry();
fe4da5cc 834
835 //Initialise geometry deposition table
cfce8870 836 sEventEnergy.Set(gMC->NofVolumes()+1);
837 sSummEnergy.Set(gMC->NofVolumes()+1);
838 sSum2Energy.Set(gMC->NofVolumes()+1);
fe4da5cc 839
fe4da5cc 840 //Compute cross-sections
b13db077 841 gVMC->BuildPhysics();
fe4da5cc 842
843 //Write Geometry object to current file.
844 fGeometry->Write();
845
846 fInitDone = kTRUE;
847}
848
849//_____________________________________________________________________________
850void AliRun::MediaTable()
851{
852 //
853 // Built media table to get from the media number to
854 // the detector id
855 //
ad51aeb0 856 Int_t kz, nz, idt, lz, i, k, ind;
857 // Int_t ibeg;
fe4da5cc 858 TObjArray &dets = *gAlice->Detectors();
8494b010 859 AliModule *det;
fe4da5cc 860 //
861 // For all detectors
862 for (kz=0;kz<fNdets;kz++) {
863 // If detector is defined
8494b010 864 if((det=(AliModule*) dets[kz])) {
ad51aeb0 865 TArrayI &idtmed = *(det->GetIdtmed());
866 for(nz=0;nz<100;nz++) {
fe4da5cc 867 // Find max and min material number
ad51aeb0 868 if((idt=idtmed[nz])) {
fe4da5cc 869 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
870 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
871 }
872 }
873 if(det->LoMedium() > det->HiMedium()) {
874 det->LoMedium() = 0;
875 det->HiMedium() = 0;
876 } else {
877 if(det->HiMedium() > fImedia->GetSize()) {
ad51aeb0 878 Error("MediaTable","Increase fImedia from %d to %d",
879 fImedia->GetSize(),det->HiMedium());
fe4da5cc 880 return;
881 }
882 // Tag all materials in rage as belonging to detector kz
883 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
884 (*fImedia)[lz]=kz;
885 }
886 }
887 }
888 }
889 //
890 // Print summary table
891 printf(" Traking media ranges:\n");
892 for(i=0;i<(fNdets-1)/6+1;i++) {
893 for(k=0;k< (6<fNdets-i*6?6:fNdets-i*6);k++) {
894 ind=i*6+k;
8494b010 895 det=(AliModule*)dets[ind];
fe4da5cc 896 if(det)
897 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
898 det->HiMedium());
899 else
900 printf(" %6s: %3d -> %3d;","NULL",0,0);
901 }
902 printf("\n");
903 }
904}
905
906//____________________________________________________________________________
907void AliRun::SetGenerator(AliGenerator *generator)
908{
909 //
910 // Load the event generator
911 //
b13db077 912 if(fGenerator)
913 Warning("SetGenerator","Replacing generator %s with %s\n",
914 fGenerator->GetName(),generator->GetName());
915 fGenerator = generator;
fe4da5cc 916}
917
918//____________________________________________________________________________
919void AliRun::SetTransPar(char* filename)
920{
921 //
922 // Read filename to set the transport parameters
923 //
924
fe4da5cc 925
926 const Int_t ncuts=10;
927 const Int_t nflags=11;
928 const Int_t npars=ncuts+nflags;
929 const char pars[npars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
930 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
931 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
932 "MULS","PAIR","PHOT","RAYL"};
933 char line[256];
ad51aeb0 934 char detName[7];
fe4da5cc 935 char* filtmp;
936 Float_t cut[ncuts];
937 Int_t flag[nflags];
938 Int_t i, itmed, iret, ktmed, kz;
939 FILE *lun;
940 //
941 // See whether the file is there
942 filtmp=gSystem->ExpandPathName(filename);
943 lun=fopen(filtmp,"r");
944 delete [] filtmp;
945 if(!lun) {
ad51aeb0 946 Warning("SetTransPar","File %s does not exist!\n",filename);
fe4da5cc 947 return;
948 }
949 //
950 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
951 printf(" *%59s\n","*");
952 printf(" * Please check carefully what you are doing!%10s\n","*");
953 printf(" *%59s\n","*");
954 //
955 while(1) {
956 // Initialise cuts and flags
957 for(i=0;i<ncuts;i++) cut[i]=-99;
958 for(i=0;i<nflags;i++) flag[i]=-99;
959 itmed=0;
960 for(i=0;i<256;i++) line[i]='\0';
961 // Read up to the end of line excluded
962 iret=fscanf(lun,"%[^\n]",line);
963 if(iret<0) {
964 //End of file
965 fclose(lun);
966 printf(" *%59s\n","*");
967 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
968 return;
969 }
970 // Read the end of line
971 fscanf(lun,"%*c");
972 if(!iret) continue;
973 if(line[0]=='*') continue;
974 // Read the numbers
ad51aeb0 975 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",
976 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
977 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
978 &flag[8],&flag[9],&flag[10]);
fe4da5cc 979 if(!iret) continue;
980 if(iret<0) {
981 //reading error
ad51aeb0 982 Warning("SetTransPar","Error reading file %s\n",filename);
fe4da5cc 983 continue;
984 }
ad51aeb0 985 // Check that the module exist
986 AliModule *mod = GetModule(detName);
987 if(mod) {
988 // Get the array of media numbers
989 TArrayI &idtmed = *mod->GetIdtmed();
990 // Check that the tracking medium code is valid
991 if(0<=itmed && itmed < 100) {
992 ktmed=idtmed[itmed];
993 if(!ktmed) {
994 Warning("SetTransPar","Invalid tracking medium code %d for %s\n",itmed,mod->GetName());
995 continue;
fe4da5cc 996 }
ad51aeb0 997 // Set energy thresholds
998 for(kz=0;kz<ncuts;kz++) {
999 if(cut[kz]>=0) {
23370b7a 1000 printf(" * %-6s set to %10.3E for tracking medium code %4d for %s\n",
ad51aeb0 1001 pars[kz],cut[kz],itmed,mod->GetName());
cfce8870 1002 gMC->Gstpar(ktmed,pars[kz],cut[kz]);
ad51aeb0 1003 }
fe4da5cc 1004 }
ad51aeb0 1005 // Set transport mechanisms
1006 for(kz=0;kz<nflags;kz++) {
1007 if(flag[kz]>=0) {
1008 printf(" * %-6s set to %10d for tracking medium code %4d for %s\n",
1009 pars[ncuts+kz],flag[kz],itmed,mod->GetName());
cfce8870 1010 gMC->Gstpar(ktmed,pars[ncuts+kz],Float_t(flag[kz]));
ad51aeb0 1011 }
1012 }
1013 } else {
1014 Warning("SetTransPar","Invalid medium code %d *\n",itmed);
1015 continue;
fe4da5cc 1016 }
1017 } else {
ad51aeb0 1018 Warning("SetTransPar","Module %s not present\n",detName);
fe4da5cc 1019 continue;
1020 }
1021 }
1022}
1023
1024//_____________________________________________________________________________
1025void AliRun::MakeTree(Option_t *option)
1026{
1027 //
1028 // Create the ROOT trees
1029 // Loop on all detectors to create the Root branch (if any)
1030 //
1031
b13db077 1032 char hname[30];
fe4da5cc 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 //
b13db077 1041 if (K && !fTreeK) {
1042 sprintf(hname,"TreeK%d",fEvent);
1043 fTreeK = new TTree(hname,"Kinematics");
1044 // Create a branch for particles
1045 fTreeK->Branch("Particles",&fParticles,4000);
1046 }
1047 if (H && !fTreeH) {
1048 sprintf(hname,"TreeH%d",fEvent);
1049 fTreeH = new TTree(hname,"Hits");
1050 fTreeH->SetAutoSave(1000000000); //no autosave
1051 }
1052 if (D && !fTreeD) {
1053 sprintf(hname,"TreeD%d",fEvent);
1054 fTreeD = new TTree(hname,"Digits");
1055 }
1056 if (R && !fTreeR) {
1057 sprintf(hname,"TreeR%d",fEvent);
1058 fTreeR = new TTree(hname,"Reconstruction");
1059 }
1060 if (E && !fTreeE) {
1061 fTreeE = new TTree("TE","Header");
1062 // Create a branch for Header
1063 fTreeE->Branch("Header","AliHeader",&header,4000);
1064 }
fe4da5cc 1065 //
1066 // Create a branch for hits/digits for each detector
1067 // Each branch is a TClonesArray. Each data member of the Hits classes
1068 // will be in turn a subbranch of the detector master branch
8494b010 1069 TIter next(fModules);
1070 AliModule *detector;
1071 while((detector = (AliModule*)next())) {
fe4da5cc 1072 if (H || D || R) detector->MakeBranch(option);
1073 }
fe4da5cc 1074}
1075
1076//_____________________________________________________________________________
1077Int_t AliRun::PurifyKine(Int_t lastSavedTrack, Int_t nofTracks)
1078{
1079 //
1080 // PurifyKine with external parameters
1081 //
1082 fHgwmk = lastSavedTrack;
1083 fNtrack = nofTracks;
1084 PurifyKine();
1085 return fHgwmk;
1086}
1087
1088//_____________________________________________________________________________
1089void AliRun::PurifyKine()
1090{
1091 //
1092 // Compress kinematic tree keeping only flagged particles
1093 // and renaming the particle id's in all the hits
1094 //
1095 TClonesArray &particles = *fParticles;
1096 int nkeep=fHgwmk+1, parent, i;
1578254f 1097 TParticle *part, *partnew, *father;
fe4da5cc 1098 int *map = new int[particles.GetEntries()];
1099
1100 // Save in Header total number of tracks before compression
1101 fHeader.SetNtrack(fHeader.GetNtrack()+fNtrack-fHgwmk);
1102
1103 // Preset map, to be removed later
1104 for(i=0; i<fNtrack; i++) {
1105 if(i<=fHgwmk) map[i]=i ; else map[i] = -99 ;}
1106 // Second pass, build map between old and new numbering
1107 for(i=fHgwmk+1; i<fNtrack; i++) {
1578254f 1108 part = (TParticle *)particles.UncheckedAt(i);
fe4da5cc 1109 if(part->TestBit(Keep_Bit)) {
1110
1111 // This particle has to be kept
1112 map[i]=nkeep;
1113 if(i!=nkeep) {
1114
1115 // Old and new are different, have to copy
1578254f 1116 partnew = (TParticle *)particles.UncheckedAt(nkeep);
fe4da5cc 1117 *partnew = *part;
1118 } else partnew = part;
1119
1120 // as the parent is always *before*, it must be already
1121 // in place. This is what we are checking anyway!
1578254f 1122 if((parent=partnew->GetFirstMother())>fHgwmk) {
fe4da5cc 1123 if(map[parent]==-99) printf("map[%d] = -99!\n",parent);
1578254f 1124 partnew->SetFirstMother(map[parent]);
fe4da5cc 1125 }
1126 nkeep++;
1127 }
1128 }
1129 fNtrack=nkeep;
1130
1578254f 1131 // Fix daughters information
ae23d366 1132 for (i=0; i<fNtrack; i++) {
1578254f 1133 part = (TParticle *)particles.UncheckedAt(i);
1134 parent = part->GetFirstMother();
ae23d366 1135 if(parent>=0) {
1136 father = (TParticle *)particles.UncheckedAt(parent);
1137 if(father->TestBit(Daughters_Bit)) {
fe4da5cc 1138
ae23d366 1139 if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
1140 if(i>father->GetLastDaughter()) father->SetLastDaughter(i);
1141 } else {
1142 // Iitialise daughters info for first pass
1143 father->SetFirstDaughter(i);
1144 father->SetLastDaughter(i);
1145 father->SetBit(Daughters_Bit);
1146 }
fe4da5cc 1147 }
1148 }
1149
1cedd08a 1150#ifdef old
fe4da5cc 1151 // Now loop on all detectors and reset the hits
b13db077 1152 AliHit *OneHit;
8494b010 1153 TIter next(fModules);
1154 AliModule *detector;
1155 while((detector = (AliModule*)next())) {
fe4da5cc 1156 if (!detector->Hits()) continue;
1157 TClonesArray &vHits=*(detector->Hits());
1158 if(vHits.GetEntries() != detector->GetNhits())
1159 printf("vHits.GetEntries()!=detector->GetNhits(): %d != %d\n",
1160 vHits.GetEntries(),detector->GetNhits());
1161 for (i=0; i<detector->GetNhits(); i++) {
1162 OneHit = (AliHit *)vHits.UncheckedAt(i);
1163 OneHit->SetTrack(map[OneHit->GetTrack()]);
1164 }
1165 }
1cedd08a 1166#else
1167
1168 // Now loop on all registered hit lists
1169 TIter next(fHitLists);
1170 TCollection *hitList;
1171 while((hitList = (TCollection*)next())) {
1172 TIter nexthit(hitList);
1173 AliHit *hit;
1174 while((hit = (AliHit*)nexthit())) {
1175 hit->SetTrack(map[hit->GetTrack()]);
1176 }
1177 }
1178#endif
fe4da5cc 1179
1180 fHgwmk=nkeep-1;
1181 particles.SetLast(fHgwmk);
1182 delete [] map;
1183}
1184
1185//_____________________________________________________________________________
1186void AliRun::Reset(Int_t run, Int_t idevent)
1187{
1188 //
1189 // Reset all Detectors & kinematics & trees
1190 //
59fe9bd2 1191 char hname[30];
1192 //
fe4da5cc 1193 ResetStack();
1194 ResetHits();
1195 ResetDigits();
1196
1197 // Initialise event header
1198 fHeader.Reset(run,idevent);
1199
59fe9bd2 1200 if(fTreeK) {
1201 fTreeK->Reset();
1202 sprintf(hname,"TreeK%d",idevent);
1203 fTreeK->SetName(hname);
1204 }
1205 if(fTreeH) {
1206 fTreeH->Reset();
1207 sprintf(hname,"TreeH%d",idevent);
1208 fTreeH->SetName(hname);
1209 }
1210 if(fTreeD) {
1211 fTreeD->Reset();
1212 sprintf(hname,"TreeD%d",idevent);
1213 fTreeD->SetName(hname);
1214 }
1215 if(fTreeR) {
1216 fTreeR->Reset();
1217 sprintf(hname,"TreeR%d",idevent);
1218 fTreeR->SetName(hname);
1219 }
fe4da5cc 1220}
1221
1222//_____________________________________________________________________________
1223void AliRun::ResetDigits()
1224{
1225 //
1226 // Reset all Detectors digits
1227 //
8494b010 1228 TIter next(fModules);
1229 AliModule *detector;
1230 while((detector = (AliModule*)next())) {
fe4da5cc 1231 detector->ResetDigits();
1232 }
1233}
1234
1235//_____________________________________________________________________________
1236void AliRun::ResetHits()
1237{
1238 //
1239 // Reset all Detectors hits
1240 //
8494b010 1241 TIter next(fModules);
1242 AliModule *detector;
1243 while((detector = (AliModule*)next())) {
fe4da5cc 1244 detector->ResetHits();
1245 }
1246}
1247
1248//_____________________________________________________________________________
1249void AliRun::ResetPoints()
1250{
1251 //
1252 // Reset all Detectors points
1253 //
8494b010 1254 TIter next(fModules);
1255 AliModule *detector;
1256 while((detector = (AliModule*)next())) {
fe4da5cc 1257 detector->ResetPoints();
1258 }
1259}
1260
1261//_____________________________________________________________________________
1262void AliRun::Run(Int_t nevent, const char *setup)
1263{
1264 //
1265 // Main function to be called to process a galice run
1266 // example
1267 // Root > gAlice.Run();
1268 // a positive number of events will cause the finish routine
1269 // to be called
1270 //
1271
1272 Int_t i, todo;
1273 // check if initialisation has been done
1274 if (!fInitDone) Init(setup);
1275
fe4da5cc 1276 // Create the Root Tree with one branch per detector
b13db077 1277 gAlice->MakeTree("KHDER");
fe4da5cc 1278
1279 todo = TMath::Abs(nevent);
1280 for (i=0; i<todo; i++) {
1281 // Process one run (one run = one event)
1282 gAlice->Reset(fRun, fEvent);
b13db077 1283 gVMC->ProcessEvent();
fe4da5cc 1284 gAlice->FinishEvent();
1285 fEvent++;
1286 }
1287
1288 // End of this run, close files
1289 if(nevent>0) gAlice->FinishRun();
1290}
1291
1292//_____________________________________________________________________________
1293void AliRun::RunLego(const char *setup,Int_t ntheta,Float_t themin,
1294 Float_t themax,Int_t nphi,Float_t phimin,Float_t phimax,
1295 Float_t rmin,Float_t rmax,Float_t zmax)
1296{
1297 //
1298 // Generates lego plots of:
1299 // - radiation length map phi vs theta
1300 // - radiation length map phi vs eta
1301 // - interaction length map
1302 // - g/cm2 length map
1303 //
1304 // ntheta bins in theta, eta
1305 // themin minimum angle in theta (degrees)
1306 // themax maximum angle in theta (degrees)
1307 // nphi bins in phi
1308 // phimin minimum angle in phi (degrees)
1309 // phimax maximum angle in phi (degrees)
1310 // rmin minimum radius
1311 // rmax maximum radius
1312 //
1313 //
1314 // The number of events generated = ntheta*nphi
1315 // run input parameters in macro setup (default="Config.C")
1316 //
1317 // Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
1318 //Begin_Html
1319 /*
1439f98e 1320 <img src="picts/AliRunLego1.gif">
fe4da5cc 1321 */
1322 //End_Html
1323 //Begin_Html
1324 /*
1439f98e 1325 <img src="picts/AliRunLego2.gif">
fe4da5cc 1326 */
1327 //End_Html
1328 //Begin_Html
1329 /*
1439f98e 1330 <img src="picts/AliRunLego3.gif">
fe4da5cc 1331 */
1332 //End_Html
1333 //
1334
1335 // check if initialisation has been done
1336 if (!fInitDone) Init(setup);
b13db077 1337
1338 //Create Lego object
1339 fLego = new AliLego("lego",ntheta,themin,themax,nphi,phimin,phimax,rmin,rmax,zmax);
1340
1341 //Run Lego Object
fe4da5cc 1342 fLego->Run();
1343
1344 // Create only the Root event Tree
1345 gAlice->MakeTree("E");
1346
1347 // End of this run, close files
1348 gAlice->FinishRun();
1349}
1350
1351//_____________________________________________________________________________
1352void AliRun::SetCurrentTrack(Int_t track)
1353{
1354 //
1355 // Set current track number
1356 //
1357 fCurrent = track;
1358}
1359
1360//_____________________________________________________________________________
1578254f 1361void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
fe4da5cc 1362 Float_t *vpos, Float_t *polar, Float_t tof,
1363 const char *mecha, Int_t &ntr, Float_t weight)
1364{
1365 //
1366 // Load a track on the stack
1367 //
1368 // done 0 if the track has to be transported
1369 // 1 if not
1370 // parent identifier of the parent track. -1 for a primary
1578254f 1371 // pdg particle code
fe4da5cc 1372 // pmom momentum GeV/c
1373 // vpos position
1374 // polar polarisation
1375 // tof time of flight in seconds
1376 // mecha production mechanism
1377 // ntr on output the number of the track stored
1378 //
1379 TClonesArray &particles = *fParticles;
1578254f 1380 TParticle *particle;
fe4da5cc 1381 Float_t mass;
1578254f 1382 const Int_t firstdaughter=-1;
1383 const Int_t lastdaughter=-1;
fe4da5cc 1384 const Int_t KS=0;
1578254f 1385 // const Float_t tlife=0;
fe4da5cc 1386
1578254f 1387 //
1388 // Here we get the static mass
1389 // For MC is ok, but a more sophisticated method could be necessary
1390 // if the calculated mass is required
1391 // also, this method is potentially dangerous if the mass
1392 // used in the MC is not the same of the PDG database
1393 //
1394 mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
fe4da5cc 1395 Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
1396 pmom[1]*pmom[1]+pmom[2]*pmom[2]);
1397
1398 //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 1399 //pname,mass,e,fNtrack,pdg,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],KS,mecha);
fe4da5cc 1400
1578254f 1401 particle=new(particles[fNtrack]) TParticle(pdg,KS,parent,-1,firstdaughter,
1402 lastdaughter,pmom[0],pmom[1],pmom[2],
1403 e,vpos[0],vpos[1],vpos[2],tof);
1404 // polar[0],polar[1],polar[2],tof,
1405 // mecha,weight);
1406 ((TParticle*)particles[fNtrack])->SetPolarisation(TVector3(polar[0],polar[1],polar[2]));
1407 ((TParticle*)particles[fNtrack])->SetWeight(weight);
fe4da5cc 1408 if(!done) particle->SetBit(Done_Bit);
1409
1410 if(parent>=0) {
1578254f 1411 particle=(TParticle*) fParticles->UncheckedAt(parent);
1412 particle->SetLastDaughter(fNtrack);
1413 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
fe4da5cc 1414 } else {
1415 //
1416 // This is a primary track. Set high water mark for this event
1417 fHgwmk=fNtrack;
1418 //
1419 // Set also number if primary tracks
1420 fHeader.SetNprimary(fHgwmk+1);
1421 fHeader.SetNtrack(fHgwmk+1);
1422 }
1423 ntr = fNtrack++;
1424}
1425
1426//_____________________________________________________________________________
1427void AliRun::KeepTrack(const Int_t track)
1428{
1429 //
1430 // flags a track to be kept
1431 //
1432 TClonesArray &particles = *fParticles;
1578254f 1433 ((TParticle*)particles[track])->SetBit(Keep_Bit);
fe4da5cc 1434}
1435
1436//_____________________________________________________________________________
1437void AliRun::StepManager(Int_t id) const
1438{
1439 //
1440 // Called at every step during transport
1441 //
1442
fe4da5cc 1443 //
1444 // --- If lego option, do it and leave
b13db077 1445 if (fLego)
fe4da5cc 1446 fLego->StepManager();
b13db077 1447 else {
1448 Int_t copy;
1449 //Update energy deposition tables
1450 sEventEnergy[gMC->CurrentVolID(copy)]+=gMC->Edep();
fe4da5cc 1451
b13db077 1452 //Call the appropriate stepping routine;
1453 AliModule *det = (AliModule*)fModules->At(id);
1454 if(det) det->StepManager();
fe4da5cc 1455 }
fe4da5cc 1456}
1457
fe4da5cc 1458//_____________________________________________________________________________
1459void AliRun::Streamer(TBuffer &R__b)
1460{
1461 //
1462 // Stream an object of class AliRun.
1463 //
1464 if (R__b.IsReading()) {
1465 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
1466 TNamed::Streamer(R__b);
1467 if (!gAlice) gAlice = this;
1468 gROOT->GetListOfBrowsables()->Add(this,"Run");
59fe9bd2 1469 fTreeE = (TTree*)gDirectory->Get("TE");
1470 if (fTreeE) fTreeE->SetBranchAddress("Header", &header);
1471 else Error("Streamer","cannot find Header Tree\n");
fe4da5cc 1472 R__b >> fNtrack;
1473 R__b >> fHgwmk;
1474 R__b >> fDebug;
1475 fHeader.Streamer(R__b);
8494b010 1476 R__b >> fModules;
fe4da5cc 1477 R__b >> fParticles;
1478 R__b >> fField;
1479 // R__b >> fMC;
1480 R__b >> fNdets;
1481 R__b >> fTrRmax;
1482 R__b >> fTrZmax;
1483 R__b >> fGenerator;
59fe9bd2 1484 if(R__v>1) {
1485 R__b >> fPDGDB; //Particle factory object!
1486 fTreeE->GetEntry(0);
1487 } else {
1488 fHeader.SetEvent(0);
1489 fPDGDB = TDatabasePDG::Instance(); //Particle factory object!
1490 }
fe4da5cc 1491 } else {
1492 R__b.WriteVersion(AliRun::IsA());
1493 TNamed::Streamer(R__b);
1494 R__b << fNtrack;
1495 R__b << fHgwmk;
1496 R__b << fDebug;
1497 fHeader.Streamer(R__b);
8494b010 1498 R__b << fModules;
fe4da5cc 1499 R__b << fParticles;
1500 R__b << fField;
1501 // R__b << fMC;
1502 R__b << fNdets;
1503 R__b << fTrRmax;
1504 R__b << fTrZmax;
1505 R__b << fGenerator;
1578254f 1506 R__b << fPDGDB; //Particle factory object!
fe4da5cc 1507 }
1508}