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

Public Member Functions

 ProtoDUNETruthBeamParticle (fhicl::ParameterSet const &pset)
 
virtual ~ProtoDUNETruthBeamParticle ()
 
void beginJob ()
 
void analyze (const art::Event &evt)
 
void reconfigure (fhicl::ParameterSet const &pset)
 
void endJob ()
 
- 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

int nbeamparticles
 
int run [NMAXBEAMPARTICLES]
 
int subrun [NMAXBEAMPARTICLES]
 
int event [NMAXBEAMPARTICLES]
 
int pdg [NMAXBEAMPARTICLES]
 
int goodparticle [NMAXBEAMPARTICLES]
 
float mom [NMAXBEAMPARTICLES]
 
float ene [NMAXBEAMPARTICLES]
 
float startpos [NMAXBEAMPARTICLES][4]
 
float p [NMAXBEAMPARTICLES][3]
 
TTree * truth
 
std::vector< int > p_pdg
 
bool p_savegoodparticle
 
bool selectall
 

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 41 of file ProtoDUNETruthBeamParticle_module.cc.

Constructor & Destructor Documentation

protoana::ProtoDUNETruthBeamParticle::ProtoDUNETruthBeamParticle ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 75 of file ProtoDUNETruthBeamParticle_module.cc.

75  : EDAnalyzer(pset){
76 
77  this->reconfigure(pset);
78 
79 }
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
void reconfigure(fhicl::ParameterSet const &pset)
protoana::ProtoDUNETruthBeamParticle::~ProtoDUNETruthBeamParticle ( )
virtual

Definition at line 82 of file ProtoDUNETruthBeamParticle_module.cc.

82 {}

Member Function Documentation

void protoana::ProtoDUNETruthBeamParticle::analyze ( const art::Event evt)

Definition at line 117 of file ProtoDUNETruthBeamParticle_module.cc.

117  {
118 
119  // Access truth info
120  std::vector< art::Ptr<simb::MCTruth> > mclist;
121  auto mctruthhandle = evt.getHandle< std::vector<simb::MCTruth> >("generator");
122  if (mctruthhandle){
123  art::fill_ptr_vector(mclist,mctruthhandle);
124  }
125  else{
126  mf::LogError("protoana::ProtoDUNETruthBeamParticle::analyze") << "Requested protoDUNE beam generator information does not exist!";
127  return;
128  }
129 
130  art::Ptr<simb::MCTruth> mctruth;
131  mctruth = mclist[0];
132 
133  nbeamparticles = 0;
134  for(int i = 0; i < mctruth->NParticles(); ++i){
135  const simb::MCParticle& part(mctruth->GetParticle(i));
136 
137  // Select partciles with the chosen pdg
138  if(!selectall){
139  bool found = false;
140  for(unsigned int j = 0; j < p_pdg.size(); j++){
141  if(part.PdgCode() == p_pdg[j]){
142  found = true;
143  break;
144  }
145  }
146  if(!found) continue;
147  }
148 
149  // Option to select only good beam particles
150  bool isgoodparticle = (part.Process() == "primary");
151  if(p_savegoodparticle && !isgoodparticle) continue;
152 
153  nbeamparticles++;
154 
155  run[i] = evt.run();
156  subrun[i] = evt.subRun();
157  event[i] = evt.id().event();
158  goodparticle[i] = (part.Process() == "primary");
159  startpos[i][0] = part.Vx();
160  startpos[i][1] = part.Vy();
161  startpos[i][2] = part.Vz();
162  startpos[i][3] = part.T();
163  p[i][0] = part.Px();
164  p[i][1] = part.Py();
165  p[i][2] = part.Pz();
166  mom[i] = part.P();
167  ene[i] = part.E();
168  pdg[i] = part.PdgCode();
169  }
170 
171  if(nbeamparticles > 0)
172  truth->Fill();
173 
174 }
Handle< PROD > getHandle(SelectorBase const &) const
Definition: DataViewImpl.h:382
int NParticles() const
Definition: MCTruth.h:75
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
SubRunNumber_t subRun() const
Definition: DataViewImpl.cc:78
RunNumber_t run() const
Definition: DataViewImpl.cc:71
const simb::MCParticle & GetParticle(int i) const
Definition: MCTruth.h:76
EventNumber_t event() const
Definition: EventID.h:116
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:297
EventID id() const
Definition: Event.cc:34
void protoana::ProtoDUNETruthBeamParticle::beginJob ( )
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 85 of file ProtoDUNETruthBeamParticle_module.cc.

85  {
86 
88 
89  // Define output tree
90  truth = tfs->make<TTree>("truth", "Beam Particle Truth Tree");
91  truth->Branch("nbeamparticles", &nbeamparticles, "nbeamparticles/I");
92  truth->Branch("run", run, "run[nbeamparticles]/I");
93  truth->Branch("subrun", subrun, "subrun[nbeamparticles]/I");
94  truth->Branch("event", event, "event[nbeamparticles]/I");
95  truth->Branch("pdg", pdg, "pdg[nbeamparticles]/I");
96  truth->Branch("goodparticle", goodparticle, "goodparticle[nbeamparticles]/I");
97  truth->Branch("mom", mom, "mom[nbeamparticles]/F");
98  truth->Branch("ene", ene, "ene[nbeamparticles]/F");
99  truth->Branch("startpos", startpos, "startpos[nbeamparticles][4]/F");
100  truth->Branch("p", p, "p[nbeamparticles][3]/F");
101 
102 }
Event finding and building.
void protoana::ProtoDUNETruthBeamParticle::endJob ( )
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 176 of file ProtoDUNETruthBeamParticle_module.cc.

176 {}
void protoana::ProtoDUNETruthBeamParticle::reconfigure ( fhicl::ParameterSet const &  pset)

Definition at line 105 of file ProtoDUNETruthBeamParticle_module.cc.

105  {
106 
107  p_savegoodparticle = pset.get<bool>("SaveOnlyGoodParticle");
108  p_pdg = pset.get< std::vector<int> >("Pdg");
109  if(p_pdg[0] == 0)
110  selectall = true;
111  else
112  selectall = false;
113 
114 }

Member Data Documentation

float protoana::ProtoDUNETruthBeamParticle::ene[NMAXBEAMPARTICLES]
private

Definition at line 61 of file ProtoDUNETruthBeamParticle_module.cc.

int protoana::ProtoDUNETruthBeamParticle::event[NMAXBEAMPARTICLES]
private

Definition at line 57 of file ProtoDUNETruthBeamParticle_module.cc.

int protoana::ProtoDUNETruthBeamParticle::goodparticle[NMAXBEAMPARTICLES]
private

Definition at line 59 of file ProtoDUNETruthBeamParticle_module.cc.

float protoana::ProtoDUNETruthBeamParticle::mom[NMAXBEAMPARTICLES]
private

Definition at line 60 of file ProtoDUNETruthBeamParticle_module.cc.

int protoana::ProtoDUNETruthBeamParticle::nbeamparticles
private

Definition at line 54 of file ProtoDUNETruthBeamParticle_module.cc.

float protoana::ProtoDUNETruthBeamParticle::p[NMAXBEAMPARTICLES][3]
private

Definition at line 63 of file ProtoDUNETruthBeamParticle_module.cc.

std::vector<int> protoana::ProtoDUNETruthBeamParticle::p_pdg
private

Definition at line 68 of file ProtoDUNETruthBeamParticle_module.cc.

bool protoana::ProtoDUNETruthBeamParticle::p_savegoodparticle
private

Definition at line 69 of file ProtoDUNETruthBeamParticle_module.cc.

int protoana::ProtoDUNETruthBeamParticle::pdg[NMAXBEAMPARTICLES]
private

Definition at line 58 of file ProtoDUNETruthBeamParticle_module.cc.

int protoana::ProtoDUNETruthBeamParticle::run[NMAXBEAMPARTICLES]
private

Definition at line 55 of file ProtoDUNETruthBeamParticle_module.cc.

bool protoana::ProtoDUNETruthBeamParticle::selectall
private

Definition at line 70 of file ProtoDUNETruthBeamParticle_module.cc.

float protoana::ProtoDUNETruthBeamParticle::startpos[NMAXBEAMPARTICLES][4]
private

Definition at line 62 of file ProtoDUNETruthBeamParticle_module.cc.

int protoana::ProtoDUNETruthBeamParticle::subrun[NMAXBEAMPARTICLES]
private

Definition at line 56 of file ProtoDUNETruthBeamParticle_module.cc.

TTree* protoana::ProtoDUNETruthBeamParticle::truth
private

Definition at line 65 of file ProtoDUNETruthBeamParticle_module.cc.


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