H3AMNuGammaPXSec.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 
12 #include "Framework/Conventions/GBuild.h"
19 
20 using namespace genie;
21 using namespace genie::constants;
22 
23 //____________________________________________________________________________
25 XSecAlgorithmI("genie::H3AMNuGammaPXSec")
26 {
27 
28 }
29 //____________________________________________________________________________
31 XSecAlgorithmI("genie::H3AMNuGammaPXSec", config)
32 {
33 
34 }
35 //____________________________________________________________________________
37 {
38 
39 }
40 //____________________________________________________________________________
42  const Interaction * interaction, KinePhaseSpace_t /*kps*/) const
43 {
44  if(! this -> ValidProcess (interaction) ) return 0.;
45  if(! this -> ValidKinematics (interaction) ) return 0.;
46 
47  LOG("AMNuGamma", pWARN)
48  << "*** No differential cross section calculation is implemented yet";
49 
50  return 1;
51 }
52 //____________________________________________________________________________
54 {
55  // Compute the cross section for a free nucleon target
56  const InitialState & init_state = interaction -> InitState();
57  const Target & target = init_state.Tgt();
58  double Ev = init_state.ProbeE(kRfHitNucRest);
59 
60  double Ecutoff = kNucleonMass / 2;
61 
62  if(Ev > Ecutoff) return 0;
63 
64  double xsec0 = 2.2E-41 * units::cm2;
65  double xsec = xsec0 * TMath::Power(Ev,6.) * TMath::Power(0.1*fGw,4.);
66 
67  LOG("AMNuGamma", pNOTICE)
68  << "*** xsec(vN->vNgamma) [free nuc](Ev="<< Ev << ") = "<< xsec;
69 
70 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
71  LOG("AMNuGamma", pDEBUG)
72  << "*** xsec(vN->vNgamma) [free nuc](Ev="<< Ev << ") = "<< xsec;
73 #endif
74 
75  // If requested return the free xsec even for nuclear target
76  if( interaction->TestBit(kIAssumeFreeNucleon) ) return xsec;
77 
78  // Scale for the number of scattering centers at the target
79  int nucpdgc = target.HitNucPdg();
80  int NNucl = (pdg::IsProton(nucpdgc)) ? target.Z() : target.N();
81  xsec*=NNucl;
82 
83  return xsec;
84 }
85 //____________________________________________________________________________
87 {
88  if(interaction->TestBit(kISkipProcessChk)) return true;
89 
90  if(interaction->ProcInfo().IsAMNuGamma()) return true;
91  return false;
92 }
93 //____________________________________________________________________________
95 {
96  if(interaction->TestBit(kISkipKinematicChk)) return true;
97  return true;
98 }
99 //____________________________________________________________________________
101 {
102  Algorithm::Configure(config);
103  this->LoadConfig();
104 }
105 //____________________________________________________________________________
107 {
108  Algorithm::Configure(config);
109  this->LoadConfig();
110 }
111 //____________________________________________________________________________
113 {
114 
115  GetParam( "AMNuGamma-Gw", fGw ) ;
116 
117 }
118 //____________________________________________________________________________
Cross Section Calculation Interface.
Basic constants.
THE MAIN GENIE PROJECT NAMESPACE
Definition: AlgCmp.h:25
static const double kNucleonMass
Definition: Constants.h:77
int HitNucPdg(void) const
Definition: Target.cxx:304
enum genie::EKinePhaseSpace KinePhaseSpace_t
Summary information for an interaction.
Definition: Interaction.h:56
bool IsProton(int pdgc)
Definition: PDGUtils.cxx:333
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool ValidKinematics(const Interaction *i) const
Is the input kinematical point a physically allowed one?
static constexpr double cm2
Definition: Units.h:69
static Config * config
Definition: config.cpp:1054
bool IsAMNuGamma(void) const
void Configure(const Registry &config)
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 Z(void) const
Definition: Target.h:68
#define pWARN
Definition: Messenger.h:60
const UInt_t kISkipKinematicChk
if set, skip kinematic validity checks
Definition: Interaction.h:48
bool ValidProcess(const Interaction *i) const
Can this cross section algorithm handle the input process?
int N(void) const
Definition: Target.h:69
double XSec(const Interaction *i, KinePhaseSpace_t k) const
Compute the cross section for the input interaction.
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:65
const UInt_t kIAssumeFreeNucleon
Definition: Interaction.h:49
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
#define pNOTICE
Definition: Messenger.h:61
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
const Target & Tgt(void) const
Definition: InitialState.h:66
double Integral(const Interaction *i) const
double ProbeE(RefFrame_t rf) const
const UInt_t kISkipProcessChk
if set, skip process validity checks
Definition: Interaction.h:47
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63