16#include "Framework/Conventions/GBuild.h"
59 <<
"Generating kinematics uniformly over the allowed phase space";
82 if(W.max <=0 || W.min>=W.max) {
83 LOG(
"RESKinematics",
pWARN) <<
"No available phase space";
86 exception.
SetReason(
"No available phase space");
91 const InitialState & init_state = interaction -> InitState();
103 double dW = W.max - W.min;
106 unsigned int iter = 0;
112 <<
"*** Could not select a valid (W,Q^2) pair after "
113 << iter <<
" iterations";
116 exception.
SetReason(
"Couldn't select kinematics");
130 gW = W.min + dW * rnd->
RndKine().Rndm();
133 if(Q2.max<=0. || Q2.min>=Q2.max)
continue;
134 gQ2 = Q2.min + (Q2.max-Q2.min) * rnd->
RndKine().Rndm();
155 gW = W.min + dW * rnd->
RndKine().Rndm();
156 gQD2 = QD2min + (QD2max - QD2min) * rnd->
RndKine().Rndm();
162 LOG(
"RESKinematics",
pINFO) <<
"Trying: W = " << gW <<
", Q2 = " <<
gQ2;
174 double t = xsec_max * rnd->
RndKine().Rndm();
186 <<
"Selected: W = " << gW <<
", Q2 = " <<
gQ2;
209 double totxsec = evrec->
XSec();
210 double wght = (vol/totxsec)*xsec;
211 LOG(
"RESKinematics",
pNOTICE) <<
"Kinematics wght = "<< wght;
215 LOG(
"RESKinematics",
pNOTICE) <<
"Current event wght = " << wght;
270 gROOT->GetListOfFunctions()->Remove(
fEnvelope);
284 double max_xsec = 0.;
286 const InitialState & init_state = interaction -> InitState();
306 Wmax = TMath::Min(Wmax,
fWcut);
308 Wmin = TMath::Max(Wmin, md-.3);
309 Wmax = TMath::Min(Wmax, md+.3);
311 if(Wmax-Wmin<0.05) { NW=1; Wmin=Wmax; }
313 double dW = (NW>1) ? (Wmax-Wmin)/(NW-1) : 0.;
315 for(
int iw=0; iw<NW; iw++)
317 double W = Wmin + iw*dW;
324 if( rQ2.
max < Q2Thres || rQ2.
min <=0 )
continue;
325 if( rQ2.
max-rQ2.
min<0.02 ) {NQ2=5; NQ2b=3;}
329 double dlogQ2 = (logQ2max - logQ2min) /(NQ2-1);
330 double xseclast = -1;
331 bool increasing =
true;
333 for(
int iq2=0; iq2<NQ2; iq2++)
335 double Q2 = TMath::Exp(logQ2min + iq2 * dlogQ2);
339 <<
"xsec(W= " << W <<
", Q2= " << Q2 <<
") = " << xsec;
340 max_xsec = TMath::Max(xsec, max_xsec);
341 increasing = xsec-xseclast>=0;
350 for(
int iq2b=0; iq2b<NQ2b; iq2b++)
352 Q2 = TMath::Exp(TMath::Log(Q2) - dlogQ2);
353 if(Q2 < rQ2.
min)
continue;
357 <<
"xsec(W= " << W <<
", Q2= " << Q2 <<
") = " << xsec;
358 max_xsec = TMath::Max(xsec, max_xsec);
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
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
Defines the EventGeneratorI interface.
virtual const XSecAlgorithmI * CrossSectionAlg(void) const =0
GENIE's GHEP MC event record.
virtual double Weight(void) const
virtual void SetDiffXSec(double xsec, KinePhaseSpace_t ps)
virtual double XSec(void) const
virtual Interaction * Summary(void) const
virtual TBits * EventFlags(void) const
virtual void SetWeight(double wght)
Initial State information.
const Target & Tgt(void) const
double ProbeE(RefFrame_t rf) const
Summary information for an interaction.
const XclsTag & ExclTag(void) const
const ProcessInfo & ProcInfo(void) const
const KPhaseSpace & PhaseSpace(void) const
Kinematics * KinePtr(void) const
Range1D_t Q2Lim_W(void) const
Q2 limits @ fixed W.
Range1D_t WLim(void) const
W limits.
Range1D_t Limits(KineVar_t kvar) const
Return the kinematical variable limits.
virtual void AssertXSecLimits(const Interaction *in, double xsec, double xsec_max) const
bool fGenerateUniformly
uniform over allowed phase space + event weight?
double fMaxXSecDiffTolerance
max{100*(xsec-maxxsec)/.5*(xsec+maxxsec)} if xsec>maxxsec
double fEMin
min E for which maxxsec is cached - forcing explicit calc.
double fSafetyFactor
ComputeMaxXSec -> ComputeMaxXSec * fSafetyFactor.
virtual double MaxXSec(GHepRecord *evrec, const int nkey=0) const
const XSecAlgorithmI * fXSecModel
void Setx(double x, bool selected=false)
void SetQ2(double Q2, bool selected=false)
void ClearRunningValues(void)
void Sety(double y, bool selected=false)
void SetW(double W, bool selected=false)
void ProcessEventRecord(GHepRecord *event_rec) const
void Configure(const Registry &config)
TF2 * fEnvelope
2-D envelope used for importance sampling
double fWcut
Wcut parameter in DIS/RES join scheme.
~RESKinematicsGenerator()
double ComputeMaxXSec(const Interaction *interaction) const
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)
const TLorentzVector & HitNucP4(void) const
Resonance_t Resonance(void) const
bool KnownResonance(void) const
An exception thrown by EventRecordVisitorI when the normal processing sequence has to be disrupted (f...
void SwitchOnFastForward(void)
void SetReason(string reason)
Misc GENIE control constants.
static const double kMinQ2Limit
static const unsigned int kRjMaxIterations
static const double kASmallNum
Simple functions for loading and reading nucleus dependent keys from config files.
static const double kMinQ2Limit
void WQ2toXY(double Ev, double M, double W, double Q2, double &x, double &y)
double QD2toQ2(double QD2)
double RESImportanceSamplingEnvelope(double *x, double *par)
double PhaseSpaceVolume(const Interaction *const i, KinePhaseSpace_t ps)
double Jacobian(const Interaction *const i, KinePhaseSpace_t f, KinePhaseSpace_t t)
double Q2toQD2(double Q2)
Baryon Resonance utilities.
double Mass(Resonance_t res)
resonance mass (GeV)
Root of GENIE utility namespaces.
THE MAIN GENIE PROJECT NAMESPACE
enum genie::EResonance Resonance_t
const UInt_t kISkipKinematicChk
if set, skip kinematic validity checks
const UInt_t kISkipProcessChk
if set, skip process validity checks