Logo ROOT  
Reference Guide
REveVSD.cxx
Go to the documentation of this file.
1// @(#)root/eve7:$Id$
2// Authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007, 2018
3
4/*************************************************************************
5 * Copyright (C) 1995-2019, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#include <ROOT/REveVSD.hxx>
13
14#include "TFile.h"
15#include "TTree.h"
16
17using namespace ROOT::Experimental;
18namespace REX = ROOT::Experimental;
19
20/** \class REveVSD
21\ingroup REve
22Visualization Summary Data - a collection of trees holding standard
23event data in experiment independent format.
24*/
25
26////////////////////////////////////////////////////////////////////////////////
27/// Constructor.
28
29REveVSD::REveVSD(const char* , const char*) :
30 TObject(),
31
32 fFile (nullptr),
33 fDirectory (nullptr),
34 fBuffSize (32000),
35 fVerbose (0),
36
37 fTreeK (nullptr),
38 fTreeH (nullptr),
39 fTreeC (nullptr),
40 fTreeR (nullptr),
41 fTreeKK (nullptr),
42 fTreeV0 (nullptr),
43 fTreeCC (nullptr),
44 fTreeGI (nullptr),
45
46 fK(), fpK (&fK),
47 fH(), fpH (&fH),
48 fC(), fpC (&fC),
49 fR(), fpR (&fR),
50 fKK(), fpKK(&fKK),
51 fV0(), fpV0(&fV0),
52 fCC(), fpCC(&fCC),
53 fGI(), fpGI(&fGI)
54{
55}
56
57////////////////////////////////////////////////////////////////////////////////
58/// Destructor.
59
61{
62}
63
64////////////////////////////////////////////////////////////////////////////////
65/// Set directory in which the trees are (or will be) created.
66
68{
69 fDirectory = dir;
70}
71
72////////////////////////////////////////////////////////////////////////////////
73/// Create internal trees.
74
76{
77 fDirectory->cd();
78 fTreeK = new TTree("Kinematics", "Simulated tracks.");
79 fTreeH = new TTree("Hits", "Combined detector hits.");
80 fTreeC = new TTree("Clusters", "Reconstructed clusters.");
81 fTreeR = new TTree("RecTracks", "Reconstructed tracks.");
82 fTreeKK = new TTree("RecKinks", "Reconstructed kinks.");
83 fTreeV0 = new TTree("RecV0s", "Reconstructed V0s.");
84 fTreeCC = new TTree("RecCascades","Reconstructed cascades.");
85 fTreeGI = new TTree("REveMCRecCrossRef", "Objects prepared for cross query.");
86}
87
88////////////////////////////////////////////////////////////////////////////////
89/// Delete internal trees.
90
92{
93 delete fTreeK; fTreeK = nullptr;
94 delete fTreeH; fTreeH = nullptr;
95 delete fTreeC; fTreeC = nullptr;
96 delete fTreeR; fTreeR = nullptr;
97 delete fTreeV0; fTreeV0 = nullptr;
98 delete fTreeKK; fTreeKK = nullptr;
99 delete fTreeGI; fTreeGI = nullptr;
100}
101
102////////////////////////////////////////////////////////////////////////////////
103/// Create internal VSD branches.
104
106{
107 if (fTreeK) fTreeK ->Branch("K", "REveMCTrack", &fpK);
108 if (fTreeH) fTreeH ->Branch("H", "REveHit", &fpH);
109 if (fTreeC) fTreeC ->Branch("C", "REveCluster", &fpC);
110 if (fTreeR) fTreeR ->Branch("R", "REveRecTrack", &fpR);
111 if (fTreeKK) fTreeKK->Branch("KK", "REveRecKink", &fpKK);
112 if (fTreeV0) fTreeV0->Branch("V0", "REveRecV0", &fpV0);
113
114 if (fTreeGI)
115 {
116 fTreeGI->Branch("GI", "REveMCRecCrossRef", &fpGI);
117 fTreeGI->Branch("K.", "REveMCTrack", &fpK);
118 fTreeGI->Branch("R.", "REveRecTrack", &fpR);
119 }
120}
121
122////////////////////////////////////////////////////////////////////////////////
123/// Set branche addresses of internal trees.
124
126{
127 if (fTreeK) fTreeK ->SetBranchAddress("K", &fpK);
128 if (fTreeH) fTreeH ->SetBranchAddress("H", &fpH);
129 if (fTreeC) fTreeC ->SetBranchAddress("C", &fpC);
130 if (fTreeR) fTreeR ->SetBranchAddress("R", &fpR);
131 if (fTreeKK) fTreeKK->SetBranchAddress("KK", &fpKK);
132 if (fTreeV0) fTreeV0->SetBranchAddress("V0", &fpV0);
133
134 if (fTreeGI)
135 {
139 }
140}
141
142////////////////////////////////////////////////////////////////////////////////
143/// Does nothing here ... reimplemented in sub-classes.
144
146{
147}
148
149////////////////////////////////////////////////////////////////////////////////
150/// Load internal trees from directory.
151
153{
154 static const REveException eH("REveVSD::LoadTrees");
155
156 if (!fDirectory)
157 throw eH + " directory not set.";
158
159 fTreeK = (TTree*) fDirectory->Get("Kinematics");
160 if (!fTreeK && fVerbose)
161 Error("REveVSD::LoadTrees","Kinematics not available in fDirectory %s.", fDirectory->GetName());
162
163 fTreeH = (TTree*) fDirectory->Get("Hits");
164 if (!fTreeH && fVerbose)
165 Error("REveVSD::LoadTrees", "Hits not available in fDirectory %s.", fDirectory->GetName());
166
167 fTreeC = (TTree*) fDirectory->Get("Clusters");
168 if (!fTreeC && fVerbose)
169 Error("REveVSD::LoadTrees", "Clusters not available in fDirectory %s.", fDirectory->GetName());
170
171 fTreeR = (TTree*) fDirectory->Get("RecTracks");
172 if (!fTreeR && fVerbose)
173 Error("REveVSD::LoadTrees", "RecTracks not available in fDirectory %s.", fDirectory->GetName());
174
175 fTreeKK = (TTree*) fDirectory->Get("RecKinks");
176 if (!fTreeKK && fVerbose)
177 Error("REveVSD::LoadTrees", "RecKinks not available in fDirectory %s.", fDirectory->GetName());
178
179 fTreeV0 = (TTree*) fDirectory->Get("RecV0s");
180 if (!fTreeV0 && fVerbose)
181 Error("REveVSD::LoadTrees", "RecV0 not available in fDirectory %s.", fDirectory->GetName());
182
183 fTreeGI = (TTree*)fDirectory->Get("REveMCRecCrossRef");
184 if(!fTreeGI && fVerbose)
185 Error("REveVSD::LoadTrees", "REveMCRecCrossRef not available in fDirectory %s.", fDirectory->GetName());
186}
187
188////////////////////////////////////////////////////////////////////////////////
189/// Disable TObject streamers for those VSD structs that inherit from
190/// TObject directly.
191
193{
194 // REveVector is not TObject
195
196 // REveMCTrack derives from TParticle
197 TParticle::Class()->IgnoreTObjectStreamer(true);
198
199 // REveHit::Class()->IgnoreTObjectStreamer(true);
200 // REveCluster::Class()->IgnoreTObjectStreamer(true);
201
202 // REveRecTrack::Class()->IgnoreTObjectStreamer(true);
203 // REveRecKink derives from REveRecTrack
204
205 // REveRecV0::Class()->IgnoreTObjectStreamer(true);
206
207 // REveMCRecCrossRef::Class()->IgnoreTObjectStreamer(true);
208}
void Class()
Definition: Class.C:29
REveException Exception-type thrown by Eve classes.
Definition: REveTypes.hxx:40
virtual void CreateBranches()
Create internal VSD branches.
Definition: REveVSD.cxx:105
TTree * fTreeR
Clusters.
Definition: REveVSD.hxx:40
TTree * fTreeH
Kinematics.
Definition: REveVSD.hxx:38
static void DisableTObjectStreamersForVSDStruct()
Disable TObject streamers for those VSD structs that inherit from TObject directly.
Definition: REveVSD.cxx:192
virtual void DeleteTrees()
Delete internal trees.
Definition: REveVSD.cxx:91
REveVSD(const REveVSD &)=delete
TTree * fTreeKK
Reconstructed tracks.
Definition: REveVSD.hxx:41
virtual void CreateTrees()
Create internal trees.
Definition: REveVSD.cxx:75
REveMCRecCrossRef * fpGI
Definition: REveVSD.hxx:53
virtual void WriteTrees()
Does nothing here ... reimplemented in sub-classes.
Definition: REveVSD.cxx:145
TTree * fTreeGI
Cascades.
Definition: REveVSD.hxx:44
virtual void SetDirectory(TDirectory *dir)
Set directory in which the trees are (or will be) created.
Definition: REveVSD.cxx:67
virtual ~REveVSD()
Destructor.
Definition: REveVSD.cxx:60
virtual void LoadTrees()
Load internal trees from directory.
Definition: REveVSD.cxx:152
virtual void SetBranchAddresses()
Set branche addresses of internal trees.
Definition: REveVSD.cxx:125
Describe directory structure in memory.
Definition: TDirectory.h:40
virtual TObject * Get(const char *namecycle)
Return pointer to object identified by namecycle.
Definition: TDirectory.cxx:806
virtual Bool_t cd(const char *path=nullptr)
Change current directory to "this" directory.
Definition: TDirectory.cxx:498
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Mother of all ROOT objects.
Definition: TObject.h:37
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:891
A TTree represents a columnar dataset.
Definition: TTree.h:78
virtual Int_t SetBranchAddress(const char *bname, void *add, TBranch **ptr=0)
Change branch address, dealing with clone trees properly.
Definition: TTree.cxx:8237
TBranch * Branch(const char *name, T *obj, Int_t bufsize=32000, Int_t splitlevel=99)
Add a new branch, and infer the data type from the type of obj being passed.
Definition: TTree.h:348