]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliMCEvent.cxx
cb81de01e5376e738293d0d8b6aebb73972d4b3c
[u/mrichter/AliRoot.git] / STEER / AliMCEvent.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, 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
18 //-------------------------------------------------------------------------
19 //     Class for Kinematic Events
20 //     Author: Andreas Morsch, CERN
21 //-------------------------------------------------------------------------
22 #include <TArrow.h>
23 #include <TMarker.h>
24 #include <TH2F.h>
25 #include <TTree.h>
26 #include <TFile.h>
27 #include <TParticle.h>
28 #include <TClonesArray.h>
29 #include <TRefArray.h>
30 #include <TList.h>
31 #include <TArrayF.h>
32
33 #include "AliLog.h"
34 #include "AliMCEvent.h"
35 #include "AliMCVertex.h"
36 #include "AliStack.h"
37 #include "AliTrackReference.h"
38 #include "AliHeader.h"
39 #include "AliGenEventHeader.h"
40
41
42 Int_t AliMCEvent::fgkBgLabelOffset(10000000);
43
44
45 AliMCEvent::AliMCEvent():
46     AliVEvent(),
47     fStack(0),
48     fMCParticles(0),
49     fMCParticleMap(0),
50     fHeader(new AliHeader()),
51     fTRBuffer(0),
52     fTrackReferences(new TClonesArray("AliTrackReference", 1000)),
53     fTreeTR(0),
54     fTmpTreeTR(0),
55     fTmpFileTR(0),
56     fNprimaries(-1),
57     fNparticles(-1),
58     fSubsidiaryEvents(0),
59     fPrimaryOffset(0),
60     fSecondaryOffset(0),
61     fExternal(0),
62     fVertex(0)
63 {
64     // Default constructor
65 }
66
67 AliMCEvent::AliMCEvent(const AliMCEvent& mcEvnt) :
68     AliVEvent(mcEvnt),
69     fStack(mcEvnt.fStack),
70     fMCParticles(mcEvnt.fMCParticles),
71     fMCParticleMap(mcEvnt.fMCParticleMap),
72     fHeader(mcEvnt.fHeader),
73     fTRBuffer(mcEvnt.fTRBuffer),
74     fTrackReferences(mcEvnt.fTrackReferences),
75     fTreeTR(mcEvnt.fTreeTR),
76     fTmpTreeTR(mcEvnt.fTmpTreeTR),
77     fTmpFileTR(mcEvnt.fTmpFileTR),
78     fNprimaries(mcEvnt.fNprimaries),
79     fNparticles(mcEvnt.fNparticles),
80     fSubsidiaryEvents(0),
81     fPrimaryOffset(0),
82     fSecondaryOffset(0),
83     fExternal(0),
84     fVertex(mcEvnt.fVertex)
85
86 // Copy constructor
87 }
88
89
90 AliMCEvent& AliMCEvent::operator=(const AliMCEvent& mcEvnt)
91 {
92     // assignment operator
93     if (this!=&mcEvnt) { 
94         AliVEvent::operator=(mcEvnt); 
95     }
96   
97     return *this; 
98 }
99
100 void AliMCEvent::ConnectTreeE (TTree* tree)
101 {
102     // Connect the event header tree
103     tree->SetBranchAddress("Header", &fHeader);
104 }
105
106 void AliMCEvent::ConnectTreeK (TTree* tree)
107 {
108     // Connect the kinematics tree to the stack
109     if (!fMCParticles) fMCParticles = new TClonesArray("AliMCParticle",1000);
110     //
111     fStack = fHeader->Stack();
112     fStack->ConnectTree(tree);
113     //
114     // Load the event
115     fStack->GetEvent();
116     fNparticles = fStack->GetNtrack();
117     fNprimaries = fStack->GetNprimary();
118     Int_t iev  = fHeader->GetEvent();
119     Int_t ievr = fHeader->GetEventNrInRun();
120     
121     AliInfo(Form("AliMCEvent# %5d %5d: Number of particles: %5d (all) %5d (primaries)\n", 
122                  iev, ievr, fNparticles, fNprimaries));
123  
124     // This is a cache for the TParticles converted to MCParticles on user request
125     if (fMCParticleMap) {
126         fMCParticleMap->Clear();
127         fMCParticles->Delete();
128         if (fNparticles>0) fMCParticleMap->Expand(fNparticles);
129     }
130     else
131         fMCParticleMap = new TRefArray(fNparticles);
132 }
133
134 void AliMCEvent::ConnectTreeTR (TTree* tree)
135 {
136     // Connect the track reference tree
137     fTreeTR = tree;
138     
139     if (fTreeTR->GetBranch("AliRun")) {
140         if (fTmpFileTR) {
141             fTmpFileTR->Close();
142             delete fTmpFileTR;
143         }
144         // This is an old format with one branch per detector not in synch with TreeK
145         ReorderAndExpandTreeTR();
146     } else {
147         // New format 
148         fTreeTR->SetBranchAddress("TrackReferences", &fTRBuffer);
149     }
150 }
151
152 Int_t AliMCEvent::GetParticleAndTR(Int_t i, TParticle*& particle, TClonesArray*& trefs)
153 {
154     // Retrieve entry i
155     if (i < 0 || i >= fNparticles) {
156         AliWarning(Form("AliMCEventHandler::GetEntry: Index out of range"));
157         particle = 0;
158         trefs    = 0;
159         return (-1);
160     }
161     particle = fStack->Particle(i);
162     if (fTreeTR) {
163         fTreeTR->GetEntry(fStack->TreeKEntry(i));
164         trefs    = fTRBuffer;
165         return trefs->GetEntries();
166     } else {
167         trefs = 0;
168         return -1;
169     }
170 }
171
172
173 void AliMCEvent::Clean()
174 {
175     // Clean-up before new trees are connected
176     delete fStack; fStack = 0;
177
178     // Clear TR
179     if (fTRBuffer) {
180         fTRBuffer->Delete();
181         delete fTRBuffer;
182         fTRBuffer = 0;
183     }
184     if (fVertex) {
185         delete fVertex;
186         fVertex = 0;
187     }
188     
189 }
190
191 void AliMCEvent::FinishEvent()
192 {
193   // Clean-up after event
194   //    
195     fStack->Reset(0);
196     fMCParticles->Delete();
197     fMCParticleMap->Clear();
198     if (fTRBuffer)
199       fTRBuffer->Delete();
200     fTrackReferences->Delete();
201     fNparticles = -1;
202     fNprimaries = -1;    
203     fStack      =  0;
204 //    fSubsidiaryEvents->Clear();
205     fSubsidiaryEvents = 0;
206 }
207
208
209
210 void AliMCEvent::DrawCheck(Int_t i, Int_t search)
211 {
212     //
213     // Simple event display for debugging
214     if (!fTreeTR) {
215         AliWarning("No Track Reference information available");
216         return;
217     } 
218     
219     if (i > -1 && i < fNparticles) {
220         fTreeTR->GetEntry(fStack->TreeKEntry(i));
221     } else {
222         AliWarning("AliMCEvent::GetEntry: Index out of range");
223     }
224     
225     Int_t nh = fTRBuffer->GetEntries();
226     
227     
228     if (search) {
229         while(nh <= search && i < fNparticles - 1) {
230             i++;
231             fTreeTR->GetEntry(fStack->TreeKEntry(i));
232             nh =  fTRBuffer->GetEntries();
233         }
234         printf("Found Hits at %5d\n", i);
235     }
236     TParticle* particle = fStack->Particle(i);
237     
238     TH2F*    h = new TH2F("", "", 100, -500, 500, 100, -500, 500);
239     Float_t x0 = particle->Vx();
240     Float_t y0 = particle->Vy();
241
242     Float_t x1 = particle->Vx() + particle->Px() * 50.;
243     Float_t y1 = particle->Vy() + particle->Py() * 50.;
244     
245     TArrow*  a = new TArrow(x0, y0, x1, y1, 0.01);
246     h->Draw();
247     a->SetLineColor(2);
248     
249     a->Draw();
250     
251     for (Int_t ih = 0; ih < nh; ih++) {
252         AliTrackReference* ref = (AliTrackReference*) fTRBuffer->At(ih);
253         TMarker* m = new TMarker(ref->X(), ref->Y(), 20);
254         m->Draw();
255         m->SetMarkerSize(0.4);
256         
257     }
258 }
259
260
261 void AliMCEvent::ReorderAndExpandTreeTR()
262 {
263 //
264 //  Reorder and expand the track reference tree in order to match the kinematics tree.
265 //  Copy the information from different branches into one
266 //
267 //  TreeTR
268
269     fTmpFileTR = new TFile("TrackRefsTmp.root", "recreate");
270     fTmpTreeTR = new TTree("TreeTR", "TrackReferences");
271     if (!fTRBuffer)  fTRBuffer = new TClonesArray("AliTrackReference", 100);
272     fTmpTreeTR->Branch("TrackReferences", "TClonesArray", &fTRBuffer, 32000, 0);
273     
274
275 //
276 //  Activate the used branches only. Otherwisw we get a bad memory leak.
277     fTreeTR->SetBranchStatus("*",        0);
278     fTreeTR->SetBranchStatus("AliRun.*", 1);
279     fTreeTR->SetBranchStatus("ITS.*",    1);
280     fTreeTR->SetBranchStatus("TPC.*",    1);
281     fTreeTR->SetBranchStatus("TRD.*",    1);
282     fTreeTR->SetBranchStatus("TOF.*",    1);
283     fTreeTR->SetBranchStatus("FRAME.*",  1);
284     fTreeTR->SetBranchStatus("MUON.*",   1);
285 //
286 //  Connect the active branches
287     TClonesArray* trefs[7];
288     for (Int_t i = 0; i < 7; i++) trefs[i] = 0;
289     if (fTreeTR){
290         // make branch for central track references
291         if (fTreeTR->GetBranch("AliRun")) fTreeTR->SetBranchAddress("AliRun", &trefs[0]);
292         if (fTreeTR->GetBranch("ITS"))    fTreeTR->SetBranchAddress("ITS",    &trefs[1]);
293         if (fTreeTR->GetBranch("TPC"))    fTreeTR->SetBranchAddress("TPC",    &trefs[2]);
294         if (fTreeTR->GetBranch("TRD"))    fTreeTR->SetBranchAddress("TRD",    &trefs[3]);
295         if (fTreeTR->GetBranch("TOF"))    fTreeTR->SetBranchAddress("TOF",    &trefs[4]);
296         if (fTreeTR->GetBranch("FRAME"))  fTreeTR->SetBranchAddress("FRAME",  &trefs[5]);
297         if (fTreeTR->GetBranch("MUON"))   fTreeTR->SetBranchAddress("MUON",   &trefs[6]);
298     }
299
300     Int_t np = fStack->GetNprimary();
301     Int_t nt = fTreeTR->GetEntries();
302     
303     //
304     // Loop over tracks and find the secondaries with the help of the kine tree
305     Int_t ifills = 0;
306     Int_t it     = 0;
307     Int_t itlast = 0;
308     TParticle* part;
309
310     for (Int_t ip = np - 1; ip > -1; ip--) {
311         part = fStack->Particle(ip);
312 //      printf("Particle %5d %5d %5d %5d %5d %5d \n", 
313 //             ip, part->GetPdgCode(), part->GetFirstMother(), part->GetFirstDaughter(), 
314 //             part->GetLastDaughter(), part->TestBit(kTransportBit));
315
316         // Determine range of secondaries produced by this primary during transport     
317         Int_t dau1  = part->GetFirstDaughter();
318         if (dau1 < np) continue;  // This particle has no secondaries produced during transport
319         Int_t dau2  = -1;
320         if (dau1 > -1) {
321             Int_t inext = ip - 1;
322             while (dau2 < 0) {
323                 if (inext >= 0) {
324                     part = fStack->Particle(inext);
325                     dau2 =  part->GetFirstDaughter();
326                     if (dau2 == -1 || dau2 < np) {
327                         dau2 = -1;
328                     } else {
329                         dau2--;
330                     }
331                 } else {
332                     dau2 = fStack->GetNtrack() - 1;
333                 }
334                 inext--;
335             } // find upper bound
336         }  // dau2 < 0
337         
338
339 //      printf("Check (1) %5d %5d %5d %5d %5d \n", ip, np, it, dau1, dau2);
340 //
341 // Loop over reference hits and find secondary label
342 // First the tricky part: find the entry in treeTR than contains the hits or
343 // make sure that no hits exist.
344 //
345         Bool_t hasHits   = kFALSE;
346         Bool_t isOutside = kFALSE;
347
348         it = itlast;
349         while (!hasHits && !isOutside && it < nt) {
350             fTreeTR->GetEntry(it++);
351             for (Int_t ib = 0; ib < 7; ib++) {
352                 if (!trefs[ib]) continue;
353                 Int_t nh = trefs[ib]->GetEntries();
354                 for (Int_t ih = 0; ih < nh; ih++) {
355                     AliTrackReference* tr = (AliTrackReference*) trefs[ib]->At(ih);
356                     Int_t label = tr->Label();
357                     if (label >= dau1 && label <= dau2) {
358                         hasHits = kTRUE;
359                         itlast = it - 1;
360                         break;
361                     }
362                     if (label > dau2 || label < ip) {
363                         isOutside = kTRUE;
364                         itlast = it - 1;
365                         break;
366                     }
367                 } // hits
368                 if (hasHits || isOutside) break;
369             } // branches
370         } // entries
371
372         if (!hasHits) {
373             // Write empty entries
374             for (Int_t id = dau1; (id <= dau2); id++) {
375                 fTmpTreeTR->Fill();
376                 ifills++;
377             } 
378         } else {
379             // Collect all hits
380             fTreeTR->GetEntry(itlast);
381             for (Int_t id = dau1; (id <= dau2) && (dau1 > -1); id++) {
382                 for (Int_t ib = 0; ib < 7; ib++) {
383                     if (!trefs[ib]) continue;
384                     Int_t nh = trefs[ib]->GetEntries();
385                     for (Int_t ih = 0; ih < nh; ih++) {
386                         AliTrackReference* tr = (AliTrackReference*) trefs[ib]->At(ih);
387                         Int_t label = tr->Label();
388                         // Skip primaries
389                         if (label == ip) continue;
390                         if (label > dau2 || label < dau1) 
391                             printf("AliMCEventHandler::Track Reference Label out of range !: %5d %5d %5d %5d \n", 
392                                    itlast, label, dau1, dau2);
393                         if (label == id) {
394                             // secondary found
395                             tr->SetDetectorId(ib-1);
396                             Int_t nref =  fTRBuffer->GetEntriesFast();
397                             TClonesArray &lref = *fTRBuffer;
398                             new(lref[nref]) AliTrackReference(*tr);
399                         }
400                     } // hits
401                 } // branches
402                 fTmpTreeTR->Fill();
403                 fTRBuffer->Delete();
404                 ifills++;
405             } // daughters
406         } // has hits
407     } // tracks
408
409     //
410     // Now loop again and write the primaries
411     //
412     it = nt - 1;
413     for (Int_t ip = 0; ip < np; ip++) {
414         Int_t labmax = -1;
415         while (labmax < ip && it > -1) {
416             fTreeTR->GetEntry(it--);
417             for (Int_t ib = 0; ib < 7; ib++) {
418                 if (!trefs[ib]) continue;
419                 Int_t nh = trefs[ib]->GetEntries();
420                 // 
421                 // Loop over reference hits and find primary labels
422                 for (Int_t ih = 0; ih < nh; ih++) {
423                     AliTrackReference* tr = (AliTrackReference*)  trefs[ib]->At(ih);
424                     Int_t label = tr->Label();
425                     if (label < np && label > labmax) {
426                         labmax = label;
427                     }
428                     
429                     if (label == ip) {
430                         tr->SetDetectorId(ib-1);
431                         Int_t nref = fTRBuffer->GetEntriesFast();
432                         TClonesArray &lref = *fTRBuffer;
433                         new(lref[nref]) AliTrackReference(*tr);
434                     }
435                 } // hits
436             } // branches
437         } // entries
438         it++;
439         fTmpTreeTR->Fill();
440         fTRBuffer->Delete();
441         ifills++;
442     } // tracks
443     // Check
444
445
446     // Clean-up
447     delete fTreeTR; fTreeTR = 0;
448     
449     for (Int_t ib = 0; ib < 7; ib++) {
450         if (trefs[ib]) {
451             trefs[ib]->Clear();
452             delete trefs[ib];
453             trefs[ib] = 0;
454         }
455     }
456
457     if (ifills != fStack->GetNtrack()) 
458         printf("AliMCEvent:Number of entries in TreeTR (%5d) unequal to TreeK (%5d) \n", 
459                ifills, fStack->GetNtrack());
460
461     fTmpTreeTR->Write();
462     fTreeTR = fTmpTreeTR;
463 }
464
465 AliVParticle* AliMCEvent::GetTrack(Int_t i) const
466 {
467     // Get MC Particle i
468     //
469
470     if (fExternal) {
471         return ((AliVParticle*) (fMCParticles->At(i)));
472     }
473     
474     //
475     // Check first if this explicitely accesses the subsidiary event
476     
477     if (i >= BgLabelOffset()) {
478         if (fSubsidiaryEvents) {
479             AliMCEvent* bgEvent = (AliMCEvent*) (fSubsidiaryEvents->At(1));
480             return (bgEvent->GetTrack(i - BgLabelOffset()));
481         } else {
482             return 0;
483         }
484     }
485     
486     //
487     AliMCParticle *mcParticle = 0;
488     TParticle     *particle   = 0;
489     TClonesArray  *trefs      = 0;
490     Int_t          ntref      = 0;
491     TRefArray     *rarray     = 0;
492
493
494
495     // Out of range check
496     if (i < 0 || i >= fNparticles) {
497         AliWarning(Form("AliMCEvent::GetEntry: Index out of range"));
498         mcParticle = 0;
499         return (mcParticle);
500     }
501
502     
503     if (fSubsidiaryEvents) {
504         AliMCEvent*   mc;
505         Int_t idx = FindIndexAndEvent(i, mc);
506         return (mc->GetTrack(idx));
507     } 
508
509     //
510     // First check If the MC Particle has been already cached
511     if(!fMCParticleMap->At(i)) {
512         // Get particle from the stack
513         particle   = fStack->Particle(i);
514         // Get track references from Tree TR
515         if (fTreeTR) {
516             fTreeTR->GetEntry(fStack->TreeKEntry(i));
517             trefs     = fTRBuffer;
518             ntref     = trefs->GetEntriesFast();
519             rarray    = new TRefArray(ntref);
520             Int_t nen = fTrackReferences->GetEntriesFast();
521             for (Int_t j = 0; j < ntref; j++) {
522                 // Save the track references in a TClonesArray
523                 AliTrackReference* ref = dynamic_cast<AliTrackReference*>((*fTRBuffer)[j]);
524                 // Save the pointer in a TRefArray
525                 new ((*fTrackReferences)[nen]) AliTrackReference(*ref);
526                 rarray->AddAt((*fTrackReferences)[nen], j);
527                 nen++;
528             } // loop over track references for entry i
529         } // if TreeTR available
530         Int_t nentries = fMCParticles->GetEntriesFast();
531         new ((*fMCParticles)[nentries]) AliMCParticle(particle, rarray, i);
532         mcParticle = dynamic_cast<AliMCParticle*>((*fMCParticles)[nentries]);
533         fMCParticleMap->AddAt(mcParticle, i);
534
535         TParticle* part = mcParticle->Particle();
536         Int_t imo  = part->GetFirstMother();
537         Int_t id1  = part->GetFirstDaughter();
538         Int_t id2  = part->GetLastDaughter();
539         if (fPrimaryOffset > 0 || fSecondaryOffset > 0) {
540             // Remapping of the mother and daughter indices
541             if (imo < fNprimaries) {
542                 mcParticle->SetMother(imo + fPrimaryOffset);
543             } else {
544                 mcParticle->SetMother(imo + fSecondaryOffset - fNprimaries);
545             }
546
547             if (id1 < fNprimaries) {
548                 mcParticle->SetFirstDaughter(id1 + fPrimaryOffset);
549                 mcParticle->SetLastDaughter (id2 + fPrimaryOffset);
550             } else {
551                 mcParticle->SetFirstDaughter(id1 + fSecondaryOffset - fNprimaries);
552                 mcParticle->SetLastDaughter (id2 + fSecondaryOffset - fNprimaries);
553             }
554
555             
556             if (i > fNprimaries) {
557                 mcParticle->SetLabel(i + fPrimaryOffset);
558             } else {
559                 mcParticle->SetLabel(i + fSecondaryOffset - fNprimaries);
560             }
561             
562             
563         } else {
564             mcParticle->SetFirstDaughter(id1);
565             mcParticle->SetLastDaughter (id2);
566             mcParticle->SetMother       (imo);
567         }
568
569     } else {
570         mcParticle = dynamic_cast<AliMCParticle*>(fMCParticleMap->At(i));
571     }
572     
573     
574     return mcParticle;
575 }
576
577 AliGenEventHeader* AliMCEvent::GenEventHeader() const {return (fHeader->GenEventHeader());}
578
579
580 void AliMCEvent::AddSubsidiaryEvent(AliMCEvent* event) 
581 {
582     // Add a subsidiary event to the list; for example merged background event.
583     if (!fSubsidiaryEvents) {
584         TList* events = new TList();
585         events->Add(new AliMCEvent(*this));
586         fSubsidiaryEvents = events;
587     }
588     
589     fSubsidiaryEvents->Add(event);
590 }
591
592 Int_t AliMCEvent::FindIndexAndEvent(Int_t oldidx, AliMCEvent*& event) const
593 {
594     // Find the index and event in case of composed events like signal + background
595     TIter next(fSubsidiaryEvents);
596     next.Reset();
597      if (oldidx < fNprimaries) {
598         while((event = (AliMCEvent*)next())) {
599             if (oldidx < (event->GetPrimaryOffset() + event->GetNumberOfPrimaries())) break;
600         }
601         return (oldidx - event->GetPrimaryOffset());
602     } else {
603         while((event = (AliMCEvent*)next())) {
604             if (oldidx < (event->GetSecondaryOffset() + (event->GetNumberOfTracks() - event->GetNumberOfPrimaries()))) break;
605         }
606         return (oldidx - event->GetSecondaryOffset() + event->GetNumberOfPrimaries());
607     }
608 }
609
610 Int_t AliMCEvent::BgLabelToIndex(Int_t label)
611 {
612     // Convert a background label to an absolute index
613     if (fSubsidiaryEvents) {
614         AliMCEvent* bgEvent = (AliMCEvent*) (fSubsidiaryEvents->At(1));
615         label -= BgLabelOffset();
616         if (label < bgEvent->GetNumberOfPrimaries()) {
617             label += bgEvent->GetPrimaryOffset();
618         } else {
619             label += (bgEvent->GetSecondaryOffset() - fNprimaries);
620         }
621     }
622     return (label);
623 }
624
625
626 Bool_t AliMCEvent::IsPhysicalPrimary(Int_t i) 
627 {
628 //
629 // Delegate to subevent if necesarry 
630
631     
632     if (!fSubsidiaryEvents) {
633         return fStack->IsPhysicalPrimary(i);
634     } else {
635         AliMCEvent* evt = 0;
636         Int_t idx = FindIndexAndEvent(i, evt);
637         return (evt->IsPhysicalPrimary(idx));
638     }
639 }
640
641 void AliMCEvent::InitEvent()
642 {
643 //
644 // Initialize the subsidiary event structure
645     if (fSubsidiaryEvents) {
646         TIter next(fSubsidiaryEvents);
647         AliMCEvent* evt;
648         fNprimaries = 0;
649         fNparticles = 0;
650         
651         while((evt = (AliMCEvent*)next())) {
652             fNprimaries += evt->GetNumberOfPrimaries(); 
653             fNparticles += evt->GetNumberOfTracks();    
654         }
655         
656         Int_t ioffp = 0;
657         Int_t ioffs = fNprimaries;
658         next.Reset();
659         
660         while((evt = (AliMCEvent*)next())) {
661             evt->SetPrimaryOffset(ioffp);
662             evt->SetSecondaryOffset(ioffs);
663             ioffp += evt->GetNumberOfPrimaries();
664             ioffs += (evt->GetNumberOfTracks() - evt->GetNumberOfPrimaries());      
665         }
666     }
667 }
668 void AliMCEvent::PreReadAll()                              
669 {
670     // Preread the MC information
671     Int_t i;
672     // secondaries
673     for (i = fStack->GetNprimary(); i < fStack->GetNtrack(); i++) 
674     {
675         GetTrack(i);
676     }
677     // primaries
678     for (i = 0; i < fStack->GetNprimary(); i++) 
679     {
680         GetTrack(i);
681     }
682     
683     
684 }
685
686
687 const AliVVertex * AliMCEvent::GetPrimaryVertex() const 
688 {
689     // Create a MCVertex object from the MCHeader information
690     if (!fVertex) {
691         TArrayF v;
692         GenEventHeader()->PrimaryVertex(v) ;
693         fVertex = new AliMCVertex(v[0], v[1], v[2]);
694     }
695     
696     return fVertex;
697 }
698
699
700 ClassImp(AliMCEvent)