Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
genie::evtlib::EvtLibPXSec Class Reference

#include <EvtLibPXSec.h>

Inheritance diagram for genie::evtlib::EvtLibPXSec:
genie::XSecAlgorithmI genie::Algorithm

Public Member Functions

 EvtLibPXSec ()
 
 EvtLibPXSec (string config)
 
virtual ~EvtLibPXSec ()
 
double XSec (const Interaction *i, KinePhaseSpace_t k) const
 Compute the cross section for the input interaction. More...
 
double Integral (const Interaction *i) const
 
bool ValidProcess (const Interaction *i) const
 Can this cross section algorithm handle the input process? More...
 
void Configure (const Registry &config)
 
void Configure (string param_set)
 
- Public Member Functions inherited from genie::XSecAlgorithmI
virtual ~XSecAlgorithmI ()
 
virtual bool ValidKinematics (const Interaction *i) const
 Is the input kinematical point a physically allowed one? More...
 
- Public Member Functions inherited from genie::Algorithm
virtual ~Algorithm ()
 
virtual void FindConfig (void)
 
virtual const RegistryGetConfig (void) const
 
RegistryGetOwnedConfig (void)
 
virtual const AlgIdId (void) const
 Get algorithm ID. More...
 
virtual AlgStatus_t GetStatus (void) const
 Get algorithm status. More...
 
virtual bool AllowReconfig (void) const
 
virtual AlgCmp_t Compare (const Algorithm *alg) const
 Compare with input algorithm. More...
 
virtual void SetId (const AlgId &id)
 Set algorithm ID. More...
 
virtual void SetId (string name, string config)
 
const AlgorithmSubAlg (const RgKey &registry_key) const
 
void AdoptConfig (void)
 
void AdoptSubstructure (void)
 
virtual void Print (ostream &stream) const
 Print algorithm info. More...
 

Protected Member Functions

TGraph * GetXSec (const Interaction *in) const
 
void LoadXSecs ()
 
void ClearXSecs ()
 
- Protected Member Functions inherited from genie::XSecAlgorithmI
 XSecAlgorithmI ()
 
 XSecAlgorithmI (string name)
 
 XSecAlgorithmI (string name, string config)
 
- Protected Member Functions inherited from genie::Algorithm
 Algorithm ()
 
 Algorithm (string name)
 
 Algorithm (string name, string config)
 
void Initialize (void)
 
void DeleteConfig (void)
 
void DeleteSubstructure (void)
 
RegistryExtractLocalConfig (const Registry &in) const
 
RegistryExtractLowerConfig (const Registry &in, const string &alg_key) const
 Split an incoming configuration Registry into a block valid for the sub-algo identified by alg_key. More...
 
template<class T >
bool GetParam (const RgKey &name, T &p, bool is_top_call=true) const
 
template<class T >
bool GetParamDef (const RgKey &name, T &p, const T &def) const
 
template<class T >
int GetParamVect (const std::string &comm_name, std::vector< T > &v, bool is_top_call=true) const
 Handle to load vectors of parameters. More...
 
int GetParamVectKeys (const std::string &comm_name, std::vector< RgKey > &k, bool is_top_call=true) const
 
int AddTopRegistry (Registry *rp, bool owns=true)
 add registry with top priority, also update ownership More...
 
int AddLowRegistry (Registry *rp, bool owns=true)
 add registry with lowest priority, also update ownership More...
 
int MergeTopRegistry (const Registry &r)
 
int AddTopRegisties (const vector< Registry * > &rs, bool owns=false)
 Add registries with top priority, also udated Ownerships. More...
 

Protected Attributes

std::map< Key, TGraph * > fXSecs
 
- Protected Attributes inherited from genie::Algorithm
bool fAllowReconfig
 
bool fOwnsSubstruc
 true if it owns its substructure (sub-algs,...) More...
 
AlgId fID
 algorithm name and configuration set More...
 
vector< Registry * > fConfVect
 
vector< boolfOwnerships
 ownership for every registry in fConfVect More...
 
AlgStatus_t fStatus
 algorithm execution status More...
 
AlgMapfOwnedSubAlgMp
 local pool for owned sub-algs (taken out of the factory pool) More...
 

Additional Inherited Members

- Static Public Member Functions inherited from genie::Algorithm
static string BuildParamVectKey (const std::string &comm_name, unsigned int i)
 
static string BuildParamVectSizeKey (const std::string &comm_name)
 

Detailed Description

Definition at line 13 of file EvtLibPXSec.h.

Constructor & Destructor Documentation

EvtLibPXSec::EvtLibPXSec ( )

