Public Member Functions | Private Attributes | Friends | List of all members
simb::MCTruth Class Reference

Event generator information. More...

#include "/cvmfs/nova.opensciencegrid.org/externals/nusimdata/v1_16_04/source/nusimdata/SimulationBase/MCTruth.h"

Public Member Functions

 MCTruth ()
 
const simb::MCGeneratorInfoGeneratorInfo () const
 
simb::Origin_t Origin () const
 
int NParticles () const
 
const simb::MCParticleGetParticle (int i) const
 
const simb::MCNeutrinoGetNeutrino () const
 
bool NeutrinoSet () const
 
void Add (simb::MCParticle &part)
 
void SetGeneratorInfo (simb::Generator_t generator, const std::string &genVersion, const std::unordered_map< std::string, std::string > &genConfig)
 
void SetOrigin (simb::Origin_t origin)
 
void SetNeutrino (int CCNC, int mode, int interactionType, int target, int nucleon, int quark, double w, double x, double y, double qsqr)
 

Private Attributes

std::vector< simb::MCParticlefPartList
 list of particles in this event More...
 
simb::MCNeutrino fMCNeutrino
 reference to neutrino info - null if not a neutrino More...
 
simb::Origin_t fOrigin
 origin for this event More...
 
simb::MCGeneratorInfo fGenInfo
 information about the generator that produced this event More...
 
bool fNeutrinoSet
 flag for whether the neutrino information has been set More...
 

Friends

std::ostream & operator<< (std::ostream &o, simb::MCTruth const &a)
 

Detailed Description

Event generator information.

Definition at line 32 of file MCTruth.h.

Constructor & Destructor Documentation

simb::MCTruth::MCTruth ( )

Definition at line 21 of file MCTruth.cxx.

22  : fPartList()
23  , fMCNeutrino()
25  , fNeutrinoSet(false)
26  {
27  }
std::vector< simb::MCParticle > fPartList
list of particles in this event
Definition: MCTruth.h:38
simb::MCNeutrino fMCNeutrino
reference to neutrino info - null if not a neutrino
Definition: MCTruth.h:39
simb::Origin_t fOrigin
origin for this event
Definition: MCTruth.h:40
bool fNeutrinoSet
flag for whether the neutrino information has been set
Definition: MCTruth.h:42

Member Function Documentation

void simb::MCTruth::Add ( simb::MCParticle part)
inline
const simb::MCGeneratorInfo & simb::MCTruth::GeneratorInfo ( ) const
inline

Definition at line 72 of file MCTruth.h.

References fGenInfo.

Referenced by caf::CAFMaker::AddMCTruthToVec(), and novarwgt::ConvertNuToolsEvent().

72 { return fGenInfo; }
simb::MCGeneratorInfo fGenInfo
information about the generator that produced this event
Definition: MCTruth.h:41
const simb::MCNeutrino & simb::MCTruth::GetNeutrino ( ) const
inline

Definition at line 76 of file MCTruth.h.

References fMCNeutrino.

Referenced by caf::CAFMaker::AddMCTruthToVec(), cheat::TestTrackIds::analyze(), tut::RecoValidationTutorial::analyze(), G4MismatchAna::analyze(), mcchk::RockAna::analyze(), mcchk::NeutrinoAna::analyze(), showere::ShowerEnergyAna::analyze(), earms::ElasticArmsValidate::analyze(), fuzz::FuzzyKValidate::analyze(), bpf::BPFCVNAna::analyze(), ncs::Xbeam::analyze(), ncs::GenieTruth::analyze(), ncs::Xeff::analyze(), slid::LIDTraining::analyze(), ncs::NCAna::analyze(), murem::MRCCAna::analyze(), mcchk::NeutrinoAna::AnalyzeNeutrinoInteraction(), murem::MRCCNeutrino::CC(), trk::KalmanTrackAna::CheckRecoTracks(), novarwgt::ConvertNuToolsEvent(), gibuu::GiBUURegen::CopyGenieEvent(), murem::MRCCNeutrino::Energy(), filter::Filter::FillTruthVariables(), fnex::AnalysisSetupBase::FillTruthVars(), showere::ShowerEnergyFilterMC::filter(), fnex::AnalysisSetupBase::FindDISWeight(), fnex::AnalysisSetupBase::FindMECWeight(), fnex::AnalysisSetupBase::FindRPAWeights(), nuesand::FillNueSandbox::GetECF(), gibuu::GiBUURegen::GetEvent(), cvn::GetInteractionType(), gibuu::GiBUURegen::GetKey(), cvn::GetMultiplicityMap(), rwgt::MCReweight::GetNOvARwgtEvRec(), cvn::GetParentParticleType(), gibuu::GiBUURegen::GetRecordList(), rwgt::MCReweight::GetWeight(), evg::MCTruthToDk2NuHackItr::match_mctruth_nuchoice(), earms::ElasticArmsHS::MCVertex(), murem::MRCCNeutrino::MRCCNeutrino(), simb::operator<<(), cvn::ParticlesSliceClassify(), MergeGenCollections::MergeGenCollections::produce(), evgen::GENIERockGen::produce(), ncid::NCNNKeras::produce(), fuzz::FuzzyKVertex::produce(), gibuu::GiBUURegen::produce(), caf::CAFMaker::produce(), MergeG4Collections::MergeG4Collections::SameMCTruth(), and cvn::SliceClassify().

76 { return fMCNeutrino; }
simb::MCNeutrino fMCNeutrino
reference to neutrino info - null if not a neutrino
Definition: MCTruth.h:39
const simb::MCParticle & simb::MCTruth::GetParticle ( int  i) const
inline
bool simb::MCTruth::NeutrinoSet ( ) const
inline
int simb::MCTruth::NParticles ( ) const
inline
simb::Origin_t simb::MCTruth::Origin ( ) const
inline

