Public Member Functions | Private Member Functions | Private Attributes | List of all members
genie::HEDISPXSec Class Reference

Computes the double differential Cross Section for HEDIS.
Is a concrete implementation of the XSecAlgorithmI interface. More...

#include <HEDISPXSec.h>

Inheritance diagram for genie::HEDISPXSec:
genie::XSecAlgorithmI genie::Algorithm

Public Member Functions

 HEDISPXSec ()
 
 HEDISPXSec (string config)
 
virtual ~HEDISPXSec ()
 
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 config)
 
- 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...
 

Private Member Functions

void LoadConfig (void)
 
double ds_dxdy (SF_xQ2 sf, double x, double y) const
 
double ds_dxdy_mass (SF_xQ2 sf, double x, double y, double e, double mt, double ml2) const
 

Private Attributes

const XSecIntegratorIfXSecIntegrator
 diff. xsec integrator More...
 
double fWmin
 Minimum value of W. More...
 
bool fMassTerms
 Account for second order effects in DDxsec. More...
 
string fSFBaseDir
 Directory were SF will be stored. More...
 
SF_info fSFinfo
 Information used to computed SF. 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)
 
- 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 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...
 

Detailed Description

Computes the double differential Cross Section for HEDIS.
Is a concrete implementation of the XSecAlgorithmI interface.

Author
Alfonso Garcia <alfonsog nikhef.nl> NIKHEF

August 28, 2019

Copyright (c) 2003-2018, The GENIE Collaboration For the full text of the license visit http://copyright.genie-mc.org or see $GENIE/LICENSE

Definition at line 30 of file HEDISPXSec.h.

Constructor & Destructor Documentation

HEDISPXSec::HEDISPXSec ( )

Definition at line 29 of file HEDISPXSec.cxx.

29  :
30 XSecAlgorithmI("genie::HEDISPXSec")
31 {
32 
33 }
HEDISPXSec::HEDISPXSec ( string  config)

Definition at line 35 of file HEDISPXSec.cxx.

35  :
36 XSecAlgorithmI("genie::HEDISPXSec", config)
37 {
38 
39 }
static Config * config
Definition: config.cpp:1054
HEDISPXSec::~HEDISPXSec ( )
virtual

Definition at line 41 of file HEDISPXSec.cxx.

42 {
43 
44 }

Member Function Documentation

void HEDISPXSec::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 178 of file HEDISPXSec.cxx.

179 {
180  Algorithm::Configure(config);
181  this->LoadConfig();
182 }
void LoadConfig(void)
Definition: HEDISPXSec.cxx:190
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void HEDISPXSec::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 184 of file HEDISPXSec.cxx.

185 {
187  this->LoadConfig();
188 }
void LoadConfig(void)
Definition: HEDISPXSec.cxx:190
static Config * config
Definition: config.cpp:1054
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
double HEDISPXSec::ds_dxdy ( SF_xQ2  sf,
double  x,
double  y 
) const
private

Definition at line 117 of file HEDISPXSec.cxx.

