15#include <Math/IFunction.h>
16#include <Math/GSLMinimizer1D.h>
17#include <Math/BrentMinimizer1D.h>
80 const InitialState & init_state = interaction -> InitState();
94 <<
"Option to generate kinematics uniformly not supported";
104 double dDNuE = DNuEnergy.
max - DNuEnergy.
min;
105 ROOT::Math::BrentMinimizer1D minimizer;
106 minimizer.SetFunction( xsec_func, DNuEnergy.
min, DNuEnergy.
max);
107 minimizer.Minimize(1000, 1, 1E-5);
108 double DNuE_for_xsec_max = minimizer.XMinimum();
109 double xsec_max = -1. * xsec_func(DNuE_for_xsec_max);
112 unsigned int iter = 0;
117 <<
"*** Could not select a valid DNuE after " << iter <<
" iterations";
118 event->EventFlags()->SetBitNumber(
kKineGenErr,
true);
120 exception.
SetReason(
"Couldn't select kinematics");
125 gDNuE = DNuEnergy.
min + dDNuE * rnd->
RndKine().Rndm();
126 LOG(
"COHDNu",
pINFO) <<
"Trying: E_N = " << gDNuE;
129 gxsec = -1. * xsec_func(gDNuE);
131 if(gxsec > xsec_max) {
132 double frac = TMath::Abs(gxsec-xsec_max)/xsec_max;
135 <<
"Current computed cross-section (" << gxsec/(
units::cm2)
136 <<
" cm2/GeV^2) exceeds the maximum cross-section ("
137 << xsec_max/(
units::cm2) <<
" beyond the specified tolerance";
145 <<
"dxsec/dQ2 = " << gxsec/(
units::cm2) <<
" cm2/GeV^2"
146 <<
"J = 1, rnd = " << t;
147 bool accept = (t<gxsec);
158 double ETimesM = E_nu * M_target;
159 double EPlusM = E_nu + M_target;
161 double p_DNu = TMath::Sqrt(gDNuE*gDNuE -
fDNuMass2);
162 double cos_theta_DNu = (gDNuE*EPlusM - ETimesM - 0.5*
fDNuMass2) / (E_nu * p_DNu);
163 double theta_DNu = TMath::ACos(cos_theta_DNu);
167 TVector3 unit_nudir = probe->
P4()->Vect().Unit();
169 TVector3 DNu_3vector = TVector3(0,0,0);
171 DNu_3vector.SetMagThetaPhi(p_DNu, theta_DNu, phi);
172 DNu_3vector.RotateUz(unit_nudir);
173 TLorentzVector P4_DNu = TLorentzVector(DNu_3vector, gDNuE);
177 double gQ2 = -(P4_DNu - *(probe->
P4())).M2();
189 const TLorentzVector & vtx = *(probe->
X4());
190 TLorentzVector x4l(vtx);
204 const TLorentzVector & p4probe = * ( probe -> P4() );
205 const TLorentzVector & p4target = * ( target -> P4() );
206 const TLorentzVector & p4fsl = * ( fsl -> P4() );
208 const TLorentzVector & p4recoil = p4probe + p4target - p4fsl;
213 const TLorentzVector & vtx = *(probe->
X4());
218 -1,-1, p4recoil, vtx);
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils.
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
virtual void Configure(const Registry &config)
bool GetParamDef(const RgKey &name, T &p, const T &def) const
void ProcessEventRecord(GHepRecord *event) const
double fMaxXSecDiffTolerance
void GenerateKinematics(GHepRecord *event) const
void AddFinalStateDarkNeutrino(GHepRecord *event) const
void Configure(const Registry &config)
void AddRecoilNucleus(GHepRecord *event) const
const XSecAlgorithmI * fXSecModel
Defines the EventGeneratorI interface.
virtual const XSecAlgorithmI * CrossSectionAlg(void) const =0
STDHEP-like event record entry that can fit a particle or a nucleus.
const TLorentzVector * P4(void) const
const TLorentzVector * X4(void) const
int FirstDaughter(void) const
GENIE's GHEP MC event record.
virtual int ProbePosition(void) const
virtual GHepParticle * TargetNucleus(void) const
virtual Interaction * Summary(void) const
virtual int TargetNucleusPosition(void) const
Initial State information.
const Target & Tgt(void) const
double ProbeE(RefFrame_t rf) const
Summary information for an interaction.
const Kinematics & Kine(void) const
const InitialState & InitState(void) const
Kinematics * KinePtr(void) const
void SetQ2(double Q2, bool selected=false)
const TLorentzVector & FSLeptonP4(void) const
void ClearRunningValues(void)
void SetFSLeptonP4(const TLorentzVector &p4)
A singleton holding random number generator classes. All random number generation in GENIE should tak...
static RandomGen * Instance()
Access instance.
TRandom3 & RndKine(void) const
rnd number generator used by kinematics generators
A simple [min,max] interval for doubles.
A registry. Provides the container for algorithm configuration parameters.
Keep info on the event generation thread currently on charge. This is used so that event generation m...
static RunningThreadInfo * Instance(void)
const EventGeneratorI * RunningThread(void)
An exception thrown by EventRecordVisitorI when the normal processing sequence has to be disrupted (f...
void SwitchOnFastForward(void)
void SetReason(string reason)
Range1D_t IntegrationRange(void) const
Misc GENIE control constants.
static const unsigned int kRjMaxIterations
static constexpr double cm2
Simple functions for loading and reading nucleus dependent keys from config files.
string P4AsString(const TLorentzVector *p)
Root of GENIE utility namespaces.
THE MAIN GENIE PROJECT NAMESPACE
const int kPdgDarkNeutrino
const int kPdgAntiDarkNeutrino
const UInt_t kISkipKinematicChk
if set, skip kinematic validity checks
const UInt_t kISkipProcessChk
if set, skip process validity checks