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