QPMDISPXSec.cxx
Go to the documentation of this file.
1 //____________________________________________________________________________
2 /*
3  Copyright (c) 2003-2020, The GENIE Collaboration
4  For the full text of the license visit http://copyright.genie-mc.org
5 
6  Costas Andreopoulos <constantinos.andreopoulos \at cern.ch>
7  University of Liverpool & STFC Rutherford Appleton Laboratory
8 */
9 //____________________________________________________________________________
10 
11 #include <sstream>
12 
13 #include <TMath.h>
14 #include <TH1D.h>
15 
20 #include "Framework/Conventions/GBuild.h"
29 #include "Framework/Utils/RunOpt.h"
32 #include "Framework/Utils/Range1.h"
34 #include "Framework/Utils/Cache.h"
36 
37 using std::ostringstream;
38 
39 using namespace genie;
40 using namespace genie::constants;
41 //using namespace genie::units;
42 
43 //____________________________________________________________________________
45 XSecAlgorithmI("genie::QPMDISPXSec")
46 {
47  fInInitPhase = true;
48 }
49 //____________________________________________________________________________
51 XSecAlgorithmI("genie::QPMDISPXSec", config)
52 {
53  fInInitPhase = true;
54 }
55 //____________________________________________________________________________
57 {
58 
59 }
60 //____________________________________________________________________________
62  const Interaction * interaction, KinePhaseSpace_t kps) const
63 {
64  if(! this -> ValidProcess (interaction) ) return 0.;
65  if(! this -> ValidKinematics (interaction) ) return 0.;
66 
67  // Get kinematical & init-state parameters
68  const Kinematics & kinematics = interaction -> Kine();
69  const InitialState & init_state = interaction -> InitState();
70  const ProcessInfo & proc_info = interaction -> ProcInfo();
71 
72  double E = init_state.ProbeE(kRfHitNucRest);
73  double ml = interaction->FSPrimLepton()->Mass();
74  double Mnuc = init_state.Tgt().HitNucMass();
75  double x = kinematics.x();
76  double y = kinematics.y();
77 
78  double E2 = E * E;
79  double ml2 = ml * ml;
80  double ml4 = ml2 * ml2;
81  double Mnuc2 = Mnuc * Mnuc;
82 
83 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
84  LOG("DISPXSec", pDEBUG)
85  << "Computing d2xsec/dxdy @ E = " << E << ", x = " << x << ", y = " << y;
86 #endif
87 
88  // One of the xsec terms changes sign for antineutrinos @ DIS/CC
89 
90  bool is_nubar_cc = pdg::IsAntiNeutrino(init_state.ProbePdg()) &&
91  proc_info.IsWeakCC();
92  int sign = (is_nubar_cc) ? -1 : 1;
93 
94  // Calculate the DIS structure functions
95  fDISSF.Calculate(interaction);
96 
97 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
98  LOG("DISPXSec", pDEBUG) << fDISSF;
99 #endif
100 
101  //
102  // Compute the differential cross section
103  //
104 
105  double g2 = kGF2;
106  // For EM interaction replace G_{Fermi} with :
107  // a_{em} * pi / ( sqrt(2) * sin^2(theta_weinberg) * Mass_{W}^2 }
108  // See C.Quigg, Gauge Theories of the Strong, Weak and E/M Interactions,
109  // ISBN 0-8053-6021-2, p.112 (6.3.57)
110  // Also, take int account that the photon propagator is 1/p^2 but the
111  // W propagator is 1/(p^2-Mass_{W}^2), so weight the EM case with
112  // Mass_{W}^4 / q^4
113  // So, overall:
114  // G_{Fermi}^2 --> a_{em}^2 * pi^2 / (2 * sin^4(theta_weinberg) * q^{4})
115  //
116  double Q2 = utils::kinematics::XYtoQ2(E,Mnuc,x,y);
117  double Q4 = Q2*Q2;
118  if(proc_info.IsEM()) {
119  g2 = kAem2 * kPi2 / (2.0 * fSin48w * Q4);
120  }
121  if (proc_info.IsWeakCC()) {
122  g2 = kGF2 * kMw2 * kMw2 / TMath::Power((Q2 + kMw2), 2);
123  } else if (proc_info.IsWeakNC()) {
124  g2 = kGF2 * kMz2 * kMz2 / TMath::Power((Q2 + kMz2), 2);
125  }
126  double front_factor = (g2*Mnuc*E) / kPi;
127 
128  // Build all dxsec/dxdy terms
129  double term1 = y * ( x*y + ml2/(2*E*Mnuc) );
130  double term2 = 1 - y - Mnuc*x*y/(2*E) - ml2/(4*E2);
131  double term3 = sign * (x*y*(1-y/2) - y*ml2/(4*Mnuc*E));
132  double term4 = x*y*ml2/(2*Mnuc*E) + ml4/(4*Mnuc2*E2);
133  double term5 = -1.*ml2/(2*Mnuc*E);
134 
135 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
136  LOG("DISPXSec", pDEBUG)
137  << "\nd2xsec/dxdy ~ (" << term1 << ")*F1+(" << term2 << ")*F2+("
138  << term3 << ")*F3+(" << term4 << ")*F4+(" << term5 << ")*F5";
139 #endif
140 
141  term1 *= fDISSF.F1();
142  term2 *= fDISSF.F2();
143  term3 *= fDISSF.F3();
144  term4 *= fDISSF.F4();
145  term5 *= fDISSF.F5();
146 
147  double xsec = front_factor * (term1 + term2 + term3 + term4 + term5);
148  xsec = TMath::Max(xsec,0.);
149 
150 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
151  LOG("DISPXSec", pINFO)
152  << "d2xsec/dxdy[FreeN] (E= " << E
153  << ", x= " << x << ", y= " << y << ") = " << xsec;
154 #endif
155 
156  // The algorithm computes d^2xsec/dxdy
157  // Check whether variable tranformation is needed
158  if(kps!=kPSxyfE) {
159  double J = utils::kinematics::Jacobian(interaction,kPSxyfE,kps);
160  xsec *= J;
161  }
162 
163  // If requested return the free nucleon xsec even for input nuclear tgt
164  if( interaction->TestBit(kIAssumeFreeNucleon) ) return xsec;
165 
166  // Compute nuclear cross section (simple scaling here, corrections must
167  // have been included in the structure functions)
168  const Target & target = init_state.Tgt();
169  int nucpdgc = target.HitNucPdg();
170  int NNucl = (pdg::IsProton(nucpdgc)) ? target.Z() : target.N();
171  xsec *= NNucl;
172 
173  // Apply scaling / if required to reach well known asymmptotic value
174  xsec *= fScale;
175 
176  // Subtract the inclusive charm production cross section
177  interaction->ExclTagPtr()->SetCharm();
178  double xsec_charm = fCharmProdModel->XSec(interaction,kps);
179  interaction->ExclTagPtr()->UnsetCharm();
180 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
181  LOG("DISPXSec", pINFO)
182  << "Subtracting charm piece: " << xsec_charm << " / out of " << xsec;
183 #endif
184  xsec = TMath::Max(0., xsec-xsec_charm);
185  return xsec;
186 }
187 //____________________________________________________________________________
189 {
190  double xsec = fXSecIntegrator->Integrate(this,interaction);
191  return xsec;
192 }
193 //____________________________________________________________________________
195 {
196  if(interaction->TestBit(kISkipProcessChk)) return true;
197 
198  const ProcessInfo & proc_info = interaction->ProcInfo();
199  if(!proc_info.IsDeepInelastic()) return false;
200 
201  const InitialState & init_state = interaction -> InitState();
202  int probe_pdg = init_state.ProbePdg();
203  if(!pdg::IsLepton(probe_pdg)) return false;
204 
205  if(! init_state.Tgt().HitNucIsSet()) return false;
206 
207  int hitnuc_pdg = init_state.Tgt().HitNucPdg();
208  if(!pdg::IsNeutronOrProton(hitnuc_pdg)) return false;
209 
210  return true;
211 }
212 //____________________________________________________________________________
214 {
215  Algorithm::Configure(config);
216  this->LoadConfig();
217 }
218 //____________________________________________________________________________
220 {
221  Algorithm::Configure(config);
222 
223  Registry r( "QPMDISPXSec_specific", false ) ;
224 
225  RgKey xdefkey = "XSecModel@genie::EventGenerator/DIS-CC-CHARM";
226  RgKey local_key = "CharmXSec" ;
227  r.Set( local_key, AlgConfigPool::Instance() -> GlobalParameterList() -> GetAlg(xdefkey) ) ;
228 
230 
231  this->LoadConfig();
232 }
233 //____________________________________________________________________________
235 {
236  // Access global defaults to use in case of missing parameters
237 
238  fDISSFModel = 0;
239  fDISSFModel =
240  dynamic_cast<const DISStructureFuncModelI *> (this->SubAlg("SFAlg"));
241  assert(fDISSFModel);
242 
243  fDISSF.SetModel(fDISSFModel); // <-- attach algorithm
244 
245  // Cross section scaling factor
246  GetParam( "DIS-XSecScale", fScale ) ;
247 
248  // sin^4(theta_weinberg)
249  double thw ;
250  GetParam( "WeinbergAngle", thw ) ;
251  fSin48w = TMath::Power( TMath::Sin(thw), 4 );
252 
253 
254  // Since this method would be called every time the current algorithm is
255  // reconfigured at run-time, remove all the data cached by this algorithm
256  // since they depend on the previous configuration
257 
258  if(!fInInitPhase) {
259  Cache * cache = Cache::Instance();
260  string keysubstr = this->Id().Key() + "/DIS-RES-Join";
261  cache->RmMatchedCacheBranches(keysubstr);
262  }
263  fInInitPhase = false;
264 
265  //-- load the differential cross section integrator
267  dynamic_cast<const XSecIntegratorI *> (this->SubAlg("XSec-Integrator"));
268  assert(fXSecIntegrator);
269 
270  // Load the charm production cross section model
271  RgKey local_key = "CharmXSec" ;
272  RgAlg xalg;
273  GetParam( local_key, xalg ) ;
274  LOG("DISXSec", pDEBUG)
275  << "Loading the cross section model: " << xalg;
276 
277  fCharmProdModel = dynamic_cast<const XSecAlgorithmI *> ( this -> SubAlg(local_key) ) ;
278  assert(fCharmProdModel);
279 }
280 //____________________________________________________________________________
void SetModel(const DISStructureFuncModelI *model)
Attach an algorithm.
Cross Section Calculation Interface.
Pure Abstract Base Class. Defines the DISStructureFuncModelI interface to be implemented by any algor...
Basic constants.
bool IsWeakCC(void) const
double F2(void) const
Get the computed structure function F2.
static const double kMw2
Definition: Constants.h:93
double J(double q0, double q3, double Enu, double ml)
Definition: MECUtils.cxx:147
THE MAIN GENIE PROJECT NAMESPACE
Definition: AlgCmp.h:25
Cross Section Integrator Interface.
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1064
int HitNucPdg(void) const
Definition: Target.cxx:304
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
const DISStructureFuncModelI * fDISSFModel
SF model.
Definition: QPMDISPXSec.h:56
double HitNucMass(void) const
Definition: Target.cxx:233
static const double kMz2
Definition: Constants.h:94
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
double x(bool selected=false) const
Definition: Kinematics.cxx:99
enum genie::EKinePhaseSpace KinePhaseSpace_t
double fSin48w
sin^4(Weingberg angle)
Definition: QPMDISPXSec.h:62
void SetCharm(int charm_pdgc=0)
Definition: XclsTag.cxx:59
double y(bool selected=false) const
Definition: Kinematics.cxx:112
double F4(void) const
Get the computed structure function F4.
Summary information for an interaction.
Definition: Interaction.h:56
double F1(void) const
Get the computed structure function F1.
virtual bool ValidKinematics(const Interaction *i) const
Is the input kinematical point a physically allowed one?
bool IsProton(int pdgc)
Definition: PDGUtils.cxx:333
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
const XSecIntegratorI * fXSecIntegrator
diff. xsec integrator
Definition: QPMDISPXSec.h:57
static Config * config
Definition: config.cpp:1054
static const double kAem2
Definition: Constants.h:57
void Configure(const Registry &config)
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
double XYtoQ2(double Ev, double M, double x, double y)
Definition: KineUtils.cxx:1195
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
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
int ProbePdg(void) const
Definition: InitialState.h:64
double fScale
cross section scaling factor
Definition: QPMDISPXSec.h:61
void LoadConfig(void)
int Z(void) const
Definition: Target.h:68
#define pINFO
Definition: Messenger.h:62
DISStructureFunc fDISSF
Definition: QPMDISPXSec.h:53
bool ValidProcess(const Interaction *i) const
Can this cross section algorithm handle the input process?
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
double F5(void) const
Get the computed structure function F5.
bool IsEM(void) const
GENIE Cache Memory.
Definition: Cache.h:38
bool IsDeepInelastic(void) const
Definition: ProcessInfo.cxx:84
XclsTag * ExclTagPtr(void) const
Definition: Interaction.h:77
const XSecAlgorithmI * fCharmProdModel
Definition: QPMDISPXSec.h:59
int sign(double val)
Definition: UtilFunc.cxx:104
double F3(void) const
Get the computed structure function F3.
virtual const AlgId & Id(void) const
Get algorithm ID.
Definition: Algorithm.h:97
int N(void) const
Definition: Target.h:69
bool HitNucIsSet(void) const
Definition: Target.cxx:283
string RgKey
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:65
const UInt_t kIAssumeFreeNucleon
Definition: Interaction.h:49
E
Definition: 018_def.c:13
bool IsNeutronOrProton(int pdgc)
Definition: PDGUtils.cxx:348
double XSec(const Interaction *i, KinePhaseSpace_t k) const
Compute the cross section for the input interaction.
Definition: QPMDISPXSec.cxx:61
Definition: 018_def.c:13
virtual double Integrate(const XSecAlgorithmI *model, const Interaction *interaction) const =0
void UnsetCharm(void)
Definition: XclsTag.cxx:65
double Jacobian(const Interaction *const i, KinePhaseSpace_t f, KinePhaseSpace_t t)
Definition: KineUtils.cxx:130
virtual ~QPMDISPXSec()
Definition: QPMDISPXSec.cxx:56
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
list x
Definition: train.py:276
void Calculate(const Interaction *interaction)
Calculate the S/F&#39;s for the input interaction using the attached algorithm.
void RmMatchedCacheBranches(string key_substring)
Definition: Cache.cxx:127
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
const Target & Tgt(void) const
Definition: InitialState.h:66
static const double kGF2
Definition: Constants.h:59
static Cache * Instance(void)
Definition: Cache.cxx:67
double Integral(const Interaction *i) const
void Set(RgIMapPair entry)
Definition: Registry.cxx:267
double ProbeE(RefFrame_t rf) const
static const double kPi
Definition: Constants.h:37
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils...
string Key(void) const
Definition: AlgId.h:46
bool IsLepton(int pdgc)
Definition: PDGUtils.cxx:83
static const double kPi2
Definition: Constants.h:38
const UInt_t kISkipProcessChk
if set, skip process validity checks
Definition: Interaction.h:47
static AlgConfigPool * Instance()
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63
const Algorithm * SubAlg(const RgKey &registry_key) const
Definition: Algorithm.cxx:345