Definition at line 73 of file MCTruth.h.

References fOrigin.

Referenced by simb::operator<<(), filter::TruthFilter::produce(), MergeGenCollections::MergeGenCollections::produce(), ifdb::IFDBSpillInfo::produce(), and ifdb::MIN::produce().

73 { return fOrigin; }
simb::Origin_t fOrigin
origin for this event
Definition: MCTruth.h:40
void simb::MCTruth::SetGeneratorInfo ( simb::Generator_t  generator,
const std::string &  genVersion,
const std::unordered_map< std::string, std::string > &  genConfig 
)
inline

Definition at line 82 of file MCTruth.h.

References fGenInfo.

Referenced by evgb::CRYHelper::Sample().

85 {
86  fGenInfo = simb::MCGeneratorInfo(generator, genVersion, genConfig);
87 }
simb::MCGeneratorInfo fGenInfo
information about the generator that produced this event
Definition: MCTruth.h:41
void simb::MCTruth::SetNeutrino ( int  CCNC,
int  mode,
int  interactionType,
int  target,
int  nucleon,
int  quark,
double  w,
double  x,
double  y,
double  qsqr 
)

Definition at line 30 of file MCTruth.cxx.

References abs(), fMCNeutrino, fNeutrinoSet, fPartList, MECModelEnuComparisons::i, art::errors::LogicError, genie::utils::res::PdgCode(), simb::MCParticle::PdgCode(), and simb::MCParticle::TrackId().

Referenced by gibuu::GiBUURegen::CopyGenieEvent(), gibuu::GiBUURegen::GetEvent(), supernova::SnovaGen::MakeMCTruth(), and MergeGenCollections::MergeGenCollections::produce().

40  {
41  if( !fNeutrinoSet ){
42  fNeutrinoSet = true;
43  // loop over the MCParticle list and get the outgoing lepton
44  // assume this is a neutral current interaction to begin with
45  // which means the outgoing lepton is the incoming neutrino
46  MCParticle nu = fPartList[0];
47  MCParticle lep = fPartList[0];
48 
49  // start at i = 1 because i = 0 is the incoming neutrino
50  for(unsigned int i = 1; i < fPartList.size(); ++i){
51  if(fPartList[i].Mother() == nu.TrackId() &&
52  (fPartList[i].PdgCode() == nu.PdgCode() ||
53  abs(fPartList[i].PdgCode()) == abs(nu.PdgCode())-1) ){
54  lep = fPartList[i];
55  break;
56  }
57  }//done looping over particles
58 
59  fMCNeutrino = simb::MCNeutrino(nu, lep,
60  CCNC, mode, interactionType,
61  target, nucleon, quark,
62  w, x, y, qsqr);
63  } // end if MCNeutrino is not already set
64  else
65  throw art::Exception(art::errors::LogicError) << "MCTruth - attempt to set neutrino when already set";
66  return;
67  }
std::vector< simb::MCParticle > fPartList
list of particles in this event
Definition: MCTruth.h:38
const XML_Char * target
Definition: expat.h:268
void abs(TH1 *hist)
simb::MCNeutrino fMCNeutrino
reference to neutrino info - null if not a neutrino
Definition: MCTruth.h:39
int PdgCode(Resonance_t res, int Q)
(resonance id, charge) -> PDG code
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
bool fNeutrinoSet
flag for whether the neutrino information has been set
Definition: MCTruth.h:42
Event generator information.
Definition: MCNeutrino.h:18
Float_t w
Definition: plot.C:20
void simb::MCTruth::SetOrigin ( simb::Origin_t  origin)
inline

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  o,
simb::MCTruth const &  a 
)
friend

Definition at line 70 of file MCTruth.cxx.

71  {
72  if(a.Origin() == kCosmicRay)
73  o << "This is a cosmic ray event" << std::endl;
74  else if(a.Origin() == kBeamNeutrino){
75  o << "This is a beam neutrino event" << std::endl;
76  o << a.GetNeutrino();
77  }
78  else if(a.Origin() == kSuperNovaNeutrino){
79  o << "This is a supernova neutrino event" << std::endl;
80  o << a.GetNeutrino();
81  }
82 
83  for (int i = 0; i < a.NParticles(); ++i)
84  o << i << " " << a.GetParticle(i) << std::endl;
85 
86  return o;
87  }
const double a
Supernova neutrinos.
Definition: MCTruth.h:25
Cosmic rays.
Definition: MCTruth.h:24
Beam neutrinos.
Definition: MCTruth.h:23

Member Data Documentation

simb::MCGeneratorInfo simb::MCTruth::fGenInfo
private

information about the generator that produced this event

Definition at line 41 of file MCTruth.h.

Referenced by GeneratorInfo(), and SetGeneratorInfo().

simb::MCNeutrino simb::MCTruth::fMCNeutrino
private

reference to neutrino info - null if not a neutrino

Definition at line 39 of file MCTruth.h.

Referenced by GetNeutrino(), and SetNeutrino().

bool simb::MCTruth::fNeutrinoSet
private

flag for whether the neutrino information has been set

Definition at line 42 of file MCTruth.h.

Referenced by NeutrinoSet(), and SetNeutrino().

simb::Origin_t simb::MCTruth::fOrigin
private

origin for this event

Definition at line 40 of file MCTruth.h.

Referenced by Origin(), and SetOrigin().

std::vector<simb::MCParticle> simb::MCTruth::fPartList
private

list of particles in this event

Definition at line 38 of file MCTruth.h.

Referenced by Add(), GetParticle(), NParticles(), and SetNeutrino().


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