GENIEGenerator
Loading...
Searching...
No Matches
QPMDISPXSec.cxx
Go to the documentation of this file.
1//____________________________________________________________________________
2/*
3 Copyright (c) 2003-2025, The GENIE Collaboration
4 For the full text of the license visit http://copyright.genie-mc.org
5
6 Costas Andreopoulos <c.andreopoulos \at cern.ch>
7 University of Liverpool
8*/
9//____________________________________________________________________________
10
11#include <sstream>
12
13#include <TMath.h>
14#include <TH1D.h>
15
20#include "Framework/Conventions/GBuild.h"
36
37using std::ostringstream;
38
39using namespace genie;
40using namespace genie::constants;
41//using namespace genie::units;
42
43//____________________________________________________________________________
45XSecAlgorithmI("genie::QPMDISPXSec")
46{
47 fInInitPhase = true;
48}
49//____________________________________________________________________________
51XSecAlgorithmI("genie::QPMDISPXSec", config)
52{
53 fInInitPhase = true;
54}
55//____________________________________________________________________________
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 if( proc_info.IsWeakCC() ) xsec *= fCCScale;
175 else if( proc_info.IsWeakNC() ) xsec *= fEMScale;
176 else if( proc_info.IsEM() ) xsec *= fEMScale;
177
178 // Subtract the inclusive charm production cross section
179 interaction->ExclTagPtr()->SetCharm();
180 double xsec_charm = fCharmProdModel->XSec(interaction,kps);
181 interaction->ExclTagPtr()->UnsetCharm();
182#ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
183 LOG("DISPXSec", pINFO)
184 << "Subtracting charm piece: " << xsec_charm << " / out of " << xsec;
185#endif
186 xsec = TMath::Max(0., xsec-xsec_charm);
187 return xsec;
188}
189//____________________________________________________________________________
190double QPMDISPXSec::Integral(const Interaction * interaction) const
191{
192 double xsec = fXSecIntegrator->Integrate(this,interaction);
193 return xsec;
194}
195//____________________________________________________________________________
196bool QPMDISPXSec::ValidProcess(const Interaction * interaction) const
197{
198 if(interaction->TestBit(kISkipProcessChk)) return true;
199
200 const ProcessInfo & proc_info = interaction->ProcInfo();
201 if(!proc_info.IsDeepInelastic()) return false;
202
203 const InitialState & init_state = interaction -> InitState();
204 int probe_pdg = init_state.ProbePdg();
205 if(!pdg::IsLepton(probe_pdg)) return false;
206
207 if(! init_state.Tgt().HitNucIsSet()) return false;
208
209 int hitnuc_pdg = init_state.Tgt().HitNucPdg();
210 if(!pdg::IsNeutronOrProton(hitnuc_pdg)) return false;
211
212 return true;
213}
214//____________________________________________________________________________
216{
217 Algorithm::Configure(config);
218 this->LoadConfig();
219}
220//____________________________________________________________________________
221void QPMDISPXSec::Configure(string config)
222{
223 Algorithm::Configure(config);
224
225 Registry r( "QPMDISPXSec_specific", false ) ;
226
227 RgKey xdefkey = "XSecModel@genie::EventGenerator/DIS-CC-CHARM";
228 RgKey local_key = "CharmXSec" ;
229 r.Set( local_key, AlgConfigPool::Instance() -> GlobalParameterList() -> GetAlg(xdefkey) ) ;
230
232
233 this->LoadConfig();
234}
235//____________________________________________________________________________
237{
238 // Access global defaults to use in case of missing parameters
239
240 fDISSFModel = 0;
242 dynamic_cast<const DISStructureFuncModelI *> (this->SubAlg("SFAlg"));
243 assert(fDISSFModel);
244
245 fDISSF.SetModel(fDISSFModel); // <-- attach algorithm
246
247 // Cross section scaling factor
248 GetParam( "DIS-CC-XSecScale", fCCScale ) ;
249 GetParam( "DIS-NC-XSecScale", fNCScale ) ;
250 GetParam( "DIS-EM-XSecScale", fEMScale ) ;
251
252 // sin^4(theta_weinberg)
253 double thw ;
254 GetParam( "WeinbergAngle", thw ) ;
255 fSin48w = TMath::Power( TMath::Sin(thw), 4 );
256
257
258 // Since this method would be called every time the current algorithm is
259 // reconfigured at run-time, remove all the data cached by this algorithm
260 // since they depend on the previous configuration
261
262 if(!fInInitPhase) {
263 Cache * cache = Cache::Instance();
264 string keysubstr = this->Id().Key() + "/DIS-RES-Join";
265 cache->RmMatchedCacheBranches(keysubstr);
266 }
267 fInInitPhase = false;
268
269 //-- load the differential cross section integrator
271 dynamic_cast<const XSecIntegratorI *> (this->SubAlg("XSec-Integrator"));
272 assert(fXSecIntegrator);
273
274 // Load the charm production cross section model
275 RgKey local_key = "CharmXSec" ;
276 RgAlg xalg;
277 GetParam( local_key, xalg ) ;
278 LOG("DISXSec", pDEBUG)
279 << "Loading the cross section model: " << xalg;
280
281 fCharmProdModel = dynamic_cast<const XSecAlgorithmI *> ( this -> SubAlg(local_key) ) ;
282 assert(fCharmProdModel);
283}
284//____________________________________________________________________________
#define pINFO
Definition Messenger.h:62
#define pDEBUG
Definition Messenger.h:63
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
Definition Messenger.h:96
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils.
string RgKey
static AlgConfigPool * Instance()
string Key(void) const
Definition AlgId.h:46
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
virtual void Configure(const Registry &config)
Definition Algorithm.cxx:62
const Algorithm * SubAlg(const RgKey &registry_key) const
virtual const AlgId & Id(void) const
Get algorithm ID.
Definition Algorithm.h:98
GENIE Cache Memory.
Definition Cache.h:39
static Cache * Instance(void)
Definition Cache.cxx:67
void RmMatchedCacheBranches(string key_substring)
Definition Cache.cxx:127
Pure Abstract Base Class. Defines the DISStructureFuncModelI interface to be implemented by any algor...
Initial State information.
const Target & Tgt(void) const
int ProbePdg(void) const
double ProbeE(RefFrame_t rf) const
Summary information for an interaction.
Definition Interaction.h:56
XclsTag * ExclTagPtr(void) const
Definition Interaction.h:77
const ProcessInfo & ProcInfo(void) const
Definition Interaction.h:70
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
Generated/set kinematical variables for an event.
Definition Kinematics.h:39
double y(bool selected=false) const
double x(bool selected=false) const
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition ProcessInfo.h:46
bool IsWeakNC(void) const
bool IsDeepInelastic(void) const
bool IsWeakCC(void) const
bool IsEM(void) const
DISStructureFunc fDISSF
Definition QPMDISPXSec.h:53
double fEMScale
cross section scaling factor
Definition QPMDISPXSec.h:63
double XSec(const Interaction *i, KinePhaseSpace_t k) const
Compute the cross section for the input interaction.
const XSecAlgorithmI * fCharmProdModel
Definition QPMDISPXSec.h:59
double Integral(const Interaction *i) const
double fCCScale
cross section scaling factor
Definition QPMDISPXSec.h:61
bool ValidProcess(const Interaction *i) const
Can this cross section algorithm handle the input process?
double fNCScale
cross section scaling factor
Definition QPMDISPXSec.h:62
void Configure(const Registry &config)
const DISStructureFuncModelI * fDISSFModel
SF model.
Definition QPMDISPXSec.h:56
double fSin48w
sin^4(Weingberg angle)
Definition QPMDISPXSec.h:64
const XSecIntegratorI * fXSecIntegrator
diff. xsec integrator
Definition QPMDISPXSec.h:57
A registry. Provides the container for algorithm configuration parameters.
Definition Registry.h:65
void Set(RgIMapPair entry)
Definition Registry.cxx:267
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition Target.h:40
int HitNucPdg(void) const
Definition Target.cxx:304
int N(void) const
Definition Target.h:69
int Z(void) const
Definition Target.h:68
double HitNucMass(void) const
Definition Target.cxx:233
bool HitNucIsSet(void) const
Definition Target.cxx:283
virtual bool ValidKinematics(const Interaction *i) const
Is the input kinematical point a physically allowed one?
Cross Section Integrator Interface.
void SetCharm(int charm_pdgc=0)
Definition XclsTag.cxx:59
void UnsetCharm(void)
Definition XclsTag.cxx:65
Basic constants.
bool IsLepton(int pdgc)
Definition PDGUtils.cxx:86
bool IsProton(int pdgc)
Definition PDGUtils.cxx:336
bool IsNeutronOrProton(int pdgc)
Definition PDGUtils.cxx:351
bool IsAntiNeutrino(int pdgc)
Definition PDGUtils.cxx:118
double Jacobian(const Interaction *const i, KinePhaseSpace_t f, KinePhaseSpace_t t)
double XYtoQ2(double Ev, double M, double x, double y)
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25
enum genie::EKinePhaseSpace KinePhaseSpace_t
@ kRfHitNucRest
Definition RefFrame.h:30
const UInt_t kISkipProcessChk
if set, skip process validity checks
Definition Interaction.h:47
const UInt_t kIAssumeFreeNucleon
Definition Interaction.h:49