118 {
119 
120  // We neglect F4 and F5 and higher order terms.
121  double term1 = y * ( x*y );
122  double term2 = ( 1 - y );
123  double term3 = ( x*y*(1-y/2) );
124 
125  LOG("HEDISPXSec", pDEBUG) << sf.F1 << " " << sf.F2 << " " << sf.F3;
126  LOG("HEDISPXSec", pDEBUG) << term1*sf.F1 + term2*sf.F2 + term3*sf.F3;
127 
128  return fmax( term1*sf.F1 + term2*sf.F2 + term3*sf.F3 , 0.);
129 
130 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
list x
Definition: train.py:276
#define pDEBUG
Definition: Messenger.h:63
double HEDISPXSec::ds_dxdy_mass ( SF_xQ2  sf,
double  x,
double  y,
double  e,
double  mt,
double  ml2 
) const
private

Definition at line 132 of file HEDISPXSec.cxx.

133 {
134 
135  double term1 = y * ( x*y + ml2/2/e/mt );
136  double term2 = ( 1 - y - mt*x*y/2/e - ml2/4/e/e );
137  double term3 = (x*y*(1-y/2) - y*ml2/4/mt/e);
138  double term4 = x*y*ml2/2/mt/e + ml2*ml2/4/mt/mt/e/e;
139  double term5 = -1.*ml2/2/mt/e;
140 
141  double F4 = 0.;
142  double F5 = sf.F2/x;
143 
144  LOG("HEDISPXSec", pDEBUG) << sf.F1 << " " << sf.F2 << " " << sf.F3;
145  LOG("HEDISPXSec", pDEBUG) << term1*sf.F1 + term2*sf.F2 + term3*sf.F3;
146 
147  return fmax( term1*sf.F1 + term2*sf.F2 + term3*sf.F3 + term4*F4 + term5*F5 , 0.);
148 
149 }
const double e
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
constexpr MyFlag_t F5
Definition: FlagSet_test.cc:44
#define F4(x, y, z)
Definition: md5.c:178
list x
Definition: train.py:276
#define pDEBUG
Definition: Messenger.h:63
double HEDISPXSec::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 151 of file HEDISPXSec.cxx.

152 {
153 
154  return fXSecIntegrator->Integrate(this,interaction);
155 
156 }
const XSecIntegratorI * fXSecIntegrator
diff. xsec integrator
Definition: HEDISPXSec.h:52
virtual double Integrate(const XSecAlgorithmI *model, const Interaction *interaction) const =0
void HEDISPXSec::LoadConfig ( void  )
private

Definition at line 190 of file HEDISPXSec.cxx.

191 {
192 
193  //-- load the differential cross section integrator
194  fXSecIntegrator = dynamic_cast<const XSecIntegratorI *> (this->SubAlg("XSec-Integrator"));
195  assert(fXSecIntegrator);
196 
197  // Minimum value of W (typically driven by hadronization limitation)
198  GetParam("Xsec-Wmin", fWmin);
199  GetParam("Mass-Terms", fMassTerms);
200 
201  GetParamDef("BaseDir", fSFBaseDir, string(""));
202 
203  // Information about Structure Functions
204  GetParam("LHAPDF-set", fSFinfo.LHAPDFset );
205  GetParam("LHAPDF-member", fSFinfo.LHAPDFmember );
206  GetParam("Is-NLO", fSFinfo.IsNLO );
207  GetParam("Scheme", fSFinfo.Scheme );
208  GetParam("Quark-Threshold", fSFinfo.QrkThrs );
209  GetParam("NGridX", fSFinfo.NGridX );
210  GetParam("NGridQ2", fSFinfo.NGridQ2 );
211  GetParam("XGrid-Min", fSFinfo.XGridMin );
212  GetParam("Q2Grid-Min", fSFinfo.Q2GridMin );
213  GetParam("Q2Grid-Max", fSFinfo.Q2GridMax );
214  GetParam("MassW", fSFinfo.MassW );
215  GetParam("MassZ", fSFinfo.MassZ );
216  GetParam("Rho", fSFinfo.Rho );
217  GetParam("Sin2ThW", fSFinfo.Sin2ThW );
218  GetParam("Mud", fSFinfo.Vud );
219  GetParam("Mus", fSFinfo.Vus );
220  GetParam("Mub", fSFinfo.Vub );
221  GetParam("Mcd", fSFinfo.Vcd );
222  GetParam("Mcs", fSFinfo.Vcs );
223  GetParam("Mcb", fSFinfo.Vcb );
224  GetParam("Mtd", fSFinfo.Vtd );
225  GetParam("Mts", fSFinfo.Vts );
226  GetParam("Mtb", fSFinfo.Vtb );
227 
228 }
Cross Section Integrator Interface.
const XSecIntegratorI * fXSecIntegrator
diff. xsec integrator
Definition: HEDISPXSec.h:52
double fWmin
Minimum value of W.
Definition: HEDISPXSec.h:54
string fSFBaseDir
Directory were SF will be stored.
Definition: HEDISPXSec.h:56
SF_info fSFinfo
Information used to computed SF.
Definition: HEDISPXSec.h:57
bool GetParamDef(const RgKey &name, T &p, const T &def) const
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
bool fMassTerms
Account for second order effects in DDxsec.
Definition: HEDISPXSec.h:55
const Algorithm * SubAlg(const RgKey &registry_key) const
Definition: Algorithm.cxx:345
bool HEDISPXSec::ValidProcess ( const Interaction i) const
virtual

Can this cross section algorithm handle the input process?

Implements genie::XSecAlgorithmI.

Definition at line 158 of file HEDISPXSec.cxx.

159 {
160 
161  if(interaction->TestBit(kISkipProcessChk)) return true;
162 
163  const ProcessInfo & proc_info = interaction->ProcInfo();
164  if(!proc_info.IsDeepInelastic()) return false;
165 
166  const InitialState & init_state = interaction -> InitState();
167  int probe_pdg = init_state.ProbePdg();
168  if(!pdg::IsLepton(probe_pdg)) return false;
169 
170  if(! init_state.Tgt().HitNucIsSet()) return false;
171 
172  int hitnuc_pdg = init_state.Tgt().HitNucPdg();
173  if(!pdg::IsNeutronOrProton(hitnuc_pdg)) return false;
174 
175  return true;
176 }
int HitNucPdg(void) const
Definition: Target.cxx:304
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
int ProbePdg(void) const
Definition: InitialState.h:64
bool IsDeepInelastic(void) const
Definition: ProcessInfo.cxx:84
bool HitNucIsSet(void) const
Definition: Target.cxx:283
bool IsNeutronOrProton(int pdgc)
Definition: PDGUtils.cxx:348
const Target & Tgt(void) const
Definition: InitialState.h:66
bool IsLepton(int pdgc)
Definition: PDGUtils.cxx:83
const UInt_t kISkipProcessChk
if set, skip process validity checks
Definition: Interaction.h:47
Initial State information.
Definition: InitialState.h:48
double HEDISPXSec::XSec ( const Interaction i,
KinePhaseSpace_t  k 
) const
virtual

Compute the cross section for the input interaction.

Implements genie::XSecAlgorithmI.

Definition at line 46 of file HEDISPXSec.cxx.

48 {
49 
50  if(! this -> ValidKinematics (interaction) ) return 0.;
51 
52  // Load SF tables
54 
55  // W limits are computed using kinematics assumption.
56  // The lower limit is tuneable because hadronization might have issues at low W (as in PYTHIA6).
57  // To be consistent the cross section must be computed in the W range where the events are generated.
58  const Kinematics & kinematics = interaction -> Kine();
59  const KPhaseSpace & ps = interaction -> PhaseSpace();
60  double W = kinematics.W();
61  Range1D_t Wl = ps.WLim();
62  Wl.min = TMath::Max(Wl.min,fWmin);
63  if ( W<Wl.min ) return 0.;
64  else if ( W>Wl.max ) return 0.;
65 
66  const InitialState & init_state = interaction -> InitState();
67 
68  double y = kinematics.y();
69  double Q2 = kinematics.Q2();
70  double x = kinematics.x();
71  double E = init_state.ProbeE(kRfLab);
72  double Mnuc = init_state.Tgt().HitNucMass();
73  double Mlep2 = TMath::Power(interaction->FSPrimLepton()->Mass(),2);
74 
75  // Get F1,F2,F3 for particular quark channel and compute differential xsec
76  SF_xQ2 sf = sf_tbl->EvalQrkSFLO( interaction, x, Q2 );
77  double xsec = (fMassTerms) ? ds_dxdy_mass( sf, x, y, E, Mnuc, Mlep2 ) : ds_dxdy( sf, x, y );
78 
79  // If NLO is enable we compute sigma_NLO/sigma_LO. Then the quark xsec
80  // is multiplied by this ratio.
81  // This is done because at NLO we can only compute the nucleon xsec. But
82  // for the hadronization we need the different quark contributions.
83  // This could be avoid if a NLO parton showering is introduced.
84  if (fSFinfo.IsNLO && xsec>0.) {
85  SF_xQ2 sflo = sf_tbl->EvalNucSFLO(interaction,x,Q2);
86  SF_xQ2 sfnlo = sf_tbl->EvalNucSFNLO(interaction,x,Q2);
87  double lo = (fMassTerms) ? ds_dxdy_mass( sflo, x, y, E, Mnuc, Mlep2 ) : ds_dxdy( sflo, x, y );
88  if (lo>0.) {
89  double nlo = (fMassTerms) ? ds_dxdy_mass( sfnlo, x, y, E, Mnuc, Mlep2 ) : ds_dxdy( sfnlo, x, y );
90  xsec *= nlo / lo;
91  }
92  }
93 
94  // Compute the front factor
95  double propagator = 0;
96  if (interaction -> ProcInfo().IsWeakCC()) propagator = TMath::Power( fSFinfo.MassW*fSFinfo.MassW/(Q2+fSFinfo.MassW*fSFinfo.MassW), 2);
97  else propagator = TMath::Power( fSFinfo.MassZ*fSFinfo.MassZ/(Q2+fSFinfo.MassZ*fSFinfo.MassZ)/(1.-fSFinfo.Rho), 2);
98 
99  xsec *= kGF2/(2*kPi*x) * propagator;
100 
101  LOG("HEDISPXSec", pINFO) << "d2xsec/dxdy[FreeN] (x= " << x << ", y= " << y << ", Q2= " << Q2 << ") = " << xsec;
102 
103  // The algorithm computes d^2xsec/dxdy. Check whether variable tranformation is needed
105 
106  // If requested return the free nucleon xsec even for input nuclear tgt
107  if( interaction->TestBit(kIAssumeFreeNucleon) ) return xsec;
108 
109  // Compute nuclear cross section (simple scaling here, corrections must have been included in the structure functions)
110  int NNucl = (pdg::IsProton(init_state.Tgt().HitNucPdg())) ? init_state.Tgt().Z() : init_state.Tgt().N();
111  xsec *= NNucl;
112 
113  return xsec;
114 
115 }
double W(bool selected=false) const
Definition: Kinematics.cxx:157
double ds_dxdy(SF_xQ2 sf, double x, double y) const
Definition: HEDISPXSec.cxx:117
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1064
int HitNucPdg(void) const
Definition: Target.cxx:304
A simple [min,max] interval for doubles.
Definition: Range1.h:42
double HitNucMass(void) const
Definition: Target.cxx:233
SF_xQ2 EvalNucSFLO(const Interaction *in, double x, double Q2)
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
double x(bool selected=false) const
Definition: Kinematics.cxx:99
double y(bool selected=false) const
Definition: Kinematics.cxx:112
virtual bool ValidKinematics(const Interaction *i) const
Is the input kinematical point a physically allowed one?
bool IsProton(int pdgc)
Definition: PDGUtils.cxx:333
double fWmin
Minimum value of W.
Definition: HEDISPXSec.h:54
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
double ds_dxdy_mass(SF_xQ2 sf, double x, double y, double e, double mt, double ml2) const
Definition: HEDISPXSec.cxx:132
Kinematical phase space.
Definition: KPhaseSpace.h:33
int Z(void) const
Definition: Target.h:68
#define pINFO
Definition: Messenger.h:62
static constexpr double ps
Definition: Units.h:99
double max
Definition: Range1.h:53
int N(void) const
Definition: Target.h:69
const UInt_t kIAssumeFreeNucleon
Definition: Interaction.h:49
string fSFBaseDir
Directory were SF will be stored.
Definition: HEDISPXSec.h:56
double Jacobian(const Interaction *const i, KinePhaseSpace_t f, KinePhaseSpace_t t)
Definition: KineUtils.cxx:130
SF_xQ2 EvalQrkSFLO(const Interaction *in, double x, double Q2)
list x
Definition: train.py:276
double min
Definition: Range1.h:52
SF_info fSFinfo
Information used to computed SF.
Definition: HEDISPXSec.h:57
double Q2(bool selected=false) const
Definition: Kinematics.cxx:125
const Target & Tgt(void) const
Definition: InitialState.h:66
static const double kGF2
Definition: Constants.h:59
bool fMassTerms
Account for second order effects in DDxsec.
Definition: HEDISPXSec.h:55
double ProbeE(RefFrame_t rf) const
static const double kPi
Definition: Constants.h:37
SF_xQ2 EvalNucSFNLO(const Interaction *in, double x, double Q2)
Range1D_t WLim(void) const
W limits.
static HEDISStrucFunc * Instance(string basedir, SF_info sfinfo)
Initial State information.
Definition: InitialState.h:48

Member Data Documentation

bool genie::HEDISPXSec::fMassTerms
private

Account for second order effects in DDxsec.

Definition at line 55 of file HEDISPXSec.h.

string genie::HEDISPXSec::fSFBaseDir
private

Directory were SF will be stored.

Definition at line 56 of file HEDISPXSec.h.

SF_info genie::HEDISPXSec::fSFinfo
private

Information used to computed SF.

Definition at line 57 of file HEDISPXSec.h.

double genie::HEDISPXSec::fWmin
private

Minimum value of W.

Definition at line 54 of file HEDISPXSec.h.

const XSecIntegratorI* genie::HEDISPXSec::fXSecIntegrator
private

diff. xsec integrator

Definition at line 52 of file HEDISPXSec.h.


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