LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
MVAPID_module.cc
Go to the documentation of this file.
1 //
3 // \file MVAPID_module.cc
4 //
5 // m.haigh@warwick.ac.uk
6 //
8 
9 // ### Generic C++ includes ###
10 
11 // ### Framework includes ###
12 #include "TTree.h"
17 #include "art_root_io/TFileService.h"
18 
19 #include "MVAAlg.h"
23 
24 namespace mvapid {
25 
26  class MVAPID : public art::EDProducer {
27 
28  public:
29  explicit MVAPID(fhicl::ParameterSet const& pset);
30  void beginJob();
31  void produce(art::Event& evt);
32 
33  private:
35  std::vector<anab::MVAPIDResult>* fResult;
36  unsigned int fRun, fSubrun, fEvent;
37  TTree* fTree;
38 
39  }; // class MVAPID
40 
41  //------------------------------------------------------------------------------
42  MVAPID::MVAPID(fhicl::ParameterSet const& pset) : EDProducer{pset}, fAlg(pset)
43  {
44  produces<std::vector<anab::MVAPIDResult>>();
45  produces<art::Assns<recob::Track, anab::MVAPIDResult, void>>();
46  produces<art::Assns<recob::Shower, anab::MVAPIDResult, void>>();
47  fResult = new std::vector<anab::MVAPIDResult>;
48  }
49 
50  // ***************** //
52  {
54  fTree =
55  tfs->make<TTree>("MVAPID", "Results");
56  fTree->Branch("run", &fRun, "run/I");
57  fTree->Branch("subrun", &fSubrun, "subrun/I");
58  fTree->Branch("event", &fEvent, "event/I");
59  fTree->Branch("MVAResult", &fResult);
62  }
63 
64  // ***************** //
66  {
67  std::unique_ptr<std::vector<anab::MVAPIDResult>> result(new std::vector<anab::MVAPIDResult>);
68  std::unique_ptr<art::Assns<recob::Track, anab::MVAPIDResult>> trackAssns(
70  std::unique_ptr<art::Assns<recob::Shower, anab::MVAPIDResult>> showerAssns(
72  fRun = evt.id().run();
73  fSubrun = evt.id().subRun();
74  fEvent = evt.id().event();
75  fAlg.RunPID(evt, *result, *trackAssns, *showerAssns);
76  *fResult = *result;
77  fTree->Fill();
78  evt.put(std::move(result));
79  evt.put(std::move(trackAssns));
80  evt.put(std::move(showerAssns));
81  }
82 
84 
85 } //namespace mvapid
MVAPID(fhicl::ParameterSet const &pset)
unsigned int fSubrun
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.cc:6
void RunPID(art::Event &evt, std::vector< anab::MVAPIDResult > &result, art::Assns< recob::Track, anab::MVAPIDResult, void > &trackAssns, art::Assns< recob::Shower, anab::MVAPIDResult, void > &showerAssns)
Definition: MVAAlg.cxx:143
RunNumber_t run() const
Definition: EventID.h:98
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Event.h:77
Definition: MVAAlg.h:41
unsigned int fRun
void GetDetectorEdges()
Definition: MVAAlg.cxx:89
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
void GetWireNormals()
Definition: MVAAlg.cxx:111
Provides recob::Track data product.
unsigned int fEvent
EventNumber_t event() const
Definition: EventID.h:116
void produce(art::Event &evt)
TCEvent evt
Definition: DataStructs.cxx:8
std::vector< anab::MVAPIDResult > * fResult
SubRunNumber_t subRun() const
Definition: EventID.h:110
EventID id() const
Definition: Event.cc:23