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

Computes, at Leading Order (LO), the differential cross section for neutrino charm production using a slow rescaling model. Is a concrete implementation of the XSecAlgorithmI interface. More...

#include <SlowRsclCharmDISPXSecLO.h>

Inheritance diagram for genie::SlowRsclCharmDISPXSecLO:
genie::XSecAlgorithmI genie::Algorithm

Public Member Functions

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

Private Member Functions

void LoadConfig (void)
 

Private Attributes

const PDFModelIfPDFModel
 
const XSecIntegratorIfXSecIntegrator
 
double fMc
 
double fVcd
 
double fVcs
 
double fMc2
 
double fVcd2
 
double fVcs2
 

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, at Leading Order (LO), the differential cross section for neutrino charm production using a slow rescaling model. Is a concrete implementation of the XSecAlgorithmI interface.

Author
Costas Andreopoulos <constantinos.andreopoulos cern.ch> University of Liverpool & STFC Rutherford Appleton Laboratory

March 17, 2005

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

Definition at line 30 of file SlowRsclCharmDISPXSecLO.h.

Constructor & Destructor Documentation

SlowRsclCharmDISPXSecLO::SlowRsclCharmDISPXSecLO ( )

Definition at line 35 of file SlowRsclCharmDISPXSecLO.cxx.

35  :
36 XSecAlgorithmI("genie::SlowRsclCharmDISPXSecLO")
37 {
38 
39 }
SlowRsclCharmDISPXSecLO::SlowRsclCharmDISPXSecLO ( string  config)

Definition at line 41 of file SlowRsclCharmDISPXSecLO.cxx.

41  :
42 XSecAlgorithmI("genie::SlowRsclCharmDISPXSecLO", config)
43 {
44 
45 }
static Config * config
Definition: config.cpp:1054
SlowRsclCharmDISPXSecLO::~SlowRsclCharmDISPXSecLO ( )
virtual

Definition at line 47 of file SlowRsclCharmDISPXSecLO.cxx.

48 {
49 
50 }

Member Function Documentation

void SlowRsclCharmDISPXSecLO::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 184 of file SlowRsclCharmDISPXSecLO.cxx.

185 {
186  Algorithm::Configure(config);
187  this->LoadConfig();
188 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void SlowRsclCharmDISPXSecLO::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 190 of file SlowRsclCharmDISPXSecLO.cxx.

191 {
192  Algorithm::Configure(param_set);
193  this->LoadConfig();
194 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
double SlowRsclCharmDISPXSecLO::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 154 of file SlowRsclCharmDISPXSecLO.cxx.

155 {
156  double xsec = fXSecIntegrator->Integrate(this,interaction);
157  return xsec;
158 }
virtual double Integrate(const XSecAlgorithmI *model, const Interaction *interaction) const =0
const XSecIntegratorI * fXSecIntegrator
void SlowRsclCharmDISPXSecLO::LoadConfig ( void  )
private

Definition at line 196 of file SlowRsclCharmDISPXSecLO.cxx.

197 {
198  // read mc, Vcd, Vcs from config or set defaults
199  GetParam( "Charm-Mass", fMc ) ;
200  GetParam( "CKM-Vcd", fVcd ) ;
201  GetParam( "CKM-Vcs", fVcs ) ;
202 
203  fMc2 = TMath::Power(fMc, 2);
204  fVcd2 = TMath::Power(fVcd, 2);
205  fVcs2 = TMath::Power(fVcs, 2);
206 
207  // load PDF set
208  fPDFModel = dynamic_cast<const PDFModelI *> (this->SubAlg("PDF-Set"));
209  assert(fPDFModel);
210 
211  // load XSec Integrator
213  dynamic_cast<const XSecIntegratorI *> (this->SubAlg("XSec-Integrator"));
214  assert(fXSecIntegrator);
215 }
Cross Section Integrator Interface.
Pure abstract base class. Defines the PDFModelI interface to be implemented by wrapper classes to exi...
Definition: PDFModelI.h:28
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
const XSecIntegratorI * fXSecIntegrator
const Algorithm * SubAlg(const RgKey &registry_key) const
Definition: Algorithm.cxx:345
bool SlowRsclCharmDISPXSecLO::ValidProcess ( const Interaction i) const
virtual

Can this cross section algorithm handle the input process?

Implements genie::XSecAlgorithmI.

Definition at line 160 of file SlowRsclCharmDISPXSecLO.cxx.

162 {
163  if(interaction->TestBit(kISkipProcessChk)) return true;
164 
165  const XclsTag & xcls = interaction->ExclTag();
166  const InitialState & init_state = interaction->InitState();
167  const ProcessInfo & proc_info = interaction->ProcInfo();
168 
169  if(!proc_info.IsDeepInelastic()) return false;
170  if(!proc_info.IsWeak()) return false;
171 
172  bool is_inclusive_charm = (xcls.IsCharmEvent() && xcls.IsInclusiveCharm());
173  if(!is_inclusive_charm) return false;
174 
175  int nu = init_state.ProbePdg();
176  int nuc = init_state.Tgt().HitNucPdg();
177 
178  if (!pdg::IsProton(nuc) && !pdg::IsNeutron(nuc)) return false;
179  if (!pdg::IsNeutrino(nu) && !pdg::IsAntiNeutrino(nu)) return false;
180 
181  return true;
182 }
bool IsWeak(void) const
bool IsNeutrino(int pdgc)
Definition: PDGUtils.cxx:107
int HitNucPdg(void) const
Definition: Target.cxx:304
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
bool IsCharmEvent(void) const
Definition: XclsTag.h:50
bool IsNeutron(int pdgc)
Definition: PDGUtils.cxx:338
bool IsProton(int pdgc)
Definition: PDGUtils.cxx:333
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
bool IsAntiNeutrino(int pdgc)
Definition: PDGUtils.cxx:115
int ProbePdg(void) const
Definition: InitialState.h:64
bool IsDeepInelastic(void) const
Definition: ProcessInfo.cxx:84
bool IsInclusiveCharm(void) const
Definition: XclsTag.cxx:54
const Target & Tgt(void) const
Definition: InitialState.h:66
const UInt_t kISkipProcessChk
if set, skip process validity checks
Definition: Interaction.h:47
Initial State information.
Definition: InitialState.h:48
double SlowRsclCharmDISPXSecLO::XSec ( const Interaction i,
KinePhaseSpace_t  k 
) const
virtual

Compute the cross section for the input interaction.

Implements genie::XSecAlgorithmI.

Definition at line 52 of file SlowRsclCharmDISPXSecLO.cxx.

54 {
55  if(! this -> ValidProcess (interaction) ) return 0.;
56  if(! this -> ValidKinematics (interaction) ) return 0.;
57 
58  if(interaction->ProcInfo().IsWeakNC()) return 0;
59 
60  //----- get kinematics & init-state parameters
61  const Kinematics & kinematics = interaction->Kine();
62  const InitialState & init_state = interaction->InitState();
63  const Target & target = init_state.Tgt();
64 
65  double Mnuc = target.HitNucMass();
66  double E = init_state.ProbeE(kRfHitNucRest);
67  double x = kinematics.x();
68  double y = kinematics.y();
69  double Q2 = 2*Mnuc*E*x*y;
70 
71  //----- get target information (hit nucleon and quark)
72  int nu = init_state.ProbePdg();
73  int nuc = target.HitNucPdg();
74  bool isP = pdg::IsProton (nuc);
75  bool isN = pdg::IsNeutron(nuc);
76  bool qset = target.HitQrkIsSet();
77  int qpdg = (qset) ? target.HitQrkPdg() : 0;
78  bool sea = (qset) ? target.HitSeaQrk() : false;
79  bool isd = (qset) ? pdg::IsDQuark (qpdg) : false;
80  bool iss = (qset) ? pdg::IsSQuark (qpdg) : false;
81  bool isdb = (qset) ? pdg::IsAntiDQuark (qpdg) : false;
82  bool issb = (qset) ? pdg::IsAntiSQuark (qpdg) : false;
83  bool isnu = pdg::IsNeutrino(nu);
84  bool isnub = pdg::IsAntiNeutrino(nu);
85 
86  //----- compute slow rescaling variable & check its value
87  double xi = x * (1 + fMc2/Q2);
88  if(xi<=0 || xi>1) return 0.;
89 
90  //----- calculate the PDFs
91  PDF pdfs;
92  pdfs.SetModel(fPDFModel); // <-- attach algorithm
93  pdfs.Calculate(xi, Q2); // <-- calculate
94 
95  //----- proton pdfs
96  double us = pdfs.UpSea()/xi;
97  double uv = pdfs.UpValence()/xi;
98  double ds = pdfs.DownSea()/xi;
99  double dv = pdfs.DownValence()/xi;
100  double s = pdfs.Strange()/xi;
101 
102  //----- if the hit nucleon is a neutron, swap u<->d
103  double tmp;
104  if (isN) {
105  tmp = uv; uv = dv; dv = tmp;
106  tmp = us; us = ds; ds = tmp;
107  }
108 
109  //----- if a hit quark has been set then switch off other contributions
110  if(qset) {
111  if(isnub) { bool pass = (isdb||issb)&&sea; if(!pass) return 0; }
112  if(isnu) { bool pass = isd||(iss&&sea); if(!pass) return 0; }
113  dv = ( isd && !sea) ? dv : 0.;
114  ds = ( (isd||isdb) && sea) ? ds : 0.;
115  s = ( (iss||issb) && sea) ? s : 0.;
116  }
117  //----- in case of a \bar{v}+N calculation (quark not set) zero the d/val contribution
118  if(isnub) dv=0;
119 
120  //----- calculate cross section
121  double Gw2 = TMath::Power( (kGF/kSqrt2)*(1+Q2/kMw2), 2);
122  double xsec_0 = Gw2 * 2*Q2/(y*kPi) * (y + xi*(1-y)/x);
123  double xsec_d = xsec_0 * fVcd2 * (dv+ds);
124  double xsec_s = xsec_0 * fVcs2 * s;
125  double xsec = xsec_d + xsec_s;
126 
127 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
128  double Mnuc2 = TMath::Power(Mnuc, 2);
129  double W2 = Mnuc2 + 2*Mnuc*E*y*(1-x);
130  double W = TMath::Max(0., TMath::Sqrt(W2));
131  LOG("DISCharmXSec", pDEBUG)
132  << "\n dxsec[DISCharm,FreeN]/dxdy (E= " << E
133  << ", x= " << x << ", y= " << y
134  << ", W= " << W << ", Q2 = " << Q2 << ") = " << xsec;
135 #endif
136 
137  //----- The algorithm computes d^2xsec/dxdy
138  // Check whether variable tranformation is needed
139  if(kps!=kPSxyfE) {
141  xsec *= J;
142  }
143 
144  //----- If requested return the free nucleon xsec even for input nuclear tgt
145  if( interaction->TestBit(kIAssumeFreeNucleon) ) return xsec;
146 
147  //----- Nuclear cross section (simple scaling here)
148  int NNucl = (isP) ? target.Z() : target.N();
149  xsec *= NNucl;
150 
151  return xsec;
152 }
bool HitSeaQrk(void) const
Definition: Target.cxx:299
static const double kSqrt2
Definition: Constants.h:115
static constexpr double us
Definition: Units.h:97
bool IsNeutrino(int pdgc)
Definition: PDGUtils.cxx:107
static const double kMw2
Definition: Constants.h:93
double J(double q0, double q3, double Enu, double ml)
Definition: MECUtils.cxx:147
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1064
int HitNucPdg(void) const
Definition: Target.cxx:304
double DownValence(void) const
Definition: PDF.h:51
int HitQrkPdg(void) const
Definition: Target.cxx:242
double HitNucMass(void) const
Definition: Target.cxx:233
A class to store PDFs.
Definition: PDF.h:37
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
double x(bool selected=false) const
Definition: Kinematics.cxx:99
bool IsSQuark(int pdgc)
Definition: PDGUtils.cxx:273
double UpSea(void) const
Definition: PDF.h:52
bool IsAntiSQuark(int pdgc)
Definition: PDGUtils.cxx:303
bool IsAntiDQuark(int pdgc)
Definition: PDGUtils.cxx:298
double y(bool selected=false) const
Definition: Kinematics.cxx:112
double DownSea(void) const
Definition: PDF.h:53
bool IsNeutron(int pdgc)
Definition: PDGUtils.cxx:338
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 Strange(void) const
Definition: PDF.h:54
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool ValidProcess(const Interaction *i) const
Can this cross section algorithm handle the input process?
HLTPathStatus const pass
bool IsAntiNeutrino(int pdgc)
Definition: PDGUtils.cxx:115
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
int ProbePdg(void) const
Definition: InitialState.h:64
void SetModel(const PDFModelI *model)
Definition: PDF.cxx:42
int Z(void) const
Definition: Target.h:68
static const double kGF
Definition: Constants.h:58
string tmp
Definition: languages.py:63
void Calculate(double x, double q2)
Definition: PDF.cxx:49
int N(void) const
Definition: Target.h:69
double UpValence(void) const
Definition: PDF.h:50
bool HitQrkIsSet(void) const
Definition: Target.cxx:292
const UInt_t kIAssumeFreeNucleon
Definition: Interaction.h:49
bool IsDQuark(int pdgc)
Definition: PDGUtils.cxx:268
double Jacobian(const Interaction *const i, KinePhaseSpace_t f, KinePhaseSpace_t t)
Definition: KineUtils.cxx:130
list x
Definition: train.py:276
const Target & Tgt(void) const
Definition: InitialState.h:66
double ProbeE(RefFrame_t rf) const
static const double kPi
Definition: Constants.h:37
static QCString * s
Definition: config.cpp:1042
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63

Member Data Documentation

double genie::SlowRsclCharmDISPXSecLO::fMc
private

Definition at line 56 of file SlowRsclCharmDISPXSecLO.h.

double genie::SlowRsclCharmDISPXSecLO::fMc2
private

Definition at line 59 of file SlowRsclCharmDISPXSecLO.h.

const PDFModelI* genie::SlowRsclCharmDISPXSecLO::fPDFModel
private

Definition at line 51 of file SlowRsclCharmDISPXSecLO.h.

double genie::SlowRsclCharmDISPXSecLO::fVcd
private

Definition at line 57 of file SlowRsclCharmDISPXSecLO.h.

double genie::SlowRsclCharmDISPXSecLO::fVcd2
private

Definition at line 60 of file SlowRsclCharmDISPXSecLO.h.

double genie::SlowRsclCharmDISPXSecLO::fVcs
private

Definition at line 58 of file SlowRsclCharmDISPXSecLO.h.

double genie::SlowRsclCharmDISPXSecLO::fVcs2
private

Definition at line 61 of file SlowRsclCharmDISPXSecLO.h.

const XSecIntegratorI* genie::SlowRsclCharmDISPXSecLO::fXSecIntegrator
private

Definition at line 52 of file SlowRsclCharmDISPXSecLO.h.


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