LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
NuGraphAnalyzer_module.cc
Go to the documentation of this file.
1 // Class: NuGraphAnalyzer
3 // Plugin Type: analyzer (Unknown Unknown)
4 // File: NuGraphAnalyzer_module.cc
5 //
6 // Generated at Mon Nov 20 13:42:17 2023 by Giuseppe Cerati using cetskelgen
7 // from version .
9 
17 #include "fhiclcpp/ParameterSet.h"
19 
20 // saving output
21 #include "TTree.h"
22 #include "art_root_io/TFileService.h"
23 
28 
29 class NuGraphAnalyzer;
30 
31 using std::vector;
32 
34 public:
35  explicit NuGraphAnalyzer(fhicl::ParameterSet const& p);
36  // The compiler-generated destructor is fine for non-base
37  // classes without bare pointers or other resource use.
38 
39  // Plugins should not be copied or assigned.
40  NuGraphAnalyzer(NuGraphAnalyzer const&) = delete;
41  NuGraphAnalyzer(NuGraphAnalyzer&&) = delete;
42  NuGraphAnalyzer& operator=(NuGraphAnalyzer const&) = delete;
44 
45  // Required functions.
46  void analyze(art::Event const& e) override;
47 
48 private:
49  // Declare member data here.
50  TTree *_treeHit, *_treeEvt;
53  float _vtx_x, _vtx_y, _vtx_z;
54 };
55 
57 // More initializers here.
58 {
59  // Call appropriate consumes<>() for any products to be retrieved by this module.
61  _treeHit = tfs->make<TTree>("NuGraphHitOutput", "NuGraphHitOutput");
62  _treeHit->Branch("run", &_run, "run/I");
63  _treeHit->Branch("subrun", &_subrun, "subrun/I");
64  _treeHit->Branch("event", &_event, "event/I");
65  _treeHit->Branch("id", &_id, "id/I");
66  _treeHit->Branch("wire", &_wire, "wire/I");
67  _treeHit->Branch("plane", &_plane, "plane/I");
68  _treeHit->Branch("x_filter", &_x_filter, "x_filter/F");
69  _treeHit->Branch("MIP", &_MIP, "MIP/F");
70  _treeHit->Branch("HIP", &_HIP, "HIP/F");
71  _treeHit->Branch("shower", &_shower, "shower/F");
72  _treeHit->Branch("michel", &_michel, "michel/F");
73  _treeHit->Branch("diffuse", &_diffuse, "diffuse/F");
74  _treeHit->Branch("time", &_time, "time/F");
75  _treeEvt = tfs->make<TTree>("NuGraphEventOutput", "NuGraphEventOutput");
76  _treeEvt->Branch("run", &_run, "run/I");
77  _treeEvt->Branch("subrun", &_subrun, "subrun/I");
78  _treeEvt->Branch("event", &_event, "event/I");
79  _treeEvt->Branch("vtx_x", &_vtx_x, "vtx_x/F");
80  _treeEvt->Branch("vtx_y", &_vtx_y, "vtx_y/F");
81  _treeEvt->Branch("vtx_z", &_vtx_z, "vtx_z/F");
82 }
83 
85 {
86 
87  auto GNNDescription = e.getHandle<anab::MVADescription<5>>(art::InputTag("NuGraph", "semantic"));
88 
89  auto const& hitsWithScores = proxy::getCollection<std::vector<recob::Hit>>(
90  e,
91  GNNDescription->dataTag(), //tag of the hit collection we ran the GNN on
92  proxy::withParallelData<anab::FeatureVector<1>>(art::InputTag("NuGraph", "filter")),
94 
95  std::cout << hitsWithScores.size() << std::endl;
96  for (auto& h : hitsWithScores) {
97  const auto& assocFilter = h.get<anab::FeatureVector<1>>();
98  const auto& assocSemantic = h.get<anab::FeatureVector<5>>();
99  _event = e.event();
100  _subrun = e.subRun();
101  _run = e.run();
102  _id = h.index();
103  _x_filter = assocFilter.at(0);
104  _MIP = assocSemantic.at(GNNDescription->getIndex("MIP"));
105  _HIP = assocSemantic.at(GNNDescription->getIndex("HIP"));
106  _shower = assocSemantic.at(GNNDescription->getIndex("shower"));
107  _michel = assocSemantic.at(GNNDescription->getIndex("michel"));
108  _diffuse = assocSemantic.at(GNNDescription->getIndex("diffuse"));
109  _treeHit->Fill();
110  }
111 
112  auto PredVertexColl = e.getHandle<std::vector<recob::Vertex>>(art::InputTag("NuGraph", "vertex"));
113  if (PredVertexColl.isValid() && PredVertexColl->size() > 0) { //there should be only one
114  _vtx_x = PredVertexColl->at(0).position().X();
115  _vtx_y = PredVertexColl->at(0).position().Y();
116  _vtx_z = PredVertexColl->at(0).position().Z();
117  _treeEvt->Fill();
118  }
119 }
120 
NuGraphAnalyzer(fhicl::ParameterSet const &p)
Declaration of signal hit object.
Base utilities for the implementation of data product facades.
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.cc:6
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
auto withParallelData(Args &&...args)
Helper function to merge an auxiliary data product into the proxy.
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
void analyze(art::Event const &e) override
Handle< PROD > getHandle(SelectorBase const &) const
float at(size_t index) const
Definition: MVAOutput.h:94
NuGraphAnalyzer & operator=(NuGraphAnalyzer const &)=delete
Float_t e
Definition: plot.C:35
size_t size() const
Definition: MVAOutput.h:199