Definition at line 17 of file EvtLibPXSec.cxx.

17  :
18 XSecAlgorithmI("genie::evtlib::EvtLibPXSec")
19 {
20 
21 }
EvtLibPXSec::EvtLibPXSec ( string  config)

Definition at line 24 of file EvtLibPXSec.cxx.

24  :
25 XSecAlgorithmI("genie::evtlib::EvtLibPXSec", config)
26 {
27 
28 }
static Config * config
Definition: config.cpp:1054
EvtLibPXSec::~EvtLibPXSec ( )
virtual

Definition at line 31 of file EvtLibPXSec.cxx.

32 {
33  ClearXSecs();
34 }

Member Function Documentation

void EvtLibPXSec::ClearXSecs ( )
protected

Definition at line 78 of file EvtLibPXSec.cxx.

79 {
80  for(auto it: fXSecs) delete it.second;
81  fXSecs.clear();
82 }
std::map< Key, TGraph * > fXSecs
Definition: EvtLibPXSec.h:35
void EvtLibPXSec::Configure ( const Registry config)
virtual

Configure the algorithm with an external registry The registry is merged with the top level registry if it is owned, Otherwise a copy of it is added with the highest priority

Reimplemented from genie::Algorithm.

Definition at line 64 of file EvtLibPXSec.cxx.

65 {
66  Algorithm::Configure(config);
67  LoadXSecs();
68 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void EvtLibPXSec::Configure ( string  config)
virtual

Configure the algorithm from the AlgoConfigPool based on param_set string given in input An algorithm contains a vector of registries coming from different xml configuration files, which are loaded according a very precise prioriy This methods will load a number registries in order of priority: 1) "Tunable" parameter set from CommonParametes. This is loaded with the highest prioriry and it is designed to be used for tuning procedure Usage not expected from the user. 2) For every string defined in "CommonParame" the corresponding parameter set will be loaded from CommonParameter.xml 3) parameter set specified by the config string and defined in the xml file of the algorithm 4) if config is not "Default" also the Default parameter set from the same xml file will be loaded Effectively this avoids the repetion of a parameter when it is not changed in the requested configuration

Reimplemented from genie::Algorithm.

Definition at line 71 of file EvtLibPXSec.cxx.

72 {
74  LoadXSecs();
75 }
static Config * config
Definition: config.cpp:1054
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
TGraph * EvtLibPXSec::GetXSec ( const Interaction in) const
protected

Definition at line 143 of file EvtLibPXSec.cxx.

144 {
145  const InitialState& init_state = in->InitState();
146 
147  const ProcessInfo& proc = in->ProcInfo();
148  if(!proc.IsWeakCC() && !proc.IsWeakNC()){
149  LOG("ELI", pINFO) << "Skipping unknown process " << proc;
150  return 0;
151  }
152 
153  int pdg = init_state.ProbePdg();
154  // Use nu_mu for NC as a convention internal to this code to index into the
155  // xsec graphs map.
156  if(proc.IsWeakNC()){
157  if(pdg > 0) pdg = kPdgNuMu; else pdg = kPdgAntiNuMu;
158  }
159 
160  const Key key(init_state.TgtPdg(), pdg, proc.IsWeakCC());
161 
162  auto it = fXSecs.find(key);
163  if(it == fXSecs.end()){
164  LOG("ELI", pINFO) << "Skipping process without xsec " << key;
165  return 0;
166  }
167 
168  return it->second;
169 }
bool IsWeakCC(void) const
const int kPdgNuMu
Definition: PDGCodes.h:30
bool IsWeakNC(void) const
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
def key(type, name=None)
Definition: graph.py:13
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
std::map< Key, TGraph * > fXSecs
Definition: EvtLibPXSec.h:35
int ProbePdg(void) const
Definition: InitialState.h:64
#define pINFO
Definition: Messenger.h:62
const int kPdgAntiNuMu
Definition: PDGCodes.h:31
int TgtPdg(void) const
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
Initial State information.
Definition: InitialState.h:48
double EvtLibPXSec::Integral ( const Interaction i) const
virtual

Integrate the model over the kinematic phase space available to the input interaction (kinematical cuts can be included)

Implements genie::XSecAlgorithmI.

Definition at line 45 of file EvtLibPXSec.cxx.

46 {
47  TGraph* g = GetXSec(in);
48  if(!g) return 0; // Reason already printed
49 
50  const InitialState& init_state = in->InitState();
51  const double E = init_state.ProbeE(kRfLab);
52 
53  // Units of the cross-section graph are expected to be 10^-38 cm^2 / nucleus
54  return g->Eval(E) * 1e-38 * genie::units::cm2;
55 }
static constexpr double g
Definition: Units.h:144
const double e
static constexpr double cm2
Definition: Units.h:69
TGraph * GetXSec(const Interaction *in) const
double ProbeE(RefFrame_t rf) const
Initial State information.
Definition: InitialState.h:48
void EvtLibPXSec::LoadXSecs ( )
protected

Definition at line 85 of file EvtLibPXSec.cxx.

86 {
87  ClearXSecs();
88 
89  std::string libPath;
90  GetParam("EventLibraryPath", libPath);
91  Expand(libPath);
92 
93  PDGLibrary* pdglib = PDGLibrary::Instance();
94 
95  TFile fin(libPath.c_str());
96  if(fin.IsZombie()) exit(1);
97 
98  TIter next(fin.GetListOfKeys());
99  while(TObject* dir = next()){
100  const std::string& tgtName = dir->GetName();
101  const TParticlePDG* tgtPart = pdglib->DBase()->GetParticle(tgtName.c_str());
102  if(!tgtPart){
103  LOG("ELI", pWARN) << "Unknown nucleus " << tgtName
104  << " found in " << libPath
105  << " -- skipping";
106  continue;
107  }
108 
109  for(int pdg: {kPdgNuE, kPdgAntiNuE,
112 
113  for(bool iscc: {true, false}){
114  // NCs should be the same for all flavours. Use nu_mu as a convention
115  // internal to this code to index into the xsecs map.
116  if(!iscc && abs(pdg) != kPdgNuMu) continue;
117 
118  std::string nuName = pdglib->Find(pdg)->GetName();
119  if(!iscc) nuName = pdg::IsAntiNeutrino(pdg) ? "nu_bar" : "nu";
120 
121  const std::string graphName =
122  TString::Format("%s/%s/%s/xsec",
123  tgtName.c_str(),
124  iscc ? "cc" : "nc",
125  nuName.c_str()).Data();
126 
127  const Key key(tgtPart->PdgCode(), pdg, iscc);
128 
129  TGraph* g = (TGraph*)fin.Get(graphName.c_str());
130  if(!g){
131  LOG("ELI", pINFO) << graphName << " not found in "
132  << libPath << " -- skipping";
133  continue;
134  }
135 
136  fXSecs[key] = new TGraph(*g);
137  } // end for iscc
138  } // end for pdg
139  } // end for dir
140 }
const int kPdgNuE
Definition: PDGCodes.h:28
Format
Definition: utils.h:7
static constexpr double g
Definition: Units.h:144
std::string string
Definition: nybbler.cc:12
const int kPdgAntiNuE
Definition: PDGCodes.h:29
void Expand(std::string &s)
Expand env vars/homedirs/wildcards in s.
Definition: Utils.cxx:8
TDatabasePDG * DBase(void)
Definition: PDGLibrary.cxx:70
const int kPdgNuMu
Definition: PDGCodes.h:30
string dir
T abs(T value)
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
def key(type, name=None)
Definition: graph.py:13
std::map< Key, TGraph * > fXSecs
Definition: EvtLibPXSec.h:35
bool IsAntiNeutrino(int pdgc)
Definition: PDGUtils.cxx:115
#define pINFO
Definition: Messenger.h:62
const int kPdgAntiNuTau
Definition: PDGCodes.h:33
const int kPdgAntiNuMu
Definition: PDGCodes.h:31
#define pWARN
Definition: Messenger.h:60
Fw2dFFT::Data Data
const int kPdgNuTau
Definition: PDGCodes.h:32
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:57
Singleton class to load & serve a TDatabasePDG.
Definition: PDGLibrary.h:32
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:75
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
bool EvtLibPXSec::ValidProcess ( const Interaction i) const
virtual

Can this cross section algorithm handle the input process?

Implements genie::XSecAlgorithmI.

Definition at line 58 of file EvtLibPXSec.cxx.

59 {
60  return true;
61 }
double EvtLibPXSec::XSec ( const Interaction i,
KinePhaseSpace_t  k 
) const
virtual

Compute the cross section for the input interaction.

Implements genie::XSecAlgorithmI.

Definition at line 37 of file EvtLibPXSec.cxx.

39 {
40  LOG("ELI", pFATAL) << "This point should not be reached" ;
41  abort();
42 }
#define pFATAL
Definition: Messenger.h:56
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96

Member Data Documentation

std::map<Key, TGraph*> genie::evtlib::EvtLibPXSec::fXSecs
protected

Definition at line 35 of file EvtLibPXSec.h.


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