Public Member Functions | Private Member Functions | Private Attributes | List of all members
rwgt::ReweightAna Class Reference

A module to check the results from the Monte Carlo generator. More...

Inheritance diagram for rwgt::ReweightAna:
art::EDAnalyzer art::detail::Analyzer art::detail::LegacyModule art::Observer art::ModuleBase

Public Member Functions

 ReweightAna (fhicl::ParameterSet const &pset)
 
virtual ~ReweightAna ()
 
void analyze (art::Event const &evt)
 
void beginSubRun (art::SubRun const &sr)
 
void beginJob ()
 
void endJob ()
 
void endSubRun (art::SubRun const &sr)
 
void reconfigure (const fhicl::ParameterSet &p)
 
- 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 ()
 
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 &)
 
std::string const & processName () const
 
bool wantAllEvents () const
 
bool wantEvent (Event const &e)
 
fhicl::ParameterSetID selectorConfig () const
 
Handle< TriggerResultsgetTriggerResults (Event const &e) 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 Member Functions

void LoadMCInfo (art::Event const &evt)
 

Private Attributes

TH1F * fEnergyNeutrino
 Total number of events. More...
 
TH1F * fNeventsSubrun
 Total number of events per subrun. More...
 
TH1F * fWgtQE [3]
 
TH1F * fWgtRES [3]
 
TH1F * fWgtDIS [3]
 
rwgt::NuReweightfGrwgt [3]
 X-sec weight calculator. More...
 
std::string fMCTruthModuleLabel
 label for module producing mc truth information More...
 
std::string fPotLabel
 Module that produced the POTSum object. More...
 

Additional Inherited Members

- Public Types inherited from art::EDAnalyzer
using WorkerType = WorkerT< EDAnalyzer >
 
using ModuleType = EDAnalyzer
 
- Protected Member Functions inherited from art::Observer
 Observer (fhicl::ParameterSet const &config)
 
 Observer (std::vector< std::string > const &paths, fhicl::ParameterSet const &config)
 
detail::ProcessAndEventSelectorsprocessAndEventSelectors ()
 
- 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

A module to check the results from the Monte Carlo generator.

Definition at line 50 of file ReweightAna_module.cc.

Constructor & Destructor Documentation

rwgt::ReweightAna::ReweightAna ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 85 of file ReweightAna_module.cc.

86  : EDAnalyzer(p)
87  {
88  this->reconfigure(p);
89  }
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:27
p
Definition: test.py:223
void reconfigure(const fhicl::ParameterSet &p)
rwgt::ReweightAna::~ReweightAna ( )
virtual

Definition at line 93 of file ReweightAna_module.cc.

93 { }

Member Function Documentation

void rwgt::ReweightAna::analyze ( art::Event const &  evt)
virtual

Implements art::EDAnalyzer.

Definition at line 138 of file ReweightAna_module.cc.

139  {
140 
141  // // Pull the MC generator information out of the event
142  mf::LogVerbatim("ReweightAna") << "Start analyze" ;
144  evt.getByLabel(fMCTruthModuleLabel, mclist);
145  if (mclist->empty()) {
146  mf::LogWarning("ReweightAna") << "Error retrieving MCTruth list" ;
147  return;
148  }
149 
151  evt.getByLabel(fMCTruthModuleLabel, gtlist);
152  if (gtlist->empty()) {
153  mf::LogWarning("ReweightAna") << "Error retrieving GTruth list" ;
154  return;
155  }
156 
157  MF_LOG_DEBUG("ReweightAna")<<"MC List sizes:" << mclist->size() << " " << gtlist->size() << "\n";
158 
159  // // Loop over neutrino interactions
160  for(size_t i_intx = 0; i_intx < mclist->size(); ++i_intx){
161  MF_LOG_DEBUG("ReweightAna") << "start loop";
162 
163  // // Link to the MCNeutrino class.
164  // // The class contains information not only about
165  // // the incoming neutrino, but about the products of the decay
166  simb::MCTruth const& truth = mclist->at(i_intx);
167  simb::GTruth const& gtruth = gtlist->at(i_intx);
168  simb::MCNeutrino const& mc_neutrino = truth.GetNeutrino();
169 
170  fEnergyNeutrino->Fill(mc_neutrino.Nu().E());
171  for(int i = 0; i < 3; i++) {
172  double wgt = fGrwgt[i]->CalcWeight(truth, gtruth);
173  //double wgt = 1.;
174  if(mc_neutrino.Mode()==0 && mc_neutrino.CCNC()==0) {
175  fWgtQE[i]->Fill(wgt);
176  }
177  else if(mc_neutrino.Mode()==1 && mc_neutrino.CCNC()==0) {
178  fWgtRES[i]->Fill(wgt);
179  }
180  else if(mc_neutrino.Mode()==2 && mc_neutrino.CCNC()==0) {
181  fWgtDIS[i]->Fill(wgt);
182  }
183  }
184 
185  MF_LOG_DEBUG("ReweightAna") << "end loop" ;
186  }//end loop over interactions
187 
188 
189  return;
190  }
double E(const int i=0) const
Definition: MCParticle.h:232
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
int CCNC() const
Definition: MCNeutrino.h:148
const simb::MCNeutrino & GetNeutrino() const
Definition: MCTruth.h:77
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:146
std::string fMCTruthModuleLabel
label for module producing mc truth information
TH1F * fEnergyNeutrino
Total number of events.
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
#define MF_LOG_DEBUG(id)
rwgt::NuReweight * fGrwgt[3]
X-sec weight calculator.
TCEvent evt
Definition: DataStructs.cxx:7
Event generator information.
Definition: MCTruth.h:32
Event generator information.
Definition: MCNeutrino.h:18
int Mode() const
Definition: MCNeutrino.h:149
double CalcWeight(const simb::MCTruth &truth, const simb::GTruth &gtruth) const
Definition: NuReweight.cxx:127
void rwgt::ReweightAna::beginJob ( )
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 96 of file ReweightAna_module.cc.

97  {
99 
100  mf::LogVerbatim("ReweightAna") << "make histograms" ;
101 
102  fEnergyNeutrino = tfs->make<TH1F>("fEnergyneutrino", "Total number of events", 50, 0., 25);
103  fNeventsSubrun = tfs->make<TH1F>("fNeventsSubrun", "Total number of events", 1, 0., 1.);
104 
105  char name[300];
106  for(int i = 0; i < 3; i++) {
107  sprintf(name, "fWgtQE_%dsigma", i+1);
108  fWgtQE[i] = tfs->make<TH1F>(name, "Evt Wgts", 100, 0., 2.0);
109  sprintf(name, "fWgtRES_%dsigma", i+1);
110  fWgtRES[i] = tfs->make<TH1F>(name, "Evt Wgts", 100, 0., 2.0);
111  sprintf(name, "fWgtDIS_%dsigma", i+1);
112  fWgtDIS[i] = tfs->make<TH1F>(name, "Evt Wgts", 100, 0., 2.0);
113 
114  double sigma = (double)(i+1);
115  fGrwgt[i] = new rwgt::NuReweight();
127 
128  fGrwgt[i]->Configure();
129  }
130  }
static QCString name
Definition: declinfo.cpp:673
void AddReweightValue(ReweightLabel_t rLabel, double value)
Change a reweight parameter. If it hasn&#39;t been added yet add it.
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
tweak the 2pi non-RES bkg in the RES region, for v+n CC
TH1F * fNeventsSubrun
Total number of events per subrun.
void Configure()
Reconfigure the weight calculators.
tweak the 2pi non-RES bkg in the RES region, for v+p CC
tweak the 2pi non-RES bkg in the RES region, for v+p NC
tweak the 2pi non-RES bkg in the RES region, for v+n NC
tweak Ma NCRES, affects d2sigma(NCRES)/dWdQ2 both in shape and normalization
tweak the 1pi non-RES bkg in the RES region, for v+n CC
tweak the 1pi non-RES bkg in the RES region, for v+p NC
TH1F * fEnergyNeutrino
Total number of events.
tweak Ma CCRES, affects d2sigma(CCRES)/dWdQ2 both in shape and normalization
tweak Ma CCQE, affects dsigma(CCQE)/dQ2 both in shape and normalization
rwgt::NuReweight * fGrwgt[3]
X-sec weight calculator.
tweak the 1pi non-RES bkg in the RES region, for v+p CC
tweak the 1pi non-RES bkg in the RES region, for v+n NC
void rwgt::ReweightAna::beginSubRun ( art::SubRun const &  sr)
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 133 of file ReweightAna_module.cc.

133  {
134 
135  }
void rwgt::ReweightAna::endJob ( )
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 216 of file ReweightAna_module.cc.

217  {
218 
219  }
void rwgt::ReweightAna::endSubRun ( art::SubRun const &  sr)
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 208 of file ReweightAna_module.cc.

209  {
210  fNeventsSubrun->Fill(sr.subRun(), cntEvent);
211  cntEvent = 0;
212 
213  }
TH1F * fNeventsSubrun
Total number of events per subrun.
static int cntEvent
static const double sr
Definition: Units.h:167
void rwgt::ReweightAna::LoadMCInfo ( art::Event const &  evt)
private

Definition at line 194 of file ReweightAna_module.cc.

195  {
196 
197 
198  }
void rwgt::ReweightAna::reconfigure ( const fhicl::ParameterSet p)

Definition at line 202 of file ReweightAna_module.cc.

203  {
204  fMCTruthModuleLabel = p.get< std::string>("MCTruthModuleLabel");
205  }
std::string string
Definition: nybbler.cc:12
std::string fMCTruthModuleLabel
label for module producing mc truth information
T get(std::string const &key) const
Definition: ParameterSet.h:231

Member Data Documentation

TH1F* rwgt::ReweightAna::fEnergyNeutrino
private

Total number of events.

Definition at line 67 of file ReweightAna_module.cc.

rwgt::NuReweight* rwgt::ReweightAna::fGrwgt[3]
private

X-sec weight calculator.

Definition at line 72 of file ReweightAna_module.cc.

std::string rwgt::ReweightAna::fMCTruthModuleLabel
private

label for module producing mc truth information

Definition at line 74 of file ReweightAna_module.cc.

TH1F* rwgt::ReweightAna::fNeventsSubrun
private

Total number of events per subrun.

Definition at line 68 of file ReweightAna_module.cc.

std::string rwgt::ReweightAna::fPotLabel
private

Module that produced the POTSum object.

Definition at line 75 of file ReweightAna_module.cc.

TH1F* rwgt::ReweightAna::fWgtDIS[3]
private

Definition at line 71 of file ReweightAna_module.cc.

TH1F* rwgt::ReweightAna::fWgtQE[3]
private

Definition at line 69 of file ReweightAna_module.cc.

TH1F* rwgt::ReweightAna::fWgtRES[3]
private

Definition at line 70 of file ReweightAna_module.cc.


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