]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliMCEvent.cxx
b31909cc0b235e6a29cef7f9e781a4b74647ea75
[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
31 #include "AliLog.h"
32 #include "AliMCEvent.h"
33 #include "AliStack.h"
34 #include "AliTrackReference.h"
35 #include "AliHeader.h"
36 #include "AliGenEventHeader.h"
37
38
39 AliMCEvent::AliMCEvent():
40     AliVEvent(),
41     fStack(0),
42     fMCParticles(new TClonesArray("AliMCParticle",1000)),
43     fMCParticleMap(0),
44     fHeader(0),
45     fTRBuffer(0),
46     fTrackReferences(new TClonesArray("AliTrackReference", 1000)),
47     fTreeTR(0),
48     fTmpTreeTR(0),
49     fTmpFileTR(0),
50     fNprimaries(-1),
51     fNparticles(-1)
52 {
53     // Default constructor
54 }
55
56 AliMCEvent::AliMCEvent(const AliMCEvent& mcEvnt) :
57     AliVEvent(mcEvnt),
58     fStack(0),
59     fMCParticles(0),
60     fMCParticleMap(0),
61     fHeader(0),
62     fTRBuffer(0),
63     fTrackReferences(0),
64     fTreeTR(0),
65     fTmpTreeTR(0),
66     fTmpFileTR(0),
67     fNprimaries(-1),
68     fNparticles(-1) 
69
70 // Copy constructor
71 }
72
73
74 AliMCEvent& AliMCEvent::operator=(const AliMCEvent& mcEvnt)
75 {
76     // assignment operator
77     if (this!=&mcEvnt) { 
78         AliVEvent::operator=(mcEvnt); 
79     }
80   
81     return *this; 
82 }
83
84 void AliMCEvent::ConnectTreeE (TTree* tree)
85 {
86     // Connect the event header tree
87     tree->SetBranchAddress("Header", &fHeader);
88 }
89
90 void AliMCEvent::ConnectTreeK (TTree* tree)
91 {
92     // Connect the kinematics tree to the stack
93     fStack = fHeader->Stack();
94     fStack->ConnectTree(tree);
95     //
96     // Load the event
97     fStack->GetEvent();
98     fNparticles = fStack->GetNtrack();
99     fNprimaries = fStack->GetNprimary();
100     Int_t iev  = fHeader->GetEvent();
101     Int_t ievr = fHeader->GetEventNrInRun();
102     
103     AliInfo(Form("AliMCEvent# %5d %5d: Number of particles: %5d (all) %5d (primaries)\n", 
104                  iev, ievr, fNparticles, fNprimaries));
105  
106     // This is a cache for the TParticles converted to MCParticles on user request
107     if (fMCParticleMap) {
108         fMCParticleMap->Clear();
109         if (fNparticles>0) fMCParticleMap->Expand(fNparticles);}
110     else
111         fMCParticleMap = new TRefArray(fNparticles);
112
113 }
114
115 void AliMCEvent::ConnectTreeTR (TTree* tree)
116 {
117     // Connect the track reference tree
118     fTreeTR = tree;
119     
120     if (fTreeTR->GetBranch("AliRun")) {
121         if (fTmpFileTR) {
122             fTmpFileTR->Close();
123             delete fTmpFileTR;
124         }
125         // This is an old format with one branch per detector not in synch with TreeK
126         ReorderAndExpandTreeTR();
127     } else {
128         // New format 
129         fTreeTR->SetBranchAddress("TrackReferences", &fTRBuffer);
130     }
131 }
132
133 Int_t AliMCEvent::GetParticleAndTR(Int_t i, TParticle*& particle, TClonesArray*& trefs)
134 {
135     // Retrieve entry i
136     if (i < 0 || i >= fNparticles) {
137         AliWarning(Form("AliMCEventHandler::GetEntry: Index out of range"));
138         particle = 0;
139         trefs    = 0;
140         return (-1);
141     }
142     particle = fStack->Particle(i);
143     if (fTreeTR) {
144         fTreeTR->GetEntry(fStack->TreeKEntry(i));
145         trefs    = fTRBuffer;
146         return trefs->GetEntries();
147     } else {
148         trefs = 0;
149         return -1;
150     }
151 }
152
153
154 void AliMCEvent::Clean()
155 {
156     // Clean-up before new trees are connected
157     
158     if (fHeader) {
159         delete fHeader;
160         fHeader = 0;
161     }
162     
163     delete fStack;
164
165     // Clear TR
166     if (fTRBuffer) {
167         fTRBuffer->Clear();
168         delete fTRBuffer;
169         fTRBuffer = 0;
170     }
171 }
172
173 void AliMCEvent::FinishEvent()
174 {
175     // Clean-up after event
176      fStack->Reset(0);
177      fMCParticles->Delete();
178      fTrackReferences->Clear();
179 }
180
181
182
183 void AliMCEvent::DrawCheck(Int_t i, Int_t search)
184 {
185     //
186     // Simple event display for debugging
187     if (!fTreeTR) {
188         AliWarning("No Track Reference information available");
189         return;
190     } 
191     
192     if (i > -1 && i < fNparticles) {
193         fTreeTR->GetEntry(fStack->TreeKEntry(i));
194     } else {
195         AliWarning("AliMCEvent::GetEntry: Index out of range");
196     }
197     
198     Int_t nh = fTRBuffer->GetEntries();
199     
200     
201     if (search) {
202         while(nh <= search && i < fNparticles - 1) {
203             i++;
204             fTreeTR->GetEntry(fStack->TreeKEntry(i));
205             nh =  fTRBuffer->GetEntries();
206         }
207         printf("Found Hits at %5d\n", i);
208     }
209     TParticle* particle = fStack->Particle(i);
210     
211     TH2F*    h = new TH2F("", "", 100, -500, 500, 100, -500, 500);
212     Float_t x0 = particle->Vx();
213     Float_t y0 = particle->Vy();
214
215     Float_t x1 = particle->Vx() + particle->Px() * 50.;
216     Float_t y1 = particle->Vy() + particle->Py() * 50.;
217     
218     TArrow*  a = new TArrow(x0, y0, x1, y1, 0.01);
219     h->Draw();
220     a->SetLineColor(2);
221     
222     a->Draw();
223     
224     for (Int_t ih = 0; ih < nh; ih++) {
225         AliTrackReference* ref = (AliTrackReference*) fTRBuffer->At(ih);
226         TMarker* m = new TMarker(ref->X(), ref->Y(), 20);
227         m->Draw();
228         m->SetMarkerSize(0.4);
229         
230     }
231 }
232
233
234 void AliMCEvent::ReorderAndExpandTreeTR()
235 {
236 //
237 //  Reorder and expand the track reference tree in order to match the kinematics tree.
238 //  Copy the information from different branches into one
239 //
240 //  TreeTR
241
242     fTmpFileTR = new TFile("TrackRefsTmp.root", "recreate");
243     fTmpTreeTR = new TTree("TreeTR", "TrackReferences");
244     if (!fTRBuffer)  fTRBuffer = new TClonesArray("AliTrackReference", 100);
245     fTmpTreeTR->Branch("TrackReferences", "TClonesArray", &fTRBuffer, 32000, 0);
246     
247
248 //
249 //  Activate the used branches only. Otherwisw we get a bad memory leak.
250     fTreeTR->SetBranchStatus("*",        0);
251     fTreeTR->SetBranchStatus("AliRun.*", 1);
252     fTreeTR->SetBranchStatus("ITS.*",    1);
253     fTreeTR->SetBranchStatus("TPC.*",    1);
254     fTreeTR->SetBranchStatus("TRD.*",    1);
255     fTreeTR->SetBranchStatus("TOF.*",    1);
256     fTreeTR->SetBranchStatus("FRAME.*",  1);
257     fTreeTR->SetBranchStatus("MUON.*",   1);
258 //
259 //  Connect the active branches
260     TClonesArray* trefs[7];
261     for (Int_t i = 0; i < 7; i++) trefs[i] = 0;
262     if (fTreeTR){
263         // make branch for central track references
264         if (fTreeTR->GetBranch("AliRun")) fTreeTR->SetBranchAddress("AliRun", &trefs[0]);
265         if (fTreeTR->GetBranch("ITS"))    fTreeTR->SetBranchAddress("ITS",    &trefs[1]);
266         if (fTreeTR->GetBranch("TPC"))    fTreeTR->SetBranchAddress("TPC",    &trefs[2]);
267         if (fTreeTR->GetBranch("TRD"))    fTreeTR->SetBranchAddress("TRD",    &trefs[3]);
268         if (fTreeTR->GetBranch("TOF"))    fTreeTR->SetBranchAddress("TOF",    &trefs[4]);
269         if (fTreeTR->GetBranch("FRAME"))  fTreeTR->SetBranchAddress("FRAME",  &trefs[5]);
270         if (fTreeTR->GetBranch("MUON"))   fTreeTR->SetBranchAddress("MUON",   &trefs[6]);
271     }
272
273     Int_t np = fStack->GetNprimary();
274     Int_t nt = fTreeTR->GetEntries();
275     
276     //
277     // Loop over tracks and find the secondaries with the help of the kine tree
278     Int_t ifills = 0;
279     Int_t it     = 0;
280     Int_t itlast = 0;
281     TParticle* part;
282
283     for (Int_t ip = np - 1; ip > -1; ip--) {
284         part = fStack->Particle(ip);
285 //      printf("Particle %5d %5d %5d %5d %5d %5d \n", 
286 //             ip, part->GetPdgCode(), part->GetFirstMother(), part->GetFirstDaughter(), 
287 //             part->GetLastDaughter(), part->TestBit(kTransportBit));
288
289         // Determine range of secondaries produced by this primary during transport     
290         Int_t dau1  = part->GetFirstDaughter();
291         if (dau1 < np) continue;  // This particle has no secondaries produced during transport
292         Int_t dau2  = -1;
293         if (dau1 > -1) {
294             Int_t inext = ip - 1;
295             while (dau2 < 0) {
296                 if (inext >= 0) {
297                     part = fStack->Particle(inext);
298                     dau2 =  part->GetFirstDaughter();
299                     if (dau2 == -1 || dau2 < np) {
300                         dau2 = -1;
301                     } else {
302                         dau2--;
303                     }
304                 } else {
305                     dau2 = fStack->GetNtrack() - 1;
306                 }
307                 inext--;
308             } // find upper bound
309         }  // dau2 < 0
310         
311
312 //      printf("Check (1) %5d %5d %5d %5d %5d \n", ip, np, it, dau1, dau2);
313 //
314 // Loop over reference hits and find secondary label
315 // First the tricky part: find the entry in treeTR than contains the hits or
316 // make sure that no hits exist.
317 //
318         Bool_t hasHits   = kFALSE;
319         Bool_t isOutside = kFALSE;
320
321         it = itlast;
322         while (!hasHits && !isOutside && it < nt) {
323             fTreeTR->GetEntry(it++);
324             for (Int_t ib = 0; ib < 7; ib++) {
325                 if (!trefs[ib]) continue;
326                 Int_t nh = trefs[ib]->GetEntries();
327                 for (Int_t ih = 0; ih < nh; ih++) {
328                     AliTrackReference* tr = (AliTrackReference*) trefs[ib]->At(ih);
329                     Int_t label = tr->Label();
330                     if (label >= dau1 && label <= dau2) {
331                         hasHits = kTRUE;
332                         itlast = it - 1;
333                         break;
334                     }
335                     if (label > dau2 || label < ip) {
336                         isOutside = kTRUE;
337                         itlast = it - 1;
338                         break;
339                     }
340                 } // hits
341                 if (hasHits || isOutside) break;
342             } // branches
343         } // entries
344
345         if (!hasHits) {
346             // Write empty entries
347             for (Int_t id = dau1; (id <= dau2); id++) {
348                 fTmpTreeTR->Fill();
349                 ifills++;
350             } 
351         } else {
352             // Collect all hits
353             fTreeTR->GetEntry(itlast);
354             for (Int_t id = dau1; (id <= dau2) && (dau1 > -1); id++) {
355                 for (Int_t ib = 0; ib < 7; ib++) {
356                     if (!trefs[ib]) continue;
357                     Int_t nh = trefs[ib]->GetEntries();
358                     for (Int_t ih = 0; ih < nh; ih++) {
359                         AliTrackReference* tr = (AliTrackReference*) trefs[ib]->At(ih);
360                         Int_t label = tr->Label();
361                         // Skip primaries
362                         if (label == ip) continue;
363                         if (label > dau2 || label < dau1) 
364                             printf("AliMCEventHandler::Track Reference Label out of range !: %5d %5d %5d %5d \n", 
365                                    itlast, label, dau1, dau2);
366                         if (label == id) {
367                             // secondary found
368                             tr->SetDetectorId(ib-1);
369                             Int_t nref =  fTRBuffer->GetEntriesFast();
370                             TClonesArray &lref = *fTRBuffer;
371                             new(lref[nref]) AliTrackReference(*tr);
372                         }
373                     } // hits
374                 } // branches
375                 fTmpTreeTR->Fill();
376                 fTRBuffer->Clear();
377                 ifills++;
378             } // daughters
379         } // has hits
380     } // tracks
381
382     //
383     // Now loop again and write the primaries
384     //
385     it = nt - 1;
386     for (Int_t ip = 0; ip < np; ip++) {
387         Int_t labmax = -1;
388         while (labmax < ip && it > -1) {
389             fTreeTR->GetEntry(it--);
390             for (Int_t ib = 0; ib < 7; ib++) {
391                 if (!trefs[ib]) continue;
392                 Int_t nh = trefs[ib]->GetEntries();
393                 // 
394                 // Loop over reference hits and find primary labels
395                 for (Int_t ih = 0; ih < nh; ih++) {
396                     AliTrackReference* tr = (AliTrackReference*)  trefs[ib]->At(ih);
397                     Int_t label = tr->Label();
398                     if (label < np && label > labmax) {
399                         labmax = label;
400                     }
401                     
402                     if (label == ip) {
403                         tr->SetDetectorId(ib-1);
404                         Int_t nref = fTRBuffer->GetEntriesFast();
405                         TClonesArray &lref = *fTRBuffer;
406                         new(lref[nref]) AliTrackReference(*tr);
407                     }
408                 } // hits
409             } // branches
410         } // entries
411         it++;
412         fTmpTreeTR->Fill();
413         fTRBuffer->Clear();
414         ifills++;
415     } // tracks
416     // Check
417
418
419     // Clean-up
420     delete fTreeTR; fTreeTR = 0;
421     
422     for (Int_t ib = 0; ib < 7; ib++) {
423         if (trefs[ib]) {
424             trefs[ib]->Clear();
425             delete trefs[ib];
426             trefs[ib] = 0;
427         }
428     }
429
430     if (ifills != fStack->GetNtrack()) 
431         printf("AliMCEvent:Number of entries in TreeTR (%5d) unequal to TreeK (%5d) \n", 
432                ifills, fStack->GetNtrack());
433
434     fTmpTreeTR->Write();
435     fTreeTR = fTmpTreeTR;
436 }
437
438 AliMCParticle* AliMCEvent::GetTrack(Int_t i) const
439 {
440     // Get MC Particle i
441     AliMCParticle *mcParticle = 0;
442     TParticle     *particle   = 0;
443     TClonesArray  *trefs      = 0;
444     Int_t          ntref      = 0;
445     TRefArray     *rarray     = 0;
446     
447     // Out of range check
448     if (i < 0 || i >= fNparticles) {
449         AliWarning(Form("AliMCEvent::GetEntry: Index out of range"));
450         mcParticle = 0;
451         return (mcParticle);
452     }
453     
454
455
456     //
457     // First check of the MC Particle has been already cached
458     if(!fMCParticleMap->At(i)) {
459         // Get particle from the stack
460         particle   = fStack->Particle(i);
461         // Get track references from Tree TR
462         if (fTreeTR) {
463             fTreeTR->GetEntry(fStack->TreeKEntry(i));
464             trefs     = fTRBuffer;
465             ntref     = trefs->GetEntriesFast();
466             rarray    = new TRefArray(ntref);
467             Int_t nen = fTrackReferences->GetEntriesFast();
468             for (Int_t j = 0; j < ntref; j++) {
469                 // Save the track references in a TClonesArray
470                 AliTrackReference* ref = dynamic_cast<AliTrackReference*>((*fTRBuffer)[j]);
471                 // Save the pointer in a TRefArray
472                 new ((*fTrackReferences)[nen]) AliTrackReference(*ref);
473                 rarray->AddAt((*fTrackReferences)[nen], j);
474                 nen++;
475             } // loop over track references for entry i
476         } // if TreeTR available
477
478         Int_t nentries = fMCParticles->GetEntriesFast();
479         new ((*fMCParticles)[nentries]) AliMCParticle(particle, rarray);
480         mcParticle = dynamic_cast<AliMCParticle*>((*fMCParticles)[nentries]);
481         fMCParticleMap->AddAt(mcParticle, i);
482     } else {
483         mcParticle = dynamic_cast<AliMCParticle*>(fMCParticleMap->At(i));
484     }
485     delete rarray;
486     return mcParticle;
487 }
488
489  AliGenEventHeader* AliMCEvent::GenEventHeader() {return (fHeader->GenEventHeader());}
490
491 ClassImp(AliMCEvent)