Public Member Functions | Private Attributes | List of all members
cvn::GCNH5 Class Reference
Inheritance diagram for cvn::GCNH5:
art::EDAnalyzer art::detail::Analyzer art::detail::LegacyModule art::Observer art::ModuleBase

Public Member Functions

 GCNH5 (fhicl::ParameterSet const &p)
 
 ~GCNH5 () noexcept
 
 GCNH5 (GCNH5 const &)=delete
 
 GCNH5 (GCNH5 &&)=delete
 
GCNH5operator= (GCNH5 const &)=delete
 
GCNH5operator= (GCNH5 &&)=delete
 
void reconfigure (fhicl::ParameterSet const &p)
 
void beginSubRun (art::SubRun const &sr) override
 
void endSubRun (art::SubRun const &sr) override
 
void analyze (art::Event const &e) override
 
- Public Member Functions inherited from art::EDAnalyzer
 EDAnalyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 EDAnalyzer (Table< Config > const &config)
 
std::string workerType () const
 
- Public Member Functions inherited from art::detail::Analyzer
virtual ~Analyzer () noexcept
 
 Analyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 Analyzer (Table< Config > const &config)
 
void doBeginJob (SharedResources const &resources)
 
void doEndJob ()
 
void doRespondToOpenInputFile (FileBlock const &fb)
 
void doRespondToCloseInputFile (FileBlock const &fb)
 
void doRespondToOpenOutputFiles (FileBlock const &fb)
 
void doRespondToCloseOutputFiles (FileBlock const &fb)
 
bool doBeginRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doEndRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doBeginSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEndSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEvent (EventPrincipal &ep, ModuleContext const &mc, std::atomic< std::size_t > &counts_run, std::atomic< std::size_t > &counts_passed, std::atomic< std::size_t > &counts_failed)
 
- Public Member Functions inherited from art::Observer
 ~Observer () noexcept
 
 Observer (Observer const &)=delete
 
 Observer (Observer &&)=delete
 
Observeroperator= (Observer const &)=delete
 
Observeroperator= (Observer &&)=delete
 
void registerProducts (ProductDescriptions &, ModuleDescription const &)
 
void fillDescriptions (ModuleDescription const &)
 
fhicl::ParameterSetID selectorConfig () const
 
- Public Member Functions inherited from art::ModuleBase
virtual ~ModuleBase () noexcept
 
 ModuleBase ()
 
ModuleDescription const & moduleDescription () const
 
void setModuleDescription (ModuleDescription const &)
 
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables () const
 
void sortConsumables (std::string const &current_process_name)
 
template<typename T , BranchType BT>
ViewToken< T > consumesView (InputTag const &tag)
 
template<typename T , BranchType BT>
ViewToken< T > mayConsumeView (InputTag const &tag)
 

Private Attributes

string fGraphModuleLabel
 Name of graph producer module. More...
 
string fGraphInstanceLabel
 Name of graph instance. More...
 
string fTruthLabel
 Name of truth producer module. More...
 
string fOutputName
 Output filename. More...
 
bool fSaveEventTruth
 Whether to save event-level truth information. More...
 
bool fSaveParticleTruth
 Whether to include particle truth information. More...
 
hep_hpc::hdf5::File fFile
 Output HDF5 file. More...
 
hep_hpc::hdf5::Ntuple< Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< int, 1 >, Column< int, 1 > > * fGraphNtuple
 graph ntuple More...
 
hep_hpc::hdf5::Ntuple< Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 > > * fEventNtuple
 Event ntuple. More...
 
hep_hpc::hdf5::Ntuple< Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< string, 1 >, Column< string, 1 > > * fParticleNtuple
 Particle ntuple. More...
 

Additional Inherited Members

- Public Types inherited from art::EDAnalyzer
using WorkerType = WorkerT< EDAnalyzer >
 
using ModuleType = EDAnalyzer
 
- Protected Member Functions inherited from art::Observer
std::string const & processName () const
 
bool wantAllEvents () const noexcept
 
bool wantEvent (ScheduleID id, Event const &e) const
 
Handle< TriggerResultsgetTriggerResults (Event const &e) const
 
 Observer (fhicl::ParameterSet const &config)
 
 Observer (std::vector< std::string > const &select_paths, std::vector< std::string > const &reject_paths, fhicl::ParameterSet const &config)
 
- Protected Member Functions inherited from art::ModuleBase
ConsumesCollectorconsumesCollector ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 

Detailed Description

Definition at line 42 of file GCNH5_module.cc.

Constructor & Destructor Documentation

cvn::GCNH5::GCNH5 ( fhicl::ParameterSet const &  p)
explicit

Definition at line 110 of file GCNH5_module.cc.

111  : EDAnalyzer{p} {
112  this->reconfigure(p);
113  }
void reconfigure(fhicl::ParameterSet const &p)
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
p
Definition: test.py:223
cvn::GCNH5::~GCNH5 ( )
inlinenoexcept

Definition at line 45 of file GCNH5_module.cc.

45 {};
cvn::GCNH5::GCNH5 ( GCNH5 const &  )
delete
cvn::GCNH5::GCNH5 ( GCNH5 &&  )
delete

Member Function Documentation

void cvn::GCNH5::analyze ( art::Event const &  e)
overridevirtual

Implements art::EDAnalyzer.

Definition at line 126 of file GCNH5_module.cc.

126  {
127 
128  // Get the graphVector
129  std::vector<art::Ptr<GCNGraph>> graphVector;
131  auto graphHandle = e.getHandle<std::vector<GCNGraph>>(itag1);
132  if (!graphHandle) {
134  << "Could not find GCNGraph vector with module label "
135  << fGraphModuleLabel << " and instance label "
136  << fGraphInstanceLabel << "!" << endl;
137  }
138  art::fill_ptr_vector(graphVector, graphHandle);
139 
140  if (graphVector.size() > 1) throw art::Exception(art::errors::LogicError)
141  << "There shouldn't be more than one GCNGraph per producer per event,"
142  << " but here there are " << graphVector.size() << "." << endl;
143 
144  if (graphVector.empty()) return;
145 
146  int run = e.id().run();
147  int subrun = e.id().subRun();
148  int event = e.id().event();
149 
150  for (size_t itNode = 0; itNode < graphVector[0]->GetNumberOfNodes(); ++itNode) {
151  GCNGraphNode node = graphVector[0]->GetNode(itNode);
152 
153  vector<float> pos = node.GetPosition();
154  vector<float> feat = node.GetFeatures();
155  vector<float> truth = node.GetGroundTruth();
156 
157  fGraphNtuple->insert(run, subrun, event, (int)round(abs(pos[0])),
158  pos[1], pos[2], feat[0], feat[1], feat[4], feat[2], feat[3],
159  (int)feat[5], truth[0]);
160  }
161 
162  // Event truth
163  if (fSaveEventTruth) {
164 
165  // Get MC truth
166  auto truthHandle = e.getHandle<vector<simb::MCTruth>>(fTruthLabel);
167  if (!truthHandle || truthHandle->size() != 1) {
169  << "Expected to find exactly one MC truth object!";
170  }
171  simb::MCNeutrino nutruth = truthHandle->at(0).GetNeutrino();
172 
173  // Fill variables
174  fEventNtuple->insert(run, subrun, event,
175  nutruth.CCNC() == simb::kCC,
176  nutruth.Nu().E(),
177  nutruth.Lepton().E(),
178  nutruth.Nu().Momentum().Vect().Unit().X(),
179  nutruth.Nu().Momentum().Vect().Unit().Y(),
180  nutruth.Nu().Momentum().Vect().Unit().Z());
181 
182  } // if saving event truth
183 
184  // // Particle tree
185  if (fSaveParticleTruth) {
186  std::vector<ptruth> ptree = GCNFeatureUtils::GetParticleTree(graphVector[0].get());
187  for (auto p : ptree) {
188  fParticleNtuple->insert(run, subrun, event,
189  std::get<0>(p), std::get<1>(p), std::get<2>(p),
190  std::get<3>(p), std::get<4>(p), std::get<5>(p),
191  std::get<6>(p), std::get<7>(p), std::get<8>(p),
192  std::get<9>(p), std::get<10>(p), std::get<11>(p));
193  }
194  }
195 
196  } // cvn::GCNH5::analyze
double E(const int i=0) const
Definition: MCParticle.h:233
int CCNC() const
Definition: MCNeutrino.h:148
static std::vector< ptruth > GetParticleTree(const cvn::GCNGraph *g)
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:146
hep_hpc::hdf5::Ntuple< Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< string, 1 >, Column< string, 1 > > * fParticleNtuple
Particle ntuple.
bool fSaveEventTruth
Whether to save event-level truth information.
Definition: GCNH5_module.cc:64
hep_hpc::hdf5::Ntuple< Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 > > * fEventNtuple
Event ntuple.
Definition: GCNH5_module.cc:90
T abs(T value)
hep_hpc::hdf5::Ntuple< Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< int, 1 >, Column< int, 1 > > * fGraphNtuple
graph ntuple
Definition: GCNH5_module.cc:80
const double e
const simb::MCParticle & Lepton() const
Definition: MCNeutrino.h:147
bool fSaveParticleTruth
Whether to include particle truth information.
Definition: GCNH5_module.cc:65
p
Definition: test.py:223
string fTruthLabel
Name of truth producer module.
Definition: GCNH5_module.cc:62
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:220
string fGraphModuleLabel
Name of graph producer module.
Definition: GCNH5_module.cc:60
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:297
Event generator information.
Definition: MCNeutrino.h:18
string fGraphInstanceLabel
Name of graph instance.
Definition: GCNH5_module.cc:61
QTextStream & endl(QTextStream &s)
Event finding and building.
void cvn::GCNH5::beginSubRun ( art::SubRun const &  sr)
overridevirtual

Reimplemented from art::EDAnalyzer.

Definition at line 198 of file GCNH5_module.cc.

198  {
199 
200  // Open HDF5 output
201  boost::uuids::random_generator generator;
202  boost::uuids::uuid uuid = generator();
203  std::ostringstream fileName;
204  fileName << fOutputName << "_r" << setfill('0') << setw(5) << sr.run()
205  << "_r" << setfill('0') << setw(5) << sr.subRun() << "_" << uuid
206  << ".h5";
207 
208  fFile = hep_hpc::hdf5::File(fileName.str(), H5F_ACC_TRUNC);
209 
210  fGraphNtuple = new hep_hpc::hdf5::Ntuple(
211  make_ntuple({fFile, "graph_table", 1000},
212  make_scalar_column<int>("run"),
213  make_scalar_column<int>("subrun"),
214  make_scalar_column<int>("event"),
215  make_scalar_column<int>("plane"),
216  make_scalar_column<float>("wire"),
217  make_scalar_column<float>("time"),
218  make_scalar_column<float>("integral"),
219  make_scalar_column<float>("rms"),
220  make_scalar_column<int>("rawplane"),
221  make_scalar_column<float>("rawwire"),
222  make_scalar_column<float>("rawtime"),
223  make_scalar_column<int>("tpc"),
224  make_scalar_column<int>("true_id")));
225 
226  if (fSaveEventTruth) {
227  fEventNtuple = new hep_hpc::hdf5::Ntuple(
228  make_ntuple({fFile, "event_table", 1000},
229  make_scalar_column<int>("run"),
230  make_scalar_column<int>("subrun"),
231  make_scalar_column<int>("event"),
232  make_scalar_column<int>("is_cc"),
233  make_scalar_column<float>("nu_energy"),
234  make_scalar_column<float>("lep_energy"),
235  make_scalar_column<float>("nu_dir_x"),
236  make_scalar_column<float>("nu_dir_y"),
237  make_scalar_column<float>("nu_dir_z")));
238  }
239 
240  if (fSaveParticleTruth) {
241  fParticleNtuple = new hep_hpc::hdf5::Ntuple(
242  make_ntuple({fFile, "particle_table", 1000},
243  make_scalar_column<int>("run"),
244  make_scalar_column<int>("subrun"),
245  make_scalar_column<int>("event"),
246  make_scalar_column<int>("id"),
247  make_scalar_column<int>("type"),
248  make_scalar_column<int>("parent_id"),
249  make_scalar_column<float>("momentum"),
250  make_scalar_column<float>("start_x"),
251  make_scalar_column<float>("start_y"),
252  make_scalar_column<float>("start_z"),
253  make_scalar_column<float>("end_x"),
254  make_scalar_column<float>("end_y"),
255  make_scalar_column<float>("end_z"),
256  make_scalar_column<string>("start_process"),
257  make_scalar_column<string>("end_process")));
258  }
259  }
hep_hpc::hdf5::Ntuple< Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< string, 1 >, Column< string, 1 > > * fParticleNtuple
Particle ntuple.
bool fSaveEventTruth
Whether to save event-level truth information.
Definition: GCNH5_module.cc:64
hep_hpc::hdf5::File fFile
Output HDF5 file.
Definition: GCNH5_module.cc:67
hep_hpc::hdf5::Ntuple< Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 > > * fEventNtuple
Event ntuple.
Definition: GCNH5_module.cc:90
hep_hpc::hdf5::Ntuple< Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< int, 1 >, Column< int, 1 > > * fGraphNtuple
graph ntuple
Definition: GCNH5_module.cc:80
fileName
Definition: dumpTree.py:9
bool fSaveParticleTruth
Whether to include particle truth information.
Definition: GCNH5_module.cc:65
generator
Definition: train.py:468
Q_EXPORT QTSManip setw(int w)
Definition: qtextstream.h:331
string fOutputName
Output filename.
Definition: GCNH5_module.cc:63
static constexpr double sr
Definition: Units.h:166
Q_EXPORT QTSManip setfill(int f)
Definition: qtextstream.h:337
void cvn::GCNH5::endSubRun ( art::SubRun const &  sr)
overridevirtual

Reimplemented from art::EDAnalyzer.

Definition at line 261 of file GCNH5_module.cc.

261  {
262  delete fGraphNtuple;
263  if (fSaveEventTruth) delete fEventNtuple;
265  fFile.close();
266  }
hep_hpc::hdf5::Ntuple< Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< string, 1 >, Column< string, 1 > > * fParticleNtuple
Particle ntuple.
bool fSaveEventTruth
Whether to save event-level truth information.
Definition: GCNH5_module.cc:64
hep_hpc::hdf5::File fFile
Output HDF5 file.
Definition: GCNH5_module.cc:67
hep_hpc::hdf5::Ntuple< Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 > > * fEventNtuple
Event ntuple.
Definition: GCNH5_module.cc:90
hep_hpc::hdf5::Ntuple< Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< float, 1 >, Column< int, 1 >, Column< float, 1 >, Column< float, 1 >, Column< int, 1 >, Column< int, 1 > > * fGraphNtuple
graph ntuple
Definition: GCNH5_module.cc:80
bool fSaveParticleTruth
Whether to include particle truth information.
Definition: GCNH5_module.cc:65
GCNH5& cvn::GCNH5::operator= ( GCNH5 const &  )
delete
GCNH5& cvn::GCNH5::operator= ( GCNH5 &&  )
delete
void cvn::GCNH5::reconfigure ( fhicl::ParameterSet const &  p)

Definition at line 115 of file GCNH5_module.cc.

115  {
116 
117  fGraphModuleLabel = p.get<string>("GraphModuleLabel");
118  fGraphInstanceLabel = p.get<string>("GraphInstanceLabel");
119  fTruthLabel = p.get<string>("TruthLabel"),
120  fOutputName = p.get<string>("OutputName");
121  fSaveEventTruth = p.get<bool>("SaveEventTruth");
122  fSaveParticleTruth = p.get<bool>("SaveParticleTruth");
123 
124  } // cvn::GCNH5::reconfigure
bool fSaveEventTruth
Whether to save event-level truth information.
Definition: GCNH5_module.cc:64
bool fSaveParticleTruth
Whether to include particle truth information.
Definition: GCNH5_module.cc:65
p
Definition: test.py:223
string fTruthLabel
Name of truth producer module.
Definition: GCNH5_module.cc:62
string fOutputName
Output filename.
Definition: GCNH5_module.cc:63
string fGraphModuleLabel
Name of graph producer module.
Definition: GCNH5_module.cc:60
string fGraphInstanceLabel
Name of graph instance.
Definition: GCNH5_module.cc:61

Member Data Documentation

hep_hpc::hdf5::Ntuple<Column<int, 1>, Column<int, 1>, Column<int, 1>, Column<int, 1>, Column<float, 1>, Column<float, 1>, Column<float, 1>, Column<float, 1>, Column<float, 1> >* cvn::GCNH5::fEventNtuple
private

Event ntuple.

Definition at line 90 of file GCNH5_module.cc.

hep_hpc::hdf5::File cvn::GCNH5::fFile
private

Output HDF5 file.

Definition at line 67 of file GCNH5_module.cc.

string cvn::GCNH5::fGraphInstanceLabel
private

Name of graph instance.

Definition at line 61 of file GCNH5_module.cc.

string cvn::GCNH5::fGraphModuleLabel
private

Name of graph producer module.

Definition at line 60 of file GCNH5_module.cc.

hep_hpc::hdf5::Ntuple<Column<int, 1>, Column<int, 1>, Column<int, 1>, Column<int, 1>, Column<float, 1>, Column<float, 1>, Column<float, 1>, Column<float, 1>, Column<int, 1>, Column<float, 1>, Column<float, 1>, Column<int, 1>, Column<int, 1> >* cvn::GCNH5::fGraphNtuple
private

graph ntuple

Definition at line 80 of file GCNH5_module.cc.

string cvn::GCNH5::fOutputName
private

Output filename.

Definition at line 63 of file GCNH5_module.cc.

hep_hpc::hdf5::Ntuple<Column<int, 1>, Column<int, 1>, Column<int, 1>, Column<int, 1>, Column<int, 1>, Column<int, 1>, Column<float, 1>, Column<float, 1>, Column<float, 1>, Column<float, 1>, Column<float, 1>, Column<float, 1>, Column<float, 1>, Column<string, 1>, Column<string, 1> >* cvn::GCNH5::fParticleNtuple
private

Particle ntuple.

Definition at line 106 of file GCNH5_module.cc.

bool cvn::GCNH5::fSaveEventTruth
private

Whether to save event-level truth information.

Definition at line 64 of file GCNH5_module.cc.

bool cvn::GCNH5::fSaveParticleTruth
private

Whether to include particle truth information.

Definition at line 65 of file GCNH5_module.cc.

string cvn::GCNH5::fTruthLabel
private

Name of truth producer module.

Definition at line 62 of file GCNH5_module.cc.


The documentation for this class was generated from the following file: