Public Member Functions | Private Attributes | List of all members
protoana::ProtoDUNETruthBeamFilter Class Reference
Inheritance diagram for protoana::ProtoDUNETruthBeamFilter:
art::EDFilter art::detail::Filter art::detail::LegacyModule art::Modifier art::ModuleBase art::ProductRegistryHelper

Public Member Functions

 ProtoDUNETruthBeamFilter (fhicl::ParameterSet const &pset)
 
virtual ~ProtoDUNETruthBeamFilter ()
 
void beginJob ()
 
bool filter (art::Event &evt)
 
void reconfigure (fhicl::ParameterSet const &pset)
 
void endJob ()
 
- Public Member Functions inherited from art::EDFilter
 EDFilter (fhicl::ParameterSet const &pset)
 
template<typename Config >
 EDFilter (Table< Config > const &config)
 
std::string workerType () const
 
- Public Member Functions inherited from art::detail::Filter
virtual ~Filter () noexcept
 
 Filter (fhicl::ParameterSet const &)
 
 Filter (Filter const &)=delete
 
 Filter (Filter &&)=delete
 
Filteroperator= (Filter const &)=delete
 
Filteroperator= (Filter &&)=delete
 
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 &)
 
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::Modifier
 ~Modifier () noexcept
 
 Modifier ()
 
 Modifier (Modifier const &)=delete
 
 Modifier (Modifier &&)=delete
 
Modifieroperator= (Modifier const &)=delete
 
Modifieroperator= (Modifier &&)=delete
 
- 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

bool fDebug
 
std::vector< int > fPDG
 Controls the verbosity. More...
 
TH1D * fSelectedEvents
 List of particle PDGs we want to keep. More...
 
TH1D * fTotalEvents
 

Additional Inherited Members

- Public Types inherited from art::EDFilter
using ModuleType = EDFilter
 
using WorkerType = WorkerT< EDFilter >
 
- Public Types inherited from art::detail::Filter
template<typename UserConfig >
using Table = Modifier::Table< UserConfig >
 
- Public Types inherited from art::Modifier
template<typename UserConfig , typename UserKeysToIgnore = void>
using Table = ProducerTable< UserConfig, detail::ModuleConfig, UserKeysToIgnore >
 
- Static Public Attributes inherited from art::detail::Filter
static constexpr bool Pass {true}
 
static constexpr bool Fail {false}
 
- 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 45 of file ProtoDUNETruthBeamFilter_module.cc.

Constructor & Destructor Documentation

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

Definition at line 72 of file ProtoDUNETruthBeamFilter_module.cc.

73 : EDFilter(pset)
74 {
75  this->reconfigure(pset);
76 
77 }
void reconfigure(fhicl::ParameterSet const &pset)
EDFilter(fhicl::ParameterSet const &pset)
Definition: EDFilter.h:21
protoana::ProtoDUNETruthBeamFilter::~ProtoDUNETruthBeamFilter ( )
virtual

Definition at line 80 of file ProtoDUNETruthBeamFilter_module.cc.

80 {}

Member Function Documentation

void protoana::ProtoDUNETruthBeamFilter::beginJob ( )
virtual

Reimplemented from art::EDFilter.

Definition at line 83 of file ProtoDUNETruthBeamFilter_module.cc.

83  {
84 
86  fSelectedEvents = tfs->make<TH1D>("fSelectedEvents", "Number of Selected Events", 3, 0, 3); //counts the number of selected events
87  fTotalEvents = tfs->make<TH1D>("fTotalEvents", "Total Events", 3, 0, 3); //counts the initial number of events in the unfiltered root input file
88 
89 }
TH1D * fSelectedEvents
List of particle PDGs we want to keep.
void protoana::ProtoDUNETruthBeamFilter::endJob ( )
virtual

Reimplemented from art::EDFilter.

Definition at line 147 of file ProtoDUNETruthBeamFilter_module.cc.

147 {}
bool protoana::ProtoDUNETruthBeamFilter::filter ( art::Event evt)
virtual

Implements art::EDFilter.

Definition at line 101 of file ProtoDUNETruthBeamFilter_module.cc.

101  {
102 
103 
104 
105  if(fDebug) std::cout << "Reading Event" << std::endl;
106 
107 
108  bool found = false;
109  // Access truth info
110  std::vector< art::Ptr<simb::MCTruth> > mclist;
111  auto mctruthhandle = evt.getHandle< std::vector<simb::MCTruth> >("generator");
112  if (mctruthhandle){
113  art::fill_ptr_vector(mclist,mctruthhandle);
114  }
115  else{
116  mf::LogError("protoana::ProtoDUNETruthBeamParticle::analyze") << "Requested protoDUNE beam generator information does not exist!";
117  return found;
118  }
119  fTotalEvents->Fill(1); //count total events
120 
121  art::Ptr<simb::MCTruth> mctruth;
122  mctruth = mclist[0];
123 
124 
125  for(int i = 0; i < mctruth->NParticles(); ++i){
126  const simb::MCParticle& part(mctruth->GetParticle(i));
127 
128  bool isgoodparticle = (part.Process() == "primary" && mctruth->Origin()==4); //primary particle and beam origin
129  if(!isgoodparticle) continue;
130  // Select partciles with the chosen pdg
131  for(unsigned int j = 0; j < fPDG.size(); j++){
132  if(part.PdgCode() == fPDG[j]){
133  found = true;
134  if(fDebug) std::cout << "this is a selected event with PDG "<< part.PdgCode() << std::endl;
135  fSelectedEvents->Fill(1); //count selected events
136  break;
137  }
138  }
139 
140  }
141 
142  return found;
143 
144 }
int PdgCode() const
Definition: MCParticle.h:212
Handle< PROD > getHandle(SelectorBase const &) const
Definition: DataViewImpl.h:382
simb::Origin_t Origin() const
Definition: MCTruth.h:74
int NParticles() const
Definition: MCTruth.h:75
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
const simb::MCParticle & GetParticle(int i) const
Definition: MCTruth.h:76
TH1D * fSelectedEvents
List of particle PDGs we want to keep.
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:297
std::vector< int > fPDG
Controls the verbosity.
QTextStream & endl(QTextStream &s)
void protoana::ProtoDUNETruthBeamFilter::reconfigure ( fhicl::ParameterSet const &  pset)

Definition at line 92 of file ProtoDUNETruthBeamFilter_module.cc.

92  {
93 
94 
95  fPDG = pset.get< std::vector<int> >("PDG");
96  fDebug = pset.get<bool>("IsVerbose");
97 
98 }
std::vector< int > fPDG
Controls the verbosity.

Member Data Documentation

bool protoana::ProtoDUNETruthBeamFilter::fDebug
private

Definition at line 58 of file ProtoDUNETruthBeamFilter_module.cc.

std::vector<int> protoana::ProtoDUNETruthBeamFilter::fPDG
private

Controls the verbosity.

Definition at line 60 of file ProtoDUNETruthBeamFilter_module.cc.

TH1D* protoana::ProtoDUNETruthBeamFilter::fSelectedEvents
private

List of particle PDGs we want to keep.

Definition at line 64 of file ProtoDUNETruthBeamFilter_module.cc.

TH1D* protoana::ProtoDUNETruthBeamFilter::fTotalEvents
private

Definition at line 65 of file ProtoDUNETruthBeamFilter_module.cc.


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