Adding protection - "fake tracks" signed with negative values
[u/mrichter/AliRoot.git] / STEER / AliMC.cxx
CommitLineData
5d12ce38 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/* $Id$ */
17
d0d4a6b3 18// This class is extracted from the AliRun class
19// and contains all the MC-related functionality
20// The number of dependencies has to be reduced...
21// Author: F.Carminati
22// Federico.Carminati@cern.ch
23
c75cfb51 24#include <RVersion.h>
5d12ce38 25#include <TBrowser.h>
7ca4655f 26#include <TClonesArray.h>
27#include <TGeoManager.h>
5d12ce38 28#include <TStopwatch.h>
29#include <TSystem.h>
30#include <TVirtualMC.h>
9ac3aec9 31#include <TParticle.h>
ced249e6 32#include <TROOT.h>
3cac0921 33#include <TFile.h>
b60e0f5e 34
21bf7095 35#include "AliLog.h"
5d12ce38 36#include "AliDetector.h"
37#include "AliGenerator.h"
38#include "AliHeader.h"
39#include "AliLego.h"
40#include "AliMC.h"
41#include "AliMCQA.h"
42#include "AliRun.h"
942a9039 43#include "AliHit.h"
5d12ce38 44#include "AliStack.h"
0561efeb 45#include "AliMagF.h"
5d12ce38 46#include "AliTrackReference.h"
39b3f8ba 47#include "AliSimulation.h"
392808e8 48#include "AliGeomManager.h"
ced249e6 49#include "AliCDBManager.h"
50#include "AliCDBStorage.h"
51#include "AliCDBEntry.h"
5d12ce38 52
53
54ClassImp(AliMC)
55
56//_______________________________________________________________________
57AliMC::AliMC() :
58 fGenerator(0),
59 fEventEnergy(0),
60 fSummEnergy(0),
61 fSum2Energy(0),
62 fTrRmax(1.e10),
63 fTrZmax(1.e10),
90e48c0c 64 fRDecayMax(1.e10),
38ca2ad6 65 fRDecayMin(-1.),
90e48c0c 66 fDecayPdg(0),
5d12ce38 67 fImedia(0),
68 fTransParName("\0"),
69 fMCQA(0),
70 fHitLists(0),
3cac0921 71 fTmpTreeTR(0),
72 fTmpFileTR(0),
73 fTrackReferences(0),
74 fTmpTrackReferences(0)
5d12ce38 75
76{
d0d4a6b3 77 //default constructor
0561efeb 78 DecayLimits();
5d12ce38 79}
80
81//_______________________________________________________________________
82AliMC::AliMC(const char *name, const char *title) :
83 TVirtualMCApplication(name, title),
84 fGenerator(0),
85 fEventEnergy(0),
86 fSummEnergy(0),
87 fSum2Energy(0),
88 fTrRmax(1.e10),
89 fTrZmax(1.e10),
90e48c0c 90 fRDecayMax(1.e10),
38ca2ad6 91 fRDecayMin(-1.),
90e48c0c 92 fDecayPdg(0),
5d12ce38 93 fImedia(new TArrayI(1000)),
94 fTransParName("\0"),
95 fMCQA(0),
96 fHitLists(new TList()),
3cac0921 97 fTmpTreeTR(0),
98 fTmpFileTR(0),
99 fTrackReferences(new TClonesArray("AliTrackReference", 100)),
100 fTmpTrackReferences(new TClonesArray("AliTrackReference", 100))
5d12ce38 101{
d0d4a6b3 102 //constructor
5d12ce38 103 // Set transport parameters
104 SetTransPar();
0561efeb 105 DecayLimits();
5d12ce38 106 // Prepare the tracking medium lists
107 for(Int_t i=0;i<1000;i++) (*fImedia)[i]=-99;
5d12ce38 108}
109
110//_______________________________________________________________________
111AliMC::AliMC(const AliMC &mc) :
112 TVirtualMCApplication(mc),
113 fGenerator(0),
114 fEventEnergy(0),
115 fSummEnergy(0),
116 fSum2Energy(0),
117 fTrRmax(1.e10),
118 fTrZmax(1.e10),
90e48c0c 119 fRDecayMax(1.e10),
38ca2ad6 120 fRDecayMin(-1.),
90e48c0c 121 fDecayPdg(0),
5d12ce38 122 fImedia(0),
123 fTransParName("\0"),
124 fMCQA(0),
125 fHitLists(0),
3cac0921 126 fTmpTreeTR(0),
127 fTmpFileTR(0),
128 fTrackReferences(0),
129 fTmpTrackReferences(0)
5d12ce38 130{
131 //
90e48c0c 132 // Copy constructor for AliMC
5d12ce38 133 //
134 mc.Copy(*this);
135}
136
137//_______________________________________________________________________
138AliMC::~AliMC()
139{
d0d4a6b3 140 //destructor
5d12ce38 141 delete fGenerator;
142 delete fImedia;
143 delete fMCQA;
144 delete fHitLists;
145 // Delete track references
146 if (fTrackReferences) {
147 fTrackReferences->Delete();
148 delete fTrackReferences;
149 fTrackReferences = 0;
150 }
151
3cac0921 152if (fTmpTrackReferences) {
153 fTmpTrackReferences->Delete();
154 delete fTmpTrackReferences;
155 fTmpTrackReferences = 0;
156 }
157
5d12ce38 158}
159
160//_______________________________________________________________________
6c4904c2 161void AliMC::Copy(TObject &) const
5d12ce38 162{
d0d4a6b3 163 //dummy Copy function
21bf7095 164 AliFatal("Not implemented!");
5d12ce38 165}
166
167//_______________________________________________________________________
168void AliMC::ConstructGeometry()
169{
170 //
4a9de4af 171 // Either load geometry from file or create it through usual
172 // loop on detectors. In the first case the method
173 // AliModule::CreateMaterials() only builds fIdtmed and is postponed
174 // at InitGeometry().
5d12ce38 175 //
176
ced249e6 177 if(gAlice->IsRootGeometry()){ //load geometry either from CDB or from file
178 if(gAlice->IsGeomFromCDB()){
179 AliInfo("Loading geometry from CDB default storage");
180 AliCDBPath path("GRP","Geometry","Data");
181 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
182 if(!entry) AliFatal("Unable to load geometry from CDB!");
183 entry->SetOwner(0);
184 gGeoManager = (TGeoManager*) entry->GetObject();
185 if (!gGeoManager) AliFatal("TGeoManager object not found in the specified CDB entry!");
4a9de4af 186 }else{
ced249e6 187 // Load geometry
188 const char *geomfilename = gAlice->GetGeometryFileName();
189 if(gSystem->ExpandPathName(geomfilename)){
190 AliInfo(Form("Loading geometry from file:\n %40s",geomfilename));
191 TGeoManager::Import(geomfilename);
192 }else{
193 AliInfo(Form("Geometry file %40s not found!\n",geomfilename));
194 return;
195 }
4a9de4af 196 }
197 }else{
198 // Create modules, materials, geometry
5d12ce38 199 TStopwatch stw;
200 TIter next(gAlice->Modules());
201 AliModule *detector;
21bf7095 202 AliDebug(1, "Geometry creation:");
5d12ce38 203 while((detector = dynamic_cast<AliModule*>(next()))) {
204 stw.Start();
205 // Initialise detector materials and geometry
206 detector->CreateMaterials();
207 detector->CreateGeometry();
21bf7095 208 AliInfo(Form("%10s R:%.2fs C:%.2fs",
209 detector->GetName(),stw.RealTime(),stw.CpuTime()));
5d12ce38 210 }
4a9de4af 211 }
212
5d12ce38 213}
214
215//_______________________________________________________________________
39b3f8ba 216Bool_t AliMC::MisalignGeometry()
217{
218// Call misalignment code if AliSimulation object was defined.
0717eed2 219
ced249e6 220 if(!gAlice->IsRootGeometry()){
221 //Set alignable volumes for the whole geometry
222 SetAllAlignableVolumes();
223 }
0717eed2 224 // Misalign geometry via AliSimulation instance
39b3f8ba 225 if (!AliSimulation::GetInstance()) return kFALSE;
392808e8 226 AliGeomManager::SetGeometry(gGeoManager);
39b3f8ba 227 return AliSimulation::GetInstance()->MisalignGeometry(gAlice->GetRunLoader());
228}
229
230//_______________________________________________________________________
661663fa 231void AliMC::ConstructOpGeometry()
232{
233 //
234 // Loop all detector modules and call DefineOpticalProperties() method
235 //
236
237 TIter next(gAlice->Modules());
238 AliModule *detector;
239 AliInfo("Optical properties definition");
240 while((detector = dynamic_cast<AliModule*>(next()))) {
241 // Initialise detector optical properties
242 detector->DefineOpticalProperties();
243 }
244}
245
246//_______________________________________________________________________
5d12ce38 247void AliMC::InitGeometry()
248{
249 //
27cf3cdc 250 // Initialize detectors
5d12ce38 251 //
252
4a9de4af 253 AliInfo("Initialisation:");
254 TStopwatch stw;
255 TIter next(gAlice->Modules());
256 AliModule *detector;
257 while((detector = dynamic_cast<AliModule*>(next()))) {
258 stw.Start();
27cf3cdc 259 // Initialise detector geometry
4a9de4af 260 if(gAlice->IsRootGeometry()) detector->CreateMaterials();
261 detector->Init();
4a9de4af 262 AliInfo(Form("%10s R:%.2fs C:%.2fs",
263 detector->GetName(),stw.RealTime(),stw.CpuTime()));
264 }
4787b401 265}
266
267//_______________________________________________________________________
268void AliMC::SetAllAlignableVolumes()
269{
270 //
271 // Add alignable volumes (TGeoPNEntries) looping on all
272 // active modules
273 //
274
275 AliInfo(Form("Setting entries for all alignable volumes of active detectors"));
276 AliModule *detector;
277 TIter next(gAlice->Modules());
278 while((detector = dynamic_cast<AliModule*>(next()))) {
279 detector->AddAlignableVolumes();
280 }
5d12ce38 281}
282
283//_______________________________________________________________________
2057aecb 284void AliMC::GeneratePrimaries()
5d12ce38 285{
286 //
287 // Generate primary particles and fill them in the stack.
288 //
289
290 Generator()->Generate();
291}
292
293//_______________________________________________________________________
294void AliMC::SetGenerator(AliGenerator *generator)
295{
296 //
297 // Load the event generator
298 //
299 if(!fGenerator) fGenerator = generator;
300}
301
302//_______________________________________________________________________
303void AliMC::ResetGenerator(AliGenerator *generator)
304{
305 //
306 // Load the event generator
307 //
9a213abc 308 if(fGenerator) {
309 if(generator) {
21bf7095 310 AliWarning(Form("Replacing generator %s with %s",
9a213abc 311 fGenerator->GetName(),generator->GetName()));
312 }
313 else {
21bf7095 314 AliWarning(Form("Replacing generator %s with NULL",
315 fGenerator->GetName()));
9a213abc 316 }
317 }
5d12ce38 318 fGenerator = generator;
319}
320
321//_______________________________________________________________________
322void AliMC::FinishRun()
323{
324 // Clean generator information
21bf7095 325 AliDebug(1, "fGenerator->FinishRun()");
5d12ce38 326 fGenerator->FinishRun();
327
328 //Output energy summary tables
21bf7095 329 AliDebug(1, "EnergySummary()");
330 ToAliDebug(1, EnergySummary());
5d12ce38 331}
332
333//_______________________________________________________________________
334void AliMC::BeginPrimary()
335{
336 //
337 // Called at the beginning of each primary track
338 //
339
340 // Reset Hits info
341 ResetHits();
342 ResetTrackReferences();
5d12ce38 343}
344
345//_______________________________________________________________________
346void AliMC::PreTrack()
347{
d0d4a6b3 348 // Actions before the track's transport
5d12ce38 349 TObjArray &dets = *gAlice->Modules();
350 AliModule *module;
351
352 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
353 if((module = dynamic_cast<AliModule*>(dets[i])))
354 module->PreTrack();
355
356 fMCQA->PreTrack();
357}
358
359//_______________________________________________________________________
360void AliMC::Stepping()
361{
362 //
363 // Called at every step during transport
364 //
38ca2ad6 365 Int_t id = DetFromMate(gMC->CurrentMedium());
5d12ce38 366 if (id < 0) return;
367
0561efeb 368
369 if ( gMC->IsNewTrack() &&
370 gMC->TrackTime() == 0. &&
38ca2ad6 371 fRDecayMin >= 0. &&
0561efeb 372 fRDecayMax > fRDecayMin &&
373 gMC->TrackPid() == fDecayPdg )
374 {
375 FixParticleDecaytime();
376 }
377
378
379
5d12ce38 380 //
381 // --- If lego option, do it and leave
382 if (gAlice->Lego())
383 gAlice->Lego()->StepManager();
384 else {
385 Int_t copy;
386 //Update energy deposition tables
387 AddEnergyDeposit(gMC->CurrentVolID(copy),gMC->Edep());
ab03c084 388 //
389 // write tracke reference for track which is dissapearing - MI
faad98df 390
391 if (gMC->IsTrackDisappeared() && !(gMC->IsTrackAlive())) {
392 if (gMC->Etot() > 0.05) AddTrackReference(GetCurrentTrackNumber(),
9ac3aec9 393 AliTrackReference::kDisappeared);
faad98df 394
395
b60e0f5e 396 }
faad98df 397
5d12ce38 398 //Call the appropriate stepping routine;
399 AliModule *det = dynamic_cast<AliModule*>(gAlice->Modules()->At(id));
400 if(det && det->StepManagerIsEnabled()) {
e569033b 401 if(AliLog::GetGlobalDebugLevel()>0) fMCQA->StepManager(id);
5d12ce38 402 det->StepManager();
403 }
404 }
405}
406
407//_______________________________________________________________________
408void AliMC::EnergySummary()
409{
9ac3aec9 410 //e
5d12ce38 411 // Print summary of deposited energy
412 //
413
414 Int_t ndep=0;
415 Float_t edtot=0;
416 Float_t ed, ed2;
417 Int_t kn, i, left, j, id;
418 const Float_t kzero=0;
419 Int_t ievent=gAlice->GetRunLoader()->GetHeader()->GetEvent()+1;
420 //
421 // Energy loss information
422 if(ievent) {
423 printf("***************** Energy Loss Information per event (GEV) *****************\n");
424 for(kn=1;kn<fEventEnergy.GetSize();kn++) {
425 ed=fSummEnergy[kn];
426 if(ed>0) {
427 fEventEnergy[ndep]=kn;
428 if(ievent>1) {
429 ed=ed/ievent;
430 ed2=fSum2Energy[kn];
431 ed2=ed2/ievent;
432 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,kzero))/ed;
433 } else
434 ed2=99;
435 fSummEnergy[ndep]=ed;
436 fSum2Energy[ndep]=TMath::Min(static_cast<Float_t>(99.),TMath::Max(ed2,kzero));
437 edtot+=ed;
438 ndep++;
439 }
440 }
441 for(kn=0;kn<(ndep-1)/3+1;kn++) {
442 left=ndep-kn*3;
443 for(i=0;i<(3<left?3:left);i++) {
444 j=kn*3+i;
445 id=Int_t (fEventEnergy[j]+0.1);
446 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),fSummEnergy[j],fSum2Energy[j]);
447 }
448 printf("\n");
449 }
450 //
451 // Relative energy loss in different detectors
452 printf("******************** Relative Energy Loss per event ********************\n");
453 printf("Total energy loss per event %10.3f GeV\n",edtot);
454 for(kn=0;kn<(ndep-1)/5+1;kn++) {
455 left=ndep-kn*5;
456 for(i=0;i<(5<left?5:left);i++) {
457 j=kn*5+i;
458 id=Int_t (fEventEnergy[j]+0.1);
459 printf(" %s %10.3f%%;",gMC->VolName(id),100*fSummEnergy[j]/edtot);
460 }
461 printf("\n");
462 }
463 for(kn=0;kn<75;kn++) printf("*");
464 printf("\n");
465 }
466 //
467 // Reset the TArray's
468 // fEventEnergy.Set(0);
469 // fSummEnergy.Set(0);
470 // fSum2Energy.Set(0);
471}
472
473//_____________________________________________________________________________
474void AliMC::BeginEvent()
475{
476 //
477 // Clean-up previous event
478 // Energy scores
21bf7095 479 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
480 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
481 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
482 AliDebug(1, " BEGINNING EVENT ");
483 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
484 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
485 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
5d12ce38 486
3cac0921 487 AliRunLoader *runloader=gAlice->GetRunLoader();
5d12ce38 488
489 /*******************************/
490 /* Clean after eventual */
491 /* previous event */
492 /*******************************/
493
494
495 //Set the next event in Run Loader -> Cleans trees (TreeK and all trees in detectors),
496 gAlice->SetEventNrInRun(gAlice->GetEventNrInRun()+1);
497 runloader->SetEventNumber(gAlice->GetEventNrInRun());// sets new files, cleans the previous event stuff, if necessary, etc.,
21bf7095 498 AliDebug(1, Form("EventNr is %d",gAlice->GetEventNrInRun()));
5d12ce38 499
500 fEventEnergy.Reset();
501 // Clean detector information
502
503 if (runloader->Stack())
3cac0921 504 runloader->Stack()->Reset();//clean stack -> tree is unloaded
5d12ce38 505 else
3cac0921 506 runloader->MakeStack();//or make a new one
504b172d 507
9ac3aec9 508
3cac0921 509 if(gAlice->Lego() == 0x0)
510 {
511 AliDebug(1, "fRunLoader->MakeTree(K)");
512 runloader->MakeTree("K");
513 }
514
21bf7095 515 AliDebug(1, "gMC->SetStack(fRunLoader->Stack())");
5d12ce38 516 gMC->SetStack(gAlice->GetRunLoader()->Stack());//Was in InitMC - but was moved here
517 //because we don't have guarantee that
518 //stack pointer is not going to change from event to event
519 //since it bellobgs to header and is obtained via RunLoader
520 //
521 // Reset all Detectors & kinematics & make/reset trees
522 //
523
8c8fbdcd 524 runloader->GetHeader()->Reset(AliCDBManager::Instance()->GetRun(),gAlice->GetEvNumber(),
3cac0921 525 gAlice->GetEventNrInRun());
5d12ce38 526// fRunLoader->WriteKinematics("OVERWRITE"); is there any reason to rewrite here since MakeTree does so
527
504b172d 528 if(gAlice->Lego())
3cac0921 529 {
530 gAlice->Lego()->BeginEvent();
531 return;
532 }
504b172d 533
5d12ce38 534
21bf7095 535 AliDebug(1, "ResetHits()");
5d12ce38 536 ResetHits();
504b172d 537
21bf7095 538 AliDebug(1, "fRunLoader->MakeTree(H)");
5d12ce38 539 runloader->MakeTree("H");
504b172d 540
5d12ce38 541
5d12ce38 542
3cac0921 543 MakeTmpTrackRefsTree();
5d12ce38 544 //create new branches and SetAdresses
545 TIter next(gAlice->Modules());
546 AliModule *detector;
547 while((detector = (AliModule*)next()))
548 {
3cac0921 549 AliDebug(2, Form("%s->MakeBranch(H)",detector->GetName()));
550 detector->MakeBranch("H");
5d12ce38 551 }
552}
553
554//_______________________________________________________________________
555void AliMC::ResetHits()
556{
557 //
558 // Reset all Detectors hits
559 //
560 TIter next(gAlice->Modules());
561 AliModule *detector;
562 while((detector = dynamic_cast<AliModule*>(next()))) {
563 detector->ResetHits();
564 }
565}
566
567//_______________________________________________________________________
568void AliMC::PostTrack()
569{
d0d4a6b3 570 // Posts tracks for each module
4a9de4af 571 TObjArray &dets = *gAlice->Modules();
572 AliModule *module;
573
574 for(Int_t i=0; i<=gAlice->GetNdets(); i++)
575 if((module = dynamic_cast<AliModule*>(dets[i])))
576 module->PostTrack();
5d12ce38 577}
578
579//_______________________________________________________________________
580void AliMC::FinishPrimary()
581{
582 //
583 // Called at the end of each primary track
584 //
585 AliRunLoader *runloader=gAlice->GetRunLoader();
586 // static Int_t count=0;
587 // const Int_t times=10;
588 // This primary is finished, purify stack
c75cfb51 589#if ROOT_VERSION_CODE > 262152
942a9039 590 if (!(gMC->SecondariesAreOrdered())) {
77db8cbf 591 if (runloader->Stack()->ReorderKine()) RemapHits();
942a9039 592 }
c75cfb51 593#endif
77db8cbf 594 if (runloader->Stack()->PurifyKine()) RemapHits();
504b172d 595
5d12ce38 596 TIter next(gAlice->Modules());
597 AliModule *detector;
598 while((detector = dynamic_cast<AliModule*>(next()))) {
599 detector->FinishPrimary();
504b172d 600 AliLoader* loader = detector->GetLoader();
601 if(loader)
602 {
603 TTree* treeH = loader->TreeH();
604 if (treeH) treeH->Fill(); //can be Lego run and treeH can not exist
605 }
5d12ce38 606 }
607
608 // Write out track references if any
3cac0921 609 if (fTmpTreeTR) fTmpTreeTR->Fill();
5d12ce38 610}
611
942a9039 612void AliMC::RemapHits()
613{
614//
615// Remaps the track labels of the hits
616 AliRunLoader *runloader=gAlice->GetRunLoader();
617 AliStack* stack = runloader->Stack();
618 TList* hitLists = GetHitLists();
619 TIter next(hitLists);
620 TCollection *hitList;
621
622 while((hitList = dynamic_cast<TCollection*>(next()))) {
623 TIter nexthit(hitList);
624 AliHit *hit;
625 while((hit = dynamic_cast<AliHit*>(nexthit()))) {
626 hit->SetTrack(stack->TrackLabel(hit->GetTrack()));
627 }
628 }
629
630 //
631 // This for detectors which have a special mapping mechanism
632 // for hits, such as TPC and TRD
633 //
ac2b85c0 634
942a9039 635
636 TObjArray* modules = gAlice->Modules();
637 TIter nextmod(modules);
ac2b85c0 638 AliModule *module;
639 while((module = (AliModule*) nextmod())) {
640 AliDetector* det = dynamic_cast<AliDetector*> (module);
641 if (det) det->RemapTrackHitIDs(stack->TrackLabelMap());
942a9039 642 }
643 //
644 RemapTrackReferencesIDs(stack->TrackLabelMap());
645}
646
5d12ce38 647//_______________________________________________________________________
648void AliMC::FinishEvent()
649{
650 //
651 // Called at the end of the event.
652 //
653
654 //
38ca2ad6 655
5d12ce38 656 if(gAlice->Lego()) gAlice->Lego()->FinishEvent();
657
658 TIter next(gAlice->Modules());
659 AliModule *detector;
660 while((detector = dynamic_cast<AliModule*>(next()))) {
661 detector->FinishEvent();
662 }
663
664 //Update the energy deposit tables
665 Int_t i;
666 for(i=0;i<fEventEnergy.GetSize();i++)
667 {
668 fSummEnergy[i]+=fEventEnergy[i];
669 fSum2Energy[i]+=fEventEnergy[i]*fEventEnergy[i];
670 }
671
672 AliRunLoader *runloader=gAlice->GetRunLoader();
673
674 AliHeader* header = runloader->GetHeader();
675 AliStack* stack = runloader->Stack();
676 if ( (header == 0x0) || (stack == 0x0) )
677 {//check if we got header and stack. If not cry and exit aliroot
21bf7095 678 AliFatal("Can not get the stack or header from LOADER");
5d12ce38 679 return;//never reached
680 }
681 // Update Header information
682 header->SetNprimary(stack->GetNprimary());
683 header->SetNtrack(stack->GetNtrack());
684
5d12ce38 685 // Write out the kinematics
38ca2ad6 686 if (!gAlice->Lego()) stack->FinishEvent();
3cac0921 687
688 // Synchronize the TreeTR with TreeK
172bc1fb 689 if (fTmpTreeTR) ReorderAndExpandTreeTR();
5d12ce38 690
691 // Write out the event Header information
692 TTree* treeE = runloader->TreeE();
693 if (treeE)
694 {
695 header->SetStack(stack);
696 treeE->Fill();
697 }
698 else
699 {
21bf7095 700 AliError("Can not get TreeE from RL");
5d12ce38 701 }
702
504b172d 703 if(gAlice->Lego() == 0x0)
704 {
705 runloader->WriteKinematics("OVERWRITE");
706 runloader->WriteTrackRefs("OVERWRITE");
707 runloader->WriteHits("OVERWRITE");
708 }
709
21bf7095 710 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
711 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
712 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
713 AliDebug(1, " FINISHING EVENT ");
714 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
715 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
716 AliDebug(1, "<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<");
5d12ce38 717}
718
719//_______________________________________________________________________
720void AliMC::Field(const Double_t* x, Double_t* b) const
721{
d0d4a6b3 722 // Calculates field "b" at point "x"
5d12ce38 723 gAlice->Field(x,b);
724}
725
726//_______________________________________________________________________
727void AliMC::Init()
728{
d0d4a6b3 729 // MC initialization
5d12ce38 730
731 //=================Create Materials and geometry
732 gMC->Init();
4e889ff4 733 // Set alignable volumes for the whole geometry (with old root)
734#if ROOT_VERSION_CODE < 331527
735 SetAllAlignableVolumes();
736#endif
5d12ce38 737 //Read the cuts for all materials
738 ReadTransPar();
739 //Build the special IMEDIA table
740 MediaTable();
741
742 //Compute cross-sections
743 gMC->BuildPhysics();
744
5d12ce38 745 //Initialise geometry deposition table
746 fEventEnergy.Set(gMC->NofVolumes()+1);
747 fSummEnergy.Set(gMC->NofVolumes()+1);
748 fSum2Energy.Set(gMC->NofVolumes()+1);
749
750 //
751 fMCQA = new AliMCQA(gAlice->GetNdets());
752
753 // Register MC in configuration
754 AliConfig::Instance()->Add(gMC);
755
756}
757
758//_______________________________________________________________________
759void AliMC::MediaTable()
760{
761 //
762 // Built media table to get from the media number to
763 // the detector id
764 //
765
766 Int_t kz, nz, idt, lz, i, k, ind;
767 // Int_t ibeg;
768 TObjArray &dets = *gAlice->Detectors();
769 AliModule *det;
770 Int_t ndets=gAlice->GetNdets();
771 //
772 // For all detectors
773 for (kz=0;kz<ndets;kz++) {
774 // If detector is defined
775 if((det=dynamic_cast<AliModule*>(dets[kz]))) {
776 TArrayI &idtmed = *(det->GetIdtmed());
777 for(nz=0;nz<100;nz++) {
0561efeb 778
5d12ce38 779 // Find max and min material number
780 if((idt=idtmed[nz])) {
781 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
782 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
783 }
784 }
785 if(det->LoMedium() > det->HiMedium()) {
786 det->LoMedium() = 0;
787 det->HiMedium() = 0;
788 } else {
789 if(det->HiMedium() > fImedia->GetSize()) {
21bf7095 790 AliError(Form("Increase fImedia from %d to %d",
791 fImedia->GetSize(),det->HiMedium()));
5d12ce38 792 return;
793 }
794 // Tag all materials in rage as belonging to detector kz
795 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
796 (*fImedia)[lz]=kz;
797 }
798 }
799 }
800 }
801 //
802 // Print summary table
21bf7095 803 AliInfo("Tracking media ranges:");
804 ToAliInfo(
5d12ce38 805 for(i=0;i<(ndets-1)/6+1;i++) {
806 for(k=0;k< (6<ndets-i*6?6:ndets-i*6);k++) {
807 ind=i*6+k;
808 det=dynamic_cast<AliModule*>(dets[ind]);
809 if(det)
810 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
811 det->HiMedium());
812 else
813 printf(" %6s: %3d -> %3d;","NULL",0,0);
814 }
815 printf("\n");
816 }
21bf7095 817 )
5d12ce38 818}
819
820//_______________________________________________________________________
821void AliMC::ReadTransPar()
822{
823 //
824 // Read filename to set the transport parameters
825 //
826
827
828 const Int_t kncuts=10;
829 const Int_t knflags=11;
830 const Int_t knpars=kncuts+knflags;
831 const char kpars[knpars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
832 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
833 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
834 "MULS","PAIR","PHOT","RAYL"};
835 char line[256];
836 char detName[7];
837 char* filtmp;
838 Float_t cut[kncuts];
839 Int_t flag[knflags];
840 Int_t i, itmed, iret, ktmed, kz;
841 FILE *lun;
842 //
843 // See whether the file is there
844 filtmp=gSystem->ExpandPathName(fTransParName.Data());
845 lun=fopen(filtmp,"r");
846 delete [] filtmp;
847 if(!lun) {
21bf7095 848 AliWarning(Form("File %s does not exist!",fTransParName.Data()));
5d12ce38 849 return;
850 }
851 //
5d12ce38 852 while(1) {
853 // Initialise cuts and flags
854 for(i=0;i<kncuts;i++) cut[i]=-99;
855 for(i=0;i<knflags;i++) flag[i]=-99;
856 itmed=0;
857 for(i=0;i<256;i++) line[i]='\0';
858 // Read up to the end of line excluded
859 iret=fscanf(lun,"%[^\n]",line);
860 if(iret<0) {
861 //End of file
862 fclose(lun);
5d12ce38 863 return;
864 }
865 // Read the end of line
866 fscanf(lun,"%*c");
867 if(!iret) continue;
868 if(line[0]=='*') continue;
869 // Read the numbers
870 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",
871 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
872 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
873 &flag[8],&flag[9],&flag[10]);
874 if(!iret) continue;
875 if(iret<0) {
876 //reading error
21bf7095 877 AliWarning(Form("Error reading file %s",fTransParName.Data()));
5d12ce38 878 continue;
879 }
880 // Check that the module exist
881 AliModule *mod = gAlice->GetModule(detName);
882 if(mod) {
883 // Get the array of media numbers
884 TArrayI &idtmed = *mod->GetIdtmed();
885 // Check that the tracking medium code is valid
886 if(0<=itmed && itmed < 100) {
887 ktmed=idtmed[itmed];
888 if(!ktmed) {
21bf7095 889 AliWarning(Form("Invalid tracking medium code %d for %s",itmed,mod->GetName()));
5d12ce38 890 continue;
891 }
892 // Set energy thresholds
893 for(kz=0;kz<kncuts;kz++) {
894 if(cut[kz]>=0) {
21bf7095 895 AliDebug(2, Form("%-6s set to %10.3E for tracking medium code %4d for %s",
896 kpars[kz],cut[kz],itmed,mod->GetName()));
5d12ce38 897 gMC->Gstpar(ktmed,kpars[kz],cut[kz]);
898 }
899 }
900 // Set transport mechanisms
901 for(kz=0;kz<knflags;kz++) {
902 if(flag[kz]>=0) {
21bf7095 903 AliDebug(2, Form("%-6s set to %10d for tracking medium code %4d for %s",
904 kpars[kncuts+kz],flag[kz],itmed,mod->GetName()));
5d12ce38 905 gMC->Gstpar(ktmed,kpars[kncuts+kz],Float_t(flag[kz]));
906 }
907 }
908 } else {
21bf7095 909 AliWarning(Form("Invalid medium code %d",itmed));
5d12ce38 910 continue;
911 }
912 } else {
21bf7095 913 AliDebug(1, Form("%s not present",detName));
5d12ce38 914 continue;
915 }
916 }
917}
918
919//_______________________________________________________________________
920void AliMC::SetTransPar(const char *filename)
921{
922 //
923 // Sets the file name for transport parameters
924 //
925 fTransParName = filename;
926}
927
928//_______________________________________________________________________
0054628d 929void AliMC::Browse(TBrowser *b)
5d12ce38 930{
931 //
932 // Called when the item "Run" is clicked on the left pane
933 // of the Root browser.
934 // It displays the Root Trees and all detectors.
935 //
936 //detectors are in folders anyway
937 b->Add(fMCQA,"AliMCQA");
938}
939
5d12ce38 940//_______________________________________________________________________
941void AliMC::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
942{
943 //
944 // Add a hit to detector id
945 //
946 TObjArray &dets = *gAlice->Modules();
947 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddHit(track,vol,hits);
948}
949
950//_______________________________________________________________________
951void AliMC::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
952{
953 //
954 // Add digit to detector id
955 //
956 TObjArray &dets = *gAlice->Modules();
957 if(dets[id]) dynamic_cast<AliModule*>(dets[id])->AddDigit(tracks,digits);
958}
959
960//_______________________________________________________________________
961Int_t AliMC::GetCurrentTrackNumber() const {
962 //
963 // Returns current track
964 //
965 return gAlice->GetRunLoader()->Stack()->GetCurrentTrackNumber();
966}
967
968//_______________________________________________________________________
969void AliMC::DumpPart (Int_t i) const
970{
971 //
972 // Dumps particle i in the stack
973 //
974 AliRunLoader * runloader = gAlice->GetRunLoader();
975 if (runloader->Stack())
976 runloader->Stack()->DumpPart(i);
977}
978
979//_______________________________________________________________________
980void AliMC::DumpPStack () const
981{
982 //
983 // Dumps the particle stack
984 //
985 AliRunLoader * runloader = gAlice->GetRunLoader();
986 if (runloader->Stack())
987 runloader->Stack()->DumpPStack();
988}
989
990//_______________________________________________________________________
991Int_t AliMC::GetNtrack() const {
992 //
993 // Returns number of tracks in stack
994 //
995 Int_t ntracks = -1;
996 AliRunLoader * runloader = gAlice->GetRunLoader();
997 if (runloader->Stack())
998 ntracks = runloader->Stack()->GetNtrack();
999 return ntracks;
1000}
1001
1002//_______________________________________________________________________
1003Int_t AliMC::GetPrimary(Int_t track) const
1004{
1005 //
1006 // return number of primary that has generated track
1007 //
1008 Int_t nprimary = -999;
1009 AliRunLoader * runloader = gAlice->GetRunLoader();
1010 if (runloader->Stack())
1011 nprimary = runloader->Stack()->GetPrimary(track);
1012 return nprimary;
1013}
1014
1015//_______________________________________________________________________
1016TParticle* AliMC::Particle(Int_t i) const
1017{
d0d4a6b3 1018 // Returns the i-th particle from the stack taking into account
1019 // the remaping done by PurifyKine
5d12ce38 1020 AliRunLoader * runloader = gAlice->GetRunLoader();
1021 if (runloader)
1022 if (runloader->Stack())
1023 return runloader->Stack()->Particle(i);
1024 return 0x0;
1025}
1026
1027//_______________________________________________________________________
1028TObjArray* AliMC::Particles() const {
1029 //
1030 // Returns pointer to Particles array
1031 //
1032 AliRunLoader * runloader = gAlice->GetRunLoader();
1033 if (runloader)
1034 if (runloader->Stack())
1035 return runloader->Stack()->Particles();
1036 return 0x0;
1037}
1038
1039//_______________________________________________________________________
1040void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
1041 Float_t *vpos, Float_t *polar, Float_t tof,
2057aecb 1042 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
5d12ce38 1043{
1044// Delegate to stack
1045//
1046 AliRunLoader * runloader = gAlice->GetRunLoader();
1047 if (runloader)
1048 if (runloader->Stack())
1049 runloader->Stack()->PushTrack(done, parent, pdg, pmom, vpos, polar, tof,
1050 mech, ntr, weight, is);
1051}
1052
1053//_______________________________________________________________________
1054void AliMC::PushTrack(Int_t done, Int_t parent, Int_t pdg,
1055 Double_t px, Double_t py, Double_t pz, Double_t e,
1056 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
1057 Double_t polx, Double_t poly, Double_t polz,
2057aecb 1058 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is) const
5d12ce38 1059{
1060 // Delegate to stack
1061 //
1062 AliRunLoader * runloader = gAlice->GetRunLoader();
1063 if (runloader)
1064 if (runloader->Stack())
1065 runloader->Stack()->PushTrack(done, parent, pdg, px, py, pz, e, vx, vy, vz, tof,
1066 polx, poly, polz, mech, ntr, weight, is);
1067}
1068
1069//_______________________________________________________________________
2057aecb 1070void AliMC::SetHighWaterMark(Int_t nt) const
5d12ce38 1071{
1072 //
1073 // Set high water mark for last track in event
1074 AliRunLoader * runloader = gAlice->GetRunLoader();
1075 if (runloader)
1076 if (runloader->Stack())
1077 runloader->Stack()->SetHighWaterMark(nt);
1078}
1079
1080//_______________________________________________________________________
2057aecb 1081void AliMC::KeepTrack(Int_t track) const
5d12ce38 1082{
1083 //
1084 // Delegate to stack
1085 //
1086 AliRunLoader * runloader = gAlice->GetRunLoader();
1087 if (runloader)
1088 if (runloader->Stack())
1089 runloader->Stack()->KeepTrack(track);
1090}
1091
1092//_______________________________________________________________________
2057aecb 1093void AliMC::FlagTrack(Int_t track) const
5d12ce38 1094{
1095 // Delegate to stack
1096 //
1097 AliRunLoader * runloader = gAlice->GetRunLoader();
1098 if (runloader)
1099 if (runloader->Stack())
1100 runloader->Stack()->FlagTrack(track);
1101}
1102
1103//_______________________________________________________________________
2057aecb 1104void AliMC::SetCurrentTrack(Int_t track) const
5d12ce38 1105{
1106 //
1107 // Set current track number
1108 //
1109 AliRunLoader * runloader = gAlice->GetRunLoader();
1110 if (runloader)
1111 if (runloader->Stack())
1112 runloader->Stack()->SetCurrentTrack(track);
1113}
1114
1115//_______________________________________________________________________
9ac3aec9 1116AliTrackReference* AliMC::AddTrackReference(Int_t label, Int_t id)
1117{
5d12ce38 1118 //
1119 // add a trackrefernce to the list
1120 if (!fTrackReferences) {
21bf7095 1121 AliError("Container trackrefernce not active");
9ac3aec9 1122 return NULL;
5d12ce38 1123 }
9ac3aec9 1124 Int_t primary = GetPrimary(label);
1125 Particle(primary)->SetBit(kKeepBit);
1126
3cac0921 1127 Int_t nref = fTmpTrackReferences->GetEntriesFast();
1128 TClonesArray &lref = *fTmpTrackReferences;
9ac3aec9 1129 return new(lref[nref]) AliTrackReference(label, id);
5d12ce38 1130}
1131
1132
1133
1134//_______________________________________________________________________
1135void AliMC::ResetTrackReferences()
1136{
1137 //
1138 // Reset all references
1139 //
3cac0921 1140 if (fTmpTrackReferences) fTmpTrackReferences->Clear();
5d12ce38 1141}
1142
1143void AliMC::RemapTrackReferencesIDs(Int_t *map)
1144{
1145 //
1146 // Remapping track reference
1147 // Called at finish primary
1148 //
3cac0921 1149 if (!fTmpTrackReferences) return;
1150 Int_t nEntries = fTmpTrackReferences->GetEntries();
4e490a96 1151 for (Int_t i=0; i < nEntries; i++){
3cac0921 1152 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTmpTrackReferences->UncheckedAt(i));
4e490a96 1153 if (ref) {
1154 Int_t newID = map[ref->GetTrack()];
1155 if (newID>=0) ref->SetTrack(newID);
1156 else {
1157 ref->SetBit(kNotDeleted,kFALSE);
3cac0921 1158 fTmpTrackReferences->RemoveAt(i);
4e490a96 1159 }
1160 } // if ref
5d12ce38 1161 }
3cac0921 1162 fTmpTrackReferences->Compress();
5d12ce38 1163}
0561efeb 1164
1165void AliMC::FixParticleDecaytime()
1166{
1167 //
1168 // Fix the particle decay time according to rmin and rmax for decays
1169 //
1170
1171 TLorentzVector p;
1172 gMC->TrackMomentum(p);
1173 Double_t tmin, tmax;
1174 Double_t b;
1175
1176 // Transverse velocity
1177 Double_t vt = p.Pt() / p.E();
1178
1179 if ((b = gAlice->Field()->SolenoidField()) > 0.) { // [kG]
1180
1181 // Radius of helix
1182
1183 Double_t rho = p.Pt() / 0.0003 / b; // [cm]
1184
1185 // Revolution frequency
1186
1187 Double_t omega = vt / rho;
1188
1189 // Maximum and minimum decay time
1190 //
1191 // Check for curlers first
1192 if (fRDecayMax * fRDecayMax / rho / rho / 2. > 1.) return;
1193
1194 //
1195
1196 tmax = TMath::ACos(1. - fRDecayMax * fRDecayMax / rho / rho / 2.) / omega; // [ct]
1197 tmin = TMath::ACos(1. - fRDecayMin * fRDecayMin / rho / rho / 2.) / omega; // [ct]
1198 } else {
1199 tmax = fRDecayMax / vt; // [ct]
1200 tmin = fRDecayMin / vt; // [ct]
1201 }
1202
1203 //
1204 // Dial t using the two limits
1205 Double_t t = tmin + (tmax - tmin) * gRandom->Rndm(); // [ct]
1206 //
1207 //
1208 // Force decay time in transport code
1209 //
fa7ca4da 1210 gMC->ForceDecayTime(t / 2.99792458e10);
0561efeb 1211}
3cac0921 1212
1213void AliMC::MakeTmpTrackRefsTree()
1214{
1215 // Make the temporary track reference tree
1216 fTmpFileTR = new TFile("TrackRefsTmp.root", "recreate");
1217 fTmpTreeTR = new TTree("TreeTR", "Track References");
1218 if (!fTmpTrackReferences) fTmpTrackReferences = new TClonesArray("AliTrackReference", 100);
1219 fTmpTreeTR->Branch("TrackReferences", "TClonesArray", &fTmpTrackReferences, 4000);
1220}
1221
1222void AliMC::ReorderAndExpandTreeTR()
1223{
1224//
1225// Reorder and expand the temporary track reference tree in order to match the kinematics tree
1226//
1227
1228 AliRunLoader *rl = gAlice->GetRunLoader();
1229//
1230// TreeTR
1231 AliDebug(1, "fRunLoader->MakeTrackRefsContainer()");
1232 rl->MakeTrackRefsContainer();
1233 TTree * treeTR = rl->TreeTR();
1234 if (treeTR){
1235 // make branch for central track references
1236 if (!fTrackReferences) fTrackReferences = new TClonesArray("AliTrackReference",0);
1237 TBranch *branch;
1238 branch = treeTR->Branch("TrackReferences",&fTrackReferences);
1239 branch->SetAddress(&fTrackReferences);
1240 }
1241
1242 AliStack* stack = rl->Stack();
1243 Int_t np = stack->GetNprimary();
1244 Int_t nt = fTmpTreeTR->GetEntries();
3cac0921 1245 //
1246 // Loop over tracks and find the secondaries with the help of the kine tree
5ce21fc4 1247 Int_t ifills = 0;
1248 Int_t it = 0;
1249 for (Int_t ip = np - 1; ip > -1; ip--) {
1250 TParticle *part = stack->Particle(ip);
5a9cf418 1251 //printf("Particle %5d %5d %5d %5d %5d \n", ip, part->GetPdgCode(), part->GetFirstMother(), part->GetFirstDaughter(), part->GetLastDaughter());
2c4e7c77 1252
5ce21fc4 1253 // Skip primaries that have not been transported
3cac0921 1254 Int_t dau1 = part->GetFirstDaughter();
1255 Int_t dau2 = -1;
f4414561 1256 if (!part->TestBit(kTransportBit)) continue;
5ce21fc4 1257 //
1258 fTmpTreeTR->GetEntry(it++);
1259 Int_t nh = fTmpTrackReferences->GetEntries();
1260 // Determine range of secondaries produced by this primary
3cac0921 1261 if (dau1 > -1) {
5ce21fc4 1262 Int_t inext = ip - 1;
3cac0921 1263 while (dau2 < 0) {
1264 if (inext >= 0) {
1265 part = stack->Particle(inext);
5ce21fc4 1266 dau2 = part->GetFirstDaughter();
5a9cf418 1267 if (!(part->TestBit(kTransportBit)) || dau2 == -1 || dau2 < np) {
1268// if (dau2 == -1 || dau2 < np) {
5ce21fc4 1269 dau2 = -1;
1270 } else {
1271 dau2--;
1272 }
3cac0921 1273 } else {
1274 dau2 = stack->GetNtrack() - 1;
1275 }
1276 inext--;
5ce21fc4 1277 } // find upper bound
1278 } // dau2 < 0
2c4e7c77 1279// printf("Check (1) %5d %5d %5d %5d %5d \n", ip, np, it, dau1, dau2);
3cac0921 1280 //
1281 // Loop over reference hits and find secondary label
1282 for (Int_t id = dau1; (id <= dau2) && (dau1 > -1); id++) {
1283 for (Int_t ih = 0; ih < nh; ih++) {
1284 AliTrackReference* tr = (AliTrackReference*) fTmpTrackReferences->At(ih);
1285 Int_t label = tr->Label();
5ce21fc4 1286 // Skip primaries
1287 if (label == ip) continue;
2c4e7c77 1288 if (label > dau2 || label < dau1)
1289 AliWarning(Form("Track Reference Label out of range !: %5d %5d %5d \n", label, dau1, dau2));
3cac0921 1290 if (label == id) {
1291 // secondary found
1292 Int_t nref = fTrackReferences->GetEntriesFast();
1293 TClonesArray &lref = *fTrackReferences;
1294 new(lref[nref]) AliTrackReference(*tr);
1295 }
1296 } // hits
1297 treeTR->Fill();
1298 fTrackReferences->Clear();
5ce21fc4 1299 ifills++;
3cac0921 1300 } // daughters
1301 } // tracks
3cac0921 1302 //
1303 // Now loop again and write the primaries
5ce21fc4 1304 it = nt - 1;
1305 for (Int_t ip = 0; ip < np; ip++) {
1306 TParticle* part = stack->Particle(ip);
f4414561 1307// if ((part->GetFirstDaughter() == -1 && part->GetStatusCode() <= 1) || part->GetFirstDaughter() >= np)
c226988d 1308 if (part->TestBit(kTransportBit))
f4414561 1309 {
5ce21fc4 1310 // Skip particles that have not been transported
1311 fTmpTreeTR->GetEntry(it--);
1312 Int_t nh = fTmpTrackReferences->GetEntries();
1313 //
1314 // Loop over reference hits and find primary labels
1315 for (Int_t ih = 0; ih < nh; ih++) {
1316 AliTrackReference* tr = (AliTrackReference*) fTmpTrackReferences->At(ih);
1317 Int_t label = tr->Label();
1318 if (label == ip) {
1319 Int_t nref = fTrackReferences->GetEntriesFast();
1320 TClonesArray &lref = *fTrackReferences;
1321 new(lref[nref]) AliTrackReference(*tr);
1322 }
1323 }
3cac0921 1324 }
1325 treeTR->Fill();
1326 fTrackReferences->Clear();
5ce21fc4 1327 ifills++;
3cac0921 1328 } // tracks
5ce21fc4 1329 // Check
1330 if (ifills != stack->GetNtrack())
1331 AliWarning(Form("Number of entries in TreeTR (%5d) unequal to TreeK (%5d) \n", ifills, stack->GetNtrack()));
3cac0921 1332//
1333// Clean-up
1334 delete fTmpTreeTR;
1335 fTmpFileTR->Close();
1336 delete fTmpFileTR;
1337 delete fTmpTrackReferences;
1338 fTmpTrackReferences = 0;
1339 gSystem->Exec("rm -rf TrackRefsTmp.root");
1340}
1341