]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveBase/AliEveVSDCreator.cxx
Extra header added to the list
[u/mrichter/AliRoot.git] / EVE / EveBase / AliEveVSDCreator.cxx
1 // $Id$
2 // Main authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
3
4 /**************************************************************************
5  * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
6  * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
7  * full copyright notice.                                                 *
8  **************************************************************************/
9
10 #include <memory>
11 #include "AliEveVSDCreator.h"
12
13 #include "AliEveEventManager.h"
14
15 #include <AliStack.h>
16 #include <AliITSLoader.h>
17 #include <AliTPCTrackHitsV2.h>
18 #include <AliPDG.h>
19 #include <AliHit.h>
20 #include <AliESDEvent.h>
21 #include <AliESDv0.h>
22 #include <AliTPCclusterMI.h>
23 #include <AliTPCClustersRow.h>
24 #include <AliITSclusterV2.h>
25 #include <AliESDkink.h>
26 #include <AliESDtrack.h>
27
28 #include <AliRun.h>
29 #include <AliRunLoader.h>
30 #include <AliTPCParam.h>
31
32 #include <TDatabasePDG.h>
33 #include <TSystem.h>
34 #include <TFile.h>
35
36 //______________________________________________________________________________
37 //
38 // Create VSD file from ALICE data.
39
40 ClassImp(AliEveVSDCreator)
41
42 AliEveVSDCreator::AliEveVSDCreator(const Text_t* name, const Text_t* title) :
43   TEveVSD(name, title),
44
45   fTPCHitRes  (2),
46   fTRDHitRes  (1),
47
48   fDebugLevel (0),
49   fRunLoader  (0),
50   fGenInfoMap ()
51 {
52   // Constructor.
53
54   // Particles not in ROOT's PDG database occuring in ALICE
55   AliPDG::AddParticlesToPdgDataBase();
56 }
57
58 /******************************************************************************/
59
60 void AliEveVSDCreator::CreateVSD(const Text_t* vsdFile)
61 {
62   // Create the VSD for current event in AliEveEventManager.
63   // Result is stored in directory "Event%04d".
64
65   static const TEveException kEH("AliEveVSDCreator::CreateVSD ");
66
67   AliEveEventManager::AssertGeometry();
68
69   fRunLoader = AliEveEventManager::AssertRunLoader();
70
71   if (fDebugLevel > 0)
72     printf("%s open seems ok. Now loading sim data.\n", kEH.Data());
73
74   fRunLoader->LoadHeader();
75   fRunLoader->LoadKinematics();
76   fRunLoader->LoadTrackRefs();
77   fRunLoader->LoadHits();
78
79   if (fDebugLevel > 0)
80     printf("%s opening output TEveVSD.\n", kEH.Data());
81
82   TFile* file = TFile::Open(vsdFile, "UPDATE", "ALICE Visualization Summary Data");
83
84   Int_t curEvent = AliEveEventManager::CurrentEventId();
85
86   TString eventDir; eventDir.Form("Event%04d", curEvent);
87
88   if (file->Get(eventDir))
89   {
90     Warning(kEH, "Event-dir '%s' existed -- replacing.", eventDir.Data());
91     file->Delete(eventDir + ";*");
92   }
93
94   fDirectory = new TDirectoryFile(eventDir, Form("Data for event %d", curEvent));
95
96   if (fDebugLevel > 0)
97     printf("%s creating trees now ...\n", kEH.Data());
98
99   CreateTrees();
100
101   if (fDebugLevel > 0)
102     printf("%s trees created, closing files.\n", kEH.Data());
103
104   file->Write();
105   file->Close();
106   delete file;
107   fDirectory = 0;
108
109   // clean after the TEveVSD data was sucessfuly written
110   fTreeK      = 0;
111   fTreeH      = 0;
112   //fTreeTR     = 0;
113   fTreeC      = 0;
114   fTreeV0     = 0;
115   fTreeKK     = 0;
116   fTreeR      = 0;
117   fTreeGI     = 0;
118
119   fRunLoader = 0;
120
121   if (fDebugLevel > 0)
122     printf("%s all done.\n", kEH.Data());
123 }
124
125 void AliEveVSDCreator::CreateTrees()
126 {
127   // Create and fill the output trees by calling all the
128   // ConvertXyzz() functions.
129   // Exceptions from individual functions are displayed as warnings.
130
131   static const TEveException kEH("AliEveVSDCreator::CreateTrees ");
132
133   if (fDirectory == 0)
134     throw kEH + "output directory not set.";
135
136   try
137   {
138     if (fDebugLevel > 1)
139       printf("%sConvertKinematics.\n", kEH.Data());
140     ConvertKinematics();
141   }
142   catch(TEveException& exc)
143   {
144     Warning(kEH, "%s", exc.Data());
145   }
146
147   try
148   {
149     if (fDebugLevel > 1)
150       printf("%sConvertHits.\n", kEH.Data());
151     ConvertHits();
152   }
153   catch(TEveException& exc)
154   {
155     Warning(kEH, "%s", exc.Data());
156   }
157
158   try
159   {
160     if (fDebugLevel > 1)
161       printf("%sConvertClusters.\n", kEH.Data());
162     ConvertClusters();
163   }
164   catch(TEveException& exc)
165   {
166     Warning(kEH, "%s", exc.Data());
167   }
168
169   try 
170   {
171     if (fDebugLevel > 1)
172       printf("%sConvertRecTracks.\n", kEH.Data());
173     ConvertRecTracks();
174   }
175   catch(TEveException& exc)
176   {
177     Warning(kEH, "Executing %s skipping AliEveV0 extraction.", exc.Data());
178     goto end_esd_processing;
179   }
180
181   // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
182
183   Warning(kEH, "Explicitly abandoning further conversion.");
184   return;
185
186   // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
187 /*
188   try
189   {
190     if (fDebugLevel > 1)
191       printf("%sConvertV0.\n", kEH.Data());
192     ConvertV0();
193   }
194   catch(TEveException& exc)
195   {
196     Warning(kEH, exc);
197   }
198
199   try
200   {
201     if (fDebugLevel > 1)
202       printf("%sConvertKinks.\n", kEH.Data());
203     ConvertKinks();
204   }
205   catch(TEveException& exc)
206   {
207     Warning(kEH, exc);
208   }
209 */
210 end_esd_processing:
211
212   try
213   {
214     if (fDebugLevel > 1)
215       printf("%sConvertGenInfo.\n", kEH.Data());
216     ConvertGenInfo();
217   }
218   catch(TEveException& exc)
219   {
220     Warning(kEH, "%s", exc.Data());
221   }
222
223   return;
224 }
225
226 /******************************************************************************/
227 // Kinematics
228 /******************************************************************************/
229
230 void AliEveVSDCreator::ConvertKinematics()
231 {
232   // Convert kinematics.
233   // Track references are not stored, they should be.
234
235   static const TEveException kEH("AliEveVSDCreator::ConvertKinematics ");
236
237   if (fTreeK != 0)
238     throw kEH + "kinematics already converted";
239
240   AliStack* stack = fRunLoader->Stack();
241   if (stack == 0)
242     throw kEH + "stack is null.";
243
244   fDirectory->cd();
245   fTreeK = new TTree("Kinematics", "Particles created during simulation.");
246
247   Int_t nentries = stack->GetNtrack();
248   std::vector<TEveMCTrack>  vmc(nentries);
249   for (Int_t idx = 0; idx < nentries; ++idx)
250   {
251     TParticle*   tp = stack->Particle(idx);
252     vmc[idx]        = *tp;
253     vmc[idx].fLabel = idx;
254   }
255
256   // read track refrences
257   // functionality now in AliEveKineTools.
258   /*
259   TTree* fTreeTR =  fRunLoader->TreeTR();
260
261   if(fTreeTR == 0) {
262     Warning(kEH, "no TrackRefs; some data will not be available.");
263   } else {
264     TClonesArray* RunArrayTR = 0;
265     fTreeTR->SetBranchAddress("AliRun", &RunArrayTR);
266
267     Int_t nPrimaries = (Int_t) fTreeTR->GetEntries();
268     for (Int_t iPrimPart = 0; iPrimPart<nPrimaries; iPrimPart++) {
269       // printf("T0 fTreeTR->GetEntry(%d) \n",iPrimPart);
270       fTreeTR->GetEntry(iPrimPart);
271       // printf("END fTreeTR->GetEntry(%d) \n",iPrimPart);
272
273       for (Int_t iTrackRef = 0; iTrackRef < RunArrayTR->GetEntriesFast(); iTrackRef++) {
274         AliTrackReference *trackRef = (AliTrackReference*)RunArrayTR->At(iTrackRef);
275         Int_t track = trackRef->GetTrack();
276         if(track < nentries && track > 0){
277           TEveMCTrack& mct = vmc[track];
278           if(trackRef->TestBit(kNotDeleted)) {
279             mct.decayed   = true;
280             mct.t_decay   = trackRef->GetTime();
281             mct.V_decay.x = trackRef->X();
282             mct.V_decay.y = trackRef->Y();
283             mct.V_decay.z = trackRef->Z();
284             mct.P_decay.x = trackRef->Px();
285             mct.P_decay.y = trackRef->Py();
286             mct.P_decay.z = trackRef->Pz();
287             if(TMath::Abs(mct.GetPdgCode()) == 11)
288               mct.decayed = false; // a bug in TreeTR
289           }
290         }
291       }
292     }
293   }
294   */
295
296   fTreeK->Branch("K", "TEveMCTrack",  &fpK, fBuffSize);
297
298   printf("sizeofvmc = %d\n", (Int_t) vmc.size());
299   for(std::vector<TEveMCTrack>::iterator k = vmc.begin(); k != vmc.end(); ++k)
300   {
301     TEveMCTrack& mct = *k;
302     fK = mct;
303
304     TParticle* m  = &mct;
305     Int_t      mi = mct.fLabel;
306     int cnt = 0;
307     while (m->GetMother(0) != -1)
308     {
309       if (cnt > 100)
310         printf("cnt %d mi=%d, mo=%d\n", cnt, mi, m->GetMother(0));
311
312       mi = m->GetMother(0);
313       m = &vmc[mi];
314       ++cnt;
315     }
316     fK.fEvaLabel = mi;
317
318     fTreeK->Fill();
319   }
320
321   fTreeK->BuildIndex("fLabel");
322 }
323
324
325 //==============================================================================
326 // Hits
327 //==============================================================================
328
329 void AliEveVSDCreator::ConvertHits()
330 {
331   // Convert MC hits.
332   // TPC hits are handled specially as they are compressed - only mayor
333   // hits are stored
334
335   static const TEveException kEH("AliEveVSDCreator::ConvertHits ");
336
337   if (fTreeH != 0)
338     throw kEH + "hits already converted.";
339
340   fDirectory->cd();
341   fTreeH = new TTree("Hits", "Simulated energy depositions in detector.");
342   fTreeH->Branch("H", "TEveHit", &fpH, fBuffSize);
343
344   try
345   {
346     ConvertAnyHits("ITS", "AliITShit", 0, 0);
347     ConvertTPCHits(                    1, fTPCHitRes*fTPCHitRes);
348     ConvertAnyHits("TRD", "AliTRDhit", 2, fTRDHitRes*fTRDHitRes);
349     ConvertAnyHits("TOF", "AliTOFhit", 3, 0);
350
351   }
352   catch(TEveException& exc)
353   {
354     Warning(kEH, "%s", exc.Data());
355   }
356 }
357
358 //------------------------------------------------------------------------------
359
360 void AliEveVSDCreator::ConvertAnyHits(const TString& detector,
361                                       const TString& clones_class,
362                                       Int_t det_id, Float_t minDistSqr)
363 {
364   // Convert hits for detector.
365
366   static const TEveException kEH("AliEveVSDCreator::ConvertAnyHits ");
367
368   Float_t x,y,z, x1,y1,z1;
369
370   TTree* treeh = fRunLoader->GetTreeH(detector, kFALSE);
371   if (treeh == 0)
372   {
373     Warning(kEH, "no hits for %s.", detector.Data());
374     return;
375   }
376
377   TClonesArray *arr = new TClonesArray(clones_class);
378   treeh->SetBranchAddress(detector, &arr);
379   Int_t np = treeh->GetEntries();
380
381   // In TreeH files hits are grouped in clones arrays
382   // each eva particle has its own clone array.
383
384   std::map<Int_t, Int_t> hmap;
385
386   for (Int_t i = 0; i < np; ++i)
387   {
388     treeh->GetEntry(i);
389     Int_t evaIdx = np - i - 1;
390     Int_t nh = arr->GetEntriesFast();
391     x = y = z = 0;
392     // printf("%d entry %d hits for primary %d \n", i, nh, evaIdx);
393     for (Int_t j = 0; j < nh; ++j)
394     {
395       AliHit* aliHit = (AliHit*)arr->UncheckedAt(j);
396
397       x1 = aliHit->X(); y1 = aliHit->Y(); z1 = aliHit->Z();
398       if (minDistSqr == 0 || (x-x1)*(x-x1) + (y-y1)*(y-y1) + (z-z1)*(z-z1) > minDistSqr)
399       {
400         fH.fDetId      = det_id;
401         fH.fSubdetId   = 0;
402         fH.fLabel      = aliHit->GetTrack();
403         fH.fEvaLabel   = evaIdx;
404         fH.fV.Set(aliHit->X(), aliHit->Y(), aliHit->Z());
405
406         fTreeH->Fill();
407
408         hmap[fH.fLabel]++;
409
410         x=x1; y=y1; z=z1;
411       }
412     }
413   }
414   delete arr;
415
416   // Set geninfo.
417   for (std::map<Int_t, Int_t>::iterator j = hmap.begin(); j != hmap.end(); ++j)
418   {
419     GetGeninfo(j->first)->fNHits += j->second;
420   }
421 }
422
423 //------------------------------------------------------------------------------
424
425 void AliEveVSDCreator::ConvertTPCHits(Int_t det_id, Float_t minDistSqr)
426 {
427   // Convert hits for TPC.
428
429   static const TEveException kEH("AliEveVSDCreator::ConvertTPCHits ");
430
431   Float_t x,y,z, x1,y1,z1;
432
433   TTree* treeh = fRunLoader->GetTreeH("TPC", false);
434   if (treeh == 0)
435   {
436     Warning(kEH, "no hits for %s.", "TPC");
437     return;
438   }
439
440   AliTPCTrackHitsV2 hv2, *hv2p = &hv2;
441   treeh->SetBranchAddress("TPC2", &hv2p);
442   Int_t np = treeh->GetEntries();
443
444   std::map<Int_t, Int_t> hmap;
445
446   for (Int_t i = 0; i < np; ++i)
447   {
448     treeh->GetEntry(i);
449     Int_t evaIdx = np -i -1;
450     if (hv2.First() == 0) continue;
451     x = y = z = 0;
452     do
453     {
454       AliHit* ah = hv2.GetHit();
455       x1 = ah->X(); y1 = ah->Y(); z1 = ah->Z();
456       if (minDistSqr == 0 || (x-x1)*(x-x1) + (y-y1)*(y-y1) + (z-z1)*(z-z1) > minDistSqr)
457       {
458         fH.fDetId    = det_id;
459         fH.fSubdetId = 0;
460         fH.fLabel    = ah->Track();
461         fH.fEvaLabel = evaIdx;
462         fH.fV.fX = x1; fH.fV.fY = y1; fH.fV.fZ = z1;
463
464         fTreeH->Fill();
465
466         hmap[fH.fLabel]++;
467
468         x = x1; y = y1; z = z1;
469       }
470     } while (hv2.Next());
471   }
472
473   // Set geninfo.
474   for (std::map<Int_t, Int_t>::iterator j = hmap.begin(); j != hmap.end(); ++j)
475   {
476     GetGeninfo(j->first)->fNHits += j->second;
477   }
478 }
479
480
481 //==============================================================================
482 // Clusters
483 //==============================================================================
484
485 void AliEveVSDCreator::ConvertClusters()
486 {
487   // Convert clusters.
488   //
489   // Only supported for ITS and TPC at the moment, see dedicated
490   // functions ConvertITSClusters() and ConvertTPCClusters().
491   //
492   // It should be possible now to do this in a general manner (with
493   // the alignment framework).
494
495   static const TEveException kEH("AliEveVSDCreator::ConvertClusters ");
496
497   if (fTreeC != 0)
498     throw kEH + "clusters already converted.";
499
500   fDirectory->cd();
501   fTreeC = new TTree("Clusters", "Reconstructed points of particle passage.");
502   fTreeC->Branch("C", "TEveCluster", &fpC, fBuffSize);
503
504   try
505   {
506     ConvertAnyClusters("ITS", "ITSRecPoints", 0);
507     ConvertTPCClusters(                       1);
508     ConvertAnyClusters("TRD", "TRDcluster",   2);
509     ConvertAnyClusters("TOF", "TOF",          3);
510
511   }
512   catch(TEveException& exc)
513   {
514     Warning(kEH, "%s", exc.Data());
515   }
516 }
517
518 //------------------------------------------------------------------------------
519
520 void AliEveVSDCreator::ConvertAnyClusters(const TString& detector,
521                                           const TString& branch_name,
522                                           Int_t det_id)
523 {
524   // Convert clusters for detector and transform them to global coordinates.
525
526   static const TEveException kEH("AliEveVSDCreator::ConvertITSClusters ");
527
528   fRunLoader->LoadRecPoints(detector);
529   TTree* tree = fRunLoader->GetTreeR(detector, false);
530   if (!tree)
531     throw kEH + "'TreeR' not found.";
532
533   TClonesArray *cl = 0;
534   TBranch *branch  = tree->GetBranch(branch_name);
535   branch->SetAddress(&cl);
536
537   Int_t nmods = branch->GetEntries();
538
539   std::map<Int_t, Int_t> cmap;
540
541   for (Int_t mod = 0; mod < nmods; ++mod)
542   {
543     branch->GetEntry(mod);
544     Int_t nc = cl->GetEntriesFast();
545     for (Int_t j = 0; j < nc; ++j)
546     {
547       AliCluster *c = (AliCluster*) cl->UncheckedAt(j);
548
549       fC.fDetId    = det_id;
550       fC.fSubdetId = mod;
551       fC.fLabel[0] = c->GetLabel(0);
552       fC.fLabel[1] = c->GetLabel(1);
553       fC.fLabel[2] = c->GetLabel(2);
554
555       c->GetGlobalXYZ(fC.fV.Arr());
556
557       fTreeC->Fill();
558
559       {
560         int i = 0;
561         while (i < 3 && fC.fLabel[i])
562           cmap[fC.fLabel[i++]]++;
563       }
564     }
565   }
566   delete cl;
567
568   for (std::map<Int_t, Int_t>::iterator j=cmap.begin(); j!=cmap.end(); ++j)
569   {
570     GetGeninfo(j->first)->fNClus += j->second;
571   }
572 }
573
574 //------------------------------------------------------------------------------
575
576 void AliEveVSDCreator::ConvertTPCClusters(Int_t det_id)
577 {
578   // Convert TPC clusters and transform them to global coordinates.
579
580   static const TEveException kEH("AliEveVSDCreator::ConvertTPCClusters ");
581
582   const Int_t kMaxCl = 100*160;
583
584   fRunLoader->LoadRecPoints("TPC");
585   TTree* tree = fRunLoader->GetTreeR("TPC", false);
586   if (!tree)
587     throw kEH + "'TreeR' not found.";
588
589   AliTPCClustersRow *clrow = new AliTPCClustersRow;
590   clrow->SetClass("AliTPCclusterMI");
591   clrow->SetArray(kMaxCl);
592   tree->SetBranchAddress("Segment", &clrow);
593
594   Int_t nEnt = tree->GetEntries();
595
596   std::map<Int_t, Int_t> cmap;
597
598   for (Int_t j = 0; j < nEnt; ++j)
599   {
600     if (!tree->GetEvent(j)) continue;
601
602     TClonesArray *cl = clrow->GetArray();
603     Int_t ncl = cl->GetEntriesFast();
604
605     while (ncl--)
606     {
607       AliCluster *c = (AliCluster*) cl->UncheckedAt(ncl);
608
609       fC.fDetId    = det_id;
610       fC.fSubdetId = j;
611       fC.fLabel[0] = c->GetLabel(0);
612       fC.fLabel[1] = c->GetLabel(1);
613       fC.fLabel[2] = c->GetLabel(2);
614
615       c->GetGlobalXYZ(fC.fV.Arr());
616
617       fTreeC->Fill();
618
619       {
620         int i = 0;
621         while (i < 3 && fC.fLabel[i])
622           cmap[fC.fLabel[i++]]++;
623       }
624     }
625
626     cl->Clear();
627   }
628   delete clrow;
629
630   //set geninfo
631   for (std::map<Int_t, Int_t>::iterator j=cmap.begin(); j!=cmap.end(); ++j)
632   {
633     GetGeninfo(j->first)->fNClus += j->second;
634   }
635 }
636
637
638 /******************************************************************************/
639 // ESD
640 /******************************************************************************/
641
642 void AliEveVSDCreator::ConvertRecTracks()
643 {
644   // Convert reconstructed tracks.
645
646   static const TEveException kEH("AliEveVSDCreator::ConvertRecTracks ");
647
648   if (fTreeR != 0)
649     throw kEH + "tracks already converted.";
650
651   AliESDEvent* esdEvent = AliEveEventManager::AssertESD();
652
653   fDirectory->cd();
654   fTreeR =  new TTree("RecTracks", "Reconstructed particle trajectories.");
655   fTreeR->Branch("R", "TEveRecTrack", &fpR, 512*1024);
656
657   // reconstructed tracks
658   AliESDtrack* esdTrack;
659   Double_t     dbuf[3];
660   for (Int_t n = 0; n < esdEvent->GetNumberOfTracks(); ++n)
661   {
662     esdTrack = esdEvent->GetTrack(n);
663
664     fR.fLabel  = esdTrack->GetLabel();
665     fR.fStatus = (Int_t) esdTrack->GetStatus();
666     fR.fSign   = (Int_t) esdTrack->GetSign();
667     esdTrack->GetXYZ(dbuf);    fR.fV.Set(dbuf);
668     esdTrack->GetPxPyPz(dbuf); fR.fP.Set(dbuf);
669     Double_t ep = esdTrack->GetP();
670     fR.fBeta = ep/TMath::Sqrt(ep*ep + TMath::C()*TMath::C()*esdTrack->GetMass()*esdTrack->GetMass());
671     fTreeR->Fill();
672   }
673
674   fTreeR->BuildIndex("fLabel");
675 }
676
677 /******************************************************************************/
678
679 void AliEveVSDCreator::ConvertV0()
680 {
681   // Convert reconstructed V0s.
682
683   static const TEveException kEH("AliEveVSDCreator::ConvertV0 ");
684
685   if (fTreeV0 != 0)
686     throw kEH + "AliEveV0 already converted.";
687
688   AliESDEvent* esdEvent = AliEveEventManager::AssertESD();
689
690   fDirectory->cd();
691   fTreeV0 =  new TTree("AliEveV0", "AliEveV0 points");
692
693   fTreeV0->Branch("AliEveV0", "TEveRecV0", &fpV0, 512*1024,1);
694
695   for (Int_t n = 0; n < esdEvent->GetNumberOfV0s(); ++n)
696   {
697     AliESDv0    *av     = esdEvent->GetV0(n);
698     AliESDtrack *trackN = esdEvent->GetTrack(av->GetNindex()); // negative daughter
699     AliESDtrack *trackP = esdEvent->GetTrack(av->GetPindex()); // positive daughter
700
701     Double_t pos[3];
702
703     fV0.fStatus = av->GetStatus();
704     // Point of closest approach
705     av->GetXYZ(pos[0],pos[1],pos[2]);
706     fV0.fVCa.fX = pos[0];
707     fV0.fVCa.fY = pos[1];
708     fV0.fVCa.fZ = pos[2];
709     // set birth vertex of neutral particle
710     av->GetXYZ(pos[0], pos[1], pos[2]);
711     fV0.fV0Birth.Set(pos);
712
713     // momentum and position of negative particle
714     av->GetParamN()->GetPxPyPz(pos);
715     fV0.fPNeg.Set(pos);
716     av->GetParamN()->GetXYZ(pos);
717     fV0.fVNeg.Set(pos);
718
719     // momentum and position of positive particle
720     av->GetParamP()->GetPxPyPz(pos);
721     fV0.fPPos.Set(pos);
722     av->GetParamP()->GetXYZ(pos);
723     fV0.fVPos.Set(pos);
724
725     fV0.fLabel = 0; // !!!! mother label unknown
726     fV0.fPdg   = av->GetPdgCode();
727
728     // daughter indices
729     fV0.fDLabel[0] = TMath::Abs(trackN->GetLabel());
730     fV0.fDLabel[1] = TMath::Abs(trackP->GetLabel());
731
732     // printf("AliEveV0 convert labels(%d,%d) index(%d,%d)\n",
733     //     fV0.d_label[0],  fV0.d_label[1],
734     //     av->GetNIndex(), av->GetPIndex());
735
736     fTreeV0->Fill();
737   }
738   // if (esdEvent->GetNumberOfV0s()) fTreeV0->BuildIndex("label");
739   delete esdEvent;
740 }
741
742 /******************************************************************************/
743
744 void AliEveVSDCreator::ConvertKinks()
745 {
746   // Convert reconstructed kinks.
747
748   static const TEveException kEH("AliEveVSDCreator::ConvertKinks ");
749
750   if (fTreeKK != 0)
751     throw kEH + "Kinks already converted.";
752
753   throw kEH + "Currently non-supported - TEveRecKink being updated.";
754
755   /*
756   AliESDEvent* esdEvent = AliEveEventManager::AssertESD();
757
758   fDirectory->cd();
759   fTreeKK =  new TTree("Kinks", "ESD Kinks");
760
761   fTreeKK->Branch("KK", "TEveRecKink", &fpKK, fBuffSize);
762
763   //  printf("CONVERT KINK Read %d entries in tree kinks \n",  esdEvent->GetNumberOfKinks());
764   for (Int_t n = 0; n < esdEvent->GetNumberOfKinks(); ++n)
765   {
766     AliESDkink* kk = esdEvent->GetKink(n);
767
768     Double_t pos[3];
769
770     fKK.fLabel  = kk->GetLabel(0);
771     fKK.fStatus = 0; // status is Char_t[12] ... have no idea how/what to extract.
772
773     // reconstructed kink position
774     fKK.fLabelSec = kk->GetLabel(1);
775     fKK.fVKink.Set(kk->GetPosition());
776
777     const AliExternalTrackParam& tpMother = kk->RefParamMother();
778     // momentum and position of mother
779     tpMother.GetPxPyPz(pos);
780     fKK.fP.Set(pos);
781     tpMother.GetXYZ(pos);
782     fKK.fV.Set(pos);
783     const Double_t* par =  tpMother.GetParameter();
784     // printf("KINK Pt %f, %f \n",1/tpMother.Pt(),par[4] );
785     fKK.fSign = (par[4] < 0) ? -1 : 1;
786
787     const AliExternalTrackParam& tpDaughter = kk->RefParamDaughter();
788     // momentum and position of daughter
789     tpDaughter.GetPxPyPz(pos);
790     fKK.fPSec.Set(pos);
791     tpDaughter.GetXYZ(pos);
792     fKK.fVEnd.Set(pos);
793
794     fTreeKK->Fill();
795   }
796   if (esdEvent->GetNumberOfKinks()) fTreeKK->BuildIndex("label");
797   delete esdEvent;
798   */
799 }
800 /******************************************************************************/
801 // TEveMCRecCrossRef
802 /******************************************************************************/
803
804 void AliEveVSDCreator::ConvertGenInfo()
805 {
806   // Build simulation-reconstruction cross-reference table.
807   // In a rather poor state at the moment.
808
809   static const TEveException kEH("AliEveVSDCreator::ConvertGenInfo ");
810
811   if (fTreeGI != 0)
812     throw kEH + "GI already converted.";
813
814   fDirectory->cd();
815   fTreeGI = new TTree("TEveMCRecCrossRef", "Objects prepared for cross querry");
816
817   TEveMCRecCrossRef::Class()->IgnoreTObjectStreamer(true);
818   fTreeGI->Branch("GI", "TEveMCRecCrossRef",  &fpGI, fBuffSize);
819   fTreeGI->Branch("K.", "TEveMCTrack",  &fpK);
820   fTreeGI->Branch("R.", "TEveRecTrack", &fpR);
821
822   for (std::map<Int_t, TEveMCRecCrossRef*>::iterator j=fGenInfoMap.begin(); j!=fGenInfoMap.end(); ++j)
823   {
824     fGI        = *(j->second);
825     fGI.fLabel = j->first;
826     fTreeK->GetEntry(j->first);
827
828     if (fTreeR) {
829       Int_t re = fTreeR->GetEntryNumberWithIndex(j->first);
830       if(re != -1)
831         fGI.fIsRec = true;
832     }
833     //    Int_t hasV0 =  fTreeV0->GetEntryNumberWithIndex(j->first);
834     //if (hasV0 != -1)
835     //  fGI.has_AliEveV0 = true;
836     if (fTreeKK)
837     {
838       Int_t hasKk =  fTreeKK->GetEntryNumberWithIndex(j->first);
839       if (hasKk != -1)
840         fGI.fHasKink = true;
841     }
842     fTreeGI->Fill();
843   }
844   fGenInfoMap.clear();
845 }
846
847 /******************************************************************************/
848 /******************************************************************************/
849 // Protected methods
850 /******************************************************************************/
851 /******************************************************************************/
852
853 TEveMCRecCrossRef* AliEveVSDCreator::GetGeninfo(Int_t label)
854 {
855   // Return the cross-reference structure for given label.
856   // If the object does not exist it is created.
857
858   // printf("get_geninfo %d\n", label);
859   TEveMCRecCrossRef* gi;
860   std::map<Int_t, TEveMCRecCrossRef*>::iterator i = fGenInfoMap.find(label);
861   if (i == fGenInfoMap.end())
862   {
863     gi =  new TEveMCRecCrossRef();
864     fGenInfoMap[label] = gi;
865   }
866   else
867   {
868     gi = i->second;
869   }
870   return gi;
871 }