16 #include "Math/Minimizer.h" 17 #include "Math/Factory.h" 20 #include "Framework/Conventions/GBuild.h" 36 using namespace genie;
63 <<
"Generating kinematics uniformly over the allowed phase space";
74 }
else if (
fXSecModel->
Id().
Name() ==
"genie::BergerSehgalFMCOHPiPXSec2015") {
81 "ProcessEventRecord >> Cannot calculate kinematics for " <<
102 double xsec_max = this->
MaxXSec(evrec);
112 const double dy = ymax - ymin;
115 const double dQ2 = Q2max - Q2min;
117 unsigned int iter = 0;
119 double xsec=-1, gy=-1,
gQ2=-1;
129 gy = ymin + dy * rnd->
RndKine().Rndm();
133 "Trying: Q^2 = " <<
gQ2 <<
", y = " << gy;
143 accept = (xsec_max * rnd->
RndKine().Rndm() < xsec);
148 <<
"Selected: Q^2 = " <<
gQ2 <<
", y = " << gy;
159 double Epi2 = TMath::Power(Epi,2);
161 double gx = interaction->
KinePtr()->
x();
163 double tA = 1. + xME - 0.5*pme2;
165 double tB = TMath::Sqrt(1.+ 2*xME) * TMath::Sqrt(1.-pme2);
166 double tmin = 2*Epi2 * (tA-tB);
167 double tmax = 2*Epi2 * (tA+tB);
168 double A = (double) init_state.
Tgt().
A();
169 double A13 = TMath::Power(A,1./3.);
171 double R2 = TMath::Power(R,2.);
172 double b = 0.33333 * R2;
173 double tsum = (TMath::Exp(-b*tmin) - TMath::Exp(-b*tmax))/b;
174 double rt = tsum * rnd->
RndKine().Rndm();
175 double gt = -1.*TMath::Log(-1.*b*rt + TMath::Exp(-1.*b*tmin))/
b;
216 double xsec_max = this->
MaxXSec(evrec);
226 const double dy = ymax - ymin;
229 const double dQ2 = Q2max - Q2min;
232 const double dt = tmax - tmin;
236 unsigned int iter = 0;
238 double xsec=-1, gy=-1,
gt=-1,
gQ2=-1;
248 gy = ymin + dy * rnd->
RndKine().Rndm();
253 "Trying: Q^2 = " <<
gQ2 <<
", y = " << gy <<
", t = " <<
gt;
263 accept = (xsec_max * rnd->
RndKine().Rndm() < xsec);
268 <<
"Selected: Q^2 = " <<
gQ2 <<
", y = " << gy <<
", t = " <<
gt;
319 const double dx = xmax - xmin;
320 const double dy = ymax - ymin;
324 unsigned int iter = 0;
326 double xsec=-1, gx=-1, gy=-1;
334 gx = xmin + dx * rnd->
RndKine().Rndm();
335 gy = ymin + dy * rnd->
RndKine().Rndm();
341 LOG(
"COHKinematics",
pNOTICE) <<
"Initializing the sampling envelope";
343 fEnvelope->SetRange(xmin,ymin,xmax,ymax);
352 LOG(
"COHKinematics",
pINFO) <<
"Trying: x = " << gx <<
", y = " << gy;
363 double t = max * rnd->
RndKine().Rndm();
366 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 368 <<
"xsec= " << xsec <<
", J= 1, Rnd= " <<
t;
378 LOG(
"COHKinematics",
pNOTICE) <<
"Selected: x = "<< gx <<
", y = "<< gy;
389 double Epi2 = TMath::Power(Epi,2);
392 double tA = 1. + xME - 0.5*pme2;
393 double tB = TMath::Sqrt(1.+ 2*xME) * TMath::Sqrt(1.-pme2);
394 double tmin = 2*Epi2 * (tA-tB);
395 double tmax = 2*Epi2 * (tA+tB);
396 double A = (double) init_state.
Tgt().
A();
397 double A13 = TMath::Power(A,1./3.);
399 double R2 = TMath::Power(R,2.);
400 double b = 0.33333 * R2;
401 double tsum = (TMath::Exp(-b*tmin) - TMath::Exp(-b*tmax))/b;
402 double rt = tsum * rnd->
RndKine().Rndm();
403 double gt = -1.*TMath::Log(-1.*b*rt + TMath::Exp(-1.*b*tmin))/
b;
406 <<
"Selected: t = "<< gt <<
", from ["<< tmin <<
", "<< tmax <<
"]";
412 double totxsec = evrec->
XSec();
413 double wght = (vol/totxsec)*xsec;
414 LOG(
"COHKinematics",
pNOTICE) <<
"Kinematics wght = "<< wght;
418 LOG(
"COHKinematics",
pNOTICE) <<
"Current event wght = " << wght;
445 LOG(
"COHKinematics",
pNOTICE) <<
"Using AlvarezRuso Model";
464 const double E_l_min = interaction->
FSPrimLepton()->Mass();
467 const double ctheta_l_min = 0.4;
470 const double ctheta_pi_min = 0.4;
471 const double ctheta_pi_max = 1.0 -
kASmallNum;
476 const double d_E_l = E_l_max - E_l_min;
477 const double d_ctheta_l = ctheta_l_max - ctheta_l_min;
478 const double d_ctheta_pi = ctheta_pi_max - ctheta_pi_min;
479 const double d_phi = phi_max - phi_min;
482 unsigned int iter = 0;
484 double xsec=-1, g_E_l=-1, g_theta_l=-1, g_phi_l=-1, g_theta_pi=-1, g_phi_pi=-1;
485 double g_ctheta_l, g_ctheta_pi;
492 g_E_l = E_l_min + d_E_l * rnd->
RndKine().Rndm();
493 g_ctheta_l = ctheta_l_min + d_ctheta_l * rnd->
RndKine().Rndm();
494 g_ctheta_pi = ctheta_pi_min + d_ctheta_pi * rnd->
RndKine().Rndm();
495 g_phi_l = phi_min + d_phi * rnd->
RndKine().Rndm();
497 g_phi_pi = g_phi_l + (phi_min + d_phi * rnd->
RndKine().Rndm());
498 g_theta_l = TMath::ACos(g_ctheta_l);
499 g_theta_pi = TMath::ACos(g_ctheta_pi);
501 LOG(
"COHKinematics",
pINFO) <<
"Trying: Lep(" <<g_E_l <<
", " <<
502 g_theta_l <<
", " << g_phi_l <<
") Pi(" <<
503 g_theta_pi <<
", " << g_phi_pi <<
")";
505 this->
SetKinematics(g_E_l, g_theta_l, g_phi_l, g_theta_pi, g_phi_pi,
506 interaction, interaction->
KinePtr());
513 double t = xsec_max * rnd->
RndKine().Rndm();
515 LOG(
"COHKinematics",
pINFO) <<
"Got: xsec = " << xsec <<
", t = " <<
516 t <<
" (max_xsec = " << xsec_max <<
")";
519 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 521 <<
"xsec= " << xsec <<
", J= 1, Rnd= " <<
t;
531 LOG(
"COHKinematics",
pNOTICE) <<
"Selected: Lepton(" <<
532 g_E_l <<
", " << g_theta_l <<
", " <<
533 g_phi_l <<
") Pion(" << g_theta_pi <<
", " << g_phi_pi <<
")";
536 double theta_l = g_theta_l;
537 double theta_pi = g_theta_pi;
538 double phi_l = g_phi_l;
539 double phi_pi = g_phi_pi;
541 double E_nu = P4_nu.E();
542 double E_pi= E_nu-E_l;
546 double p_l = TMath::Sqrt(E_l*E_l - m_l*m_l);
547 TVector3 lepton_3vector = TVector3(0,0,0);
548 lepton_3vector.SetMagThetaPhi(p_l,theta_l,phi_l);
549 TLorentzVector P4_lep = TLorentzVector(lepton_3vector , E_l );
551 double p_pi = TMath::Sqrt(E_pi*E_pi - m_pi*m_pi);
552 TVector3 pion_3vector = TVector3(0,0,0);
553 pion_3vector.SetMagThetaPhi(p_pi,theta_pi,phi_pi);
554 TLorentzVector P4_pion = TLorentzVector(pion_3vector , E_pi);
556 TLorentzVector q = P4_nu - P4_lep;
557 double Q2 = -q.Mag2();
559 double y = E_pi/E_nu;
561 double t = TMath::Abs( (q - P4_pion).Mag2() );
567 double vol = d_E_l*d_ctheta_l*d_phi*d_ctheta_pi*d_phi;
568 double totxsec = evrec->
XSec();
569 double wght = (vol/totxsec)*xsec;
570 LOG(
"COHKinematics",
pNOTICE) <<
"Kinematics wght = "<< wght;
574 LOG(
"COHKinematics",
pNOTICE) <<
"Current event wght = " << wght;
596 const double theta_l,
598 const double theta_pi,
604 double E_nu = P4_nu.E();
605 double E_pi= E_nu-E_l;
615 p_l = TMath::Sqrt(E_l*E_l - m_l*m_l);
617 TVector3 lepton_3vector = TVector3(0,0,0);
618 lepton_3vector.SetMagThetaPhi(p_l,theta_l,phi_l);
619 TLorentzVector P4_lep = TLorentzVector(lepton_3vector , E_l );
623 p_pi = TMath::Sqrt(E_pi*E_pi - m_pi*m_pi);
625 TVector3 pion_3vector = TVector3(0,0,0);
626 pion_3vector.SetMagThetaPhi(p_pi,theta_pi,phi_pi);
627 TLorentzVector P4_pion = TLorentzVector(pion_3vector , E_pi);
629 double Q2 = -(P4_nu-P4_lep).Mag2();
631 double y = E_pi/E_nu;
649 double E_nu = P4_nu.E();
650 double E_pi= E_nu-E_l;
675 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 677 <<
"Scanning the allowed phase space {K} for the max(dxsec/d{K})";
679 double max_xsec = 0.;
682 }
else if ((
fXSecModel->
Id().
Name() ==
"genie::BergerSehgalCOHPiPXSec2015")) {
684 }
else if ((
fXSecModel->
Id().
Name() ==
"genie::BergerSehgalFMCOHPiPXSec2015")) {
691 "ComputeMaxXSec >> Cannot calculate max cross-section for " <<
699 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 701 SLOG(
"COHKinematics",
pDEBUG) <<
"Max xsec in phase space = " << max_xsec;
719 const double logQ2max = TMath::Log10(Q2r.
max);
720 const double dlogQ2 = (logQ2max - logQ2min) /(NQ2-1);
722 for(
int i=0; i<NQ2; i++) {
723 double Q2 = TMath::Power(10, logQ2min + i * dlogQ2);
728 (yr.
max > 1) || (yr.
min < 0)) {
731 const double logymin = TMath::Log10(yr.
min);
732 const double logymax = TMath::Log10(yr.
max);
733 const double dlogy = (logymax - logymin) /(Ny-1);
735 for(
int j=0; j<Ny; j++) {
736 double gy = TMath::Power(10, logymin + j * dlogy);
744 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 746 <<
"xsec(Q2= " << Q2 <<
", y= " << gy <<
", t = " <<
gt <<
") = " << xsec;
748 max_xsec = TMath::Max(max_xsec, xsec);
768 const double logQ2max = TMath::Log10(Q2r.
max);
769 const double logtmin = TMath::Log10(
kASmallNum);
771 const double dlogQ2 = (logQ2max - logQ2min) /(NQ2-1);
772 const double dlogt = (logtmax - logtmin) /(Nt-1);
774 for(
int i=0; i<NQ2; i++) {
775 double Q2 = TMath::Power(10, logQ2min + i * dlogQ2);
780 (yr.
max > 1) || (yr.
min < 0)) {
783 const double logymin = TMath::Log10(yr.
min);
784 const double logymax = TMath::Log10(yr.
max);
785 const double dlogy = (logymax - logymin) /(Ny-1);
787 for(
int j=0; j<Ny; j++) {
788 double gy = TMath::Power(10, logymin + j * dlogy);
790 for(
int k=0;
k<Nt;
k++) {
791 double gt = TMath::Power(10, logtmin +
k * dlogt);
797 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 799 <<
"xsec(Q2= " << Q2 <<
", y= " << gy <<
", t = " << gt <<
") = " << xsec;
801 max_xsec = TMath::Max(max_xsec, xsec);
820 const double logxmin = TMath::Log10(1
E-5);
821 const double logxmax = TMath::Log10(1.0);
822 const double logymin = TMath::Log10(y.
min);
823 const double logymax = TMath::Log10(y.
max);
825 const double dlogx = (logxmax - logxmin) /(Nx-1);
826 const double dlogy = (logymax - logymin) /(Ny-1);
828 for(
int i=0; i<Nx; i++) {
829 double gx = TMath::Power(10, logxmin + i * dlogx);
830 for(
int j=0; j<Ny; j++) {
831 double gy = TMath::Power(10, logymin + j * dlogy);
834 if(Ev>1.0 && Q2>0.01)
continue;
840 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 842 <<
"xsec(x= " << gx <<
", y= " << gy <<
") = " << xsec;
844 max_xsec = TMath::Max(max_xsec, xsec);
858 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 860 <<
"Scanning the allowed phase space {K} for the max(dxsec/d{K})";
862 double max_xsec = 0.;
868 ROOT::Math::Minimizer *
min = ROOT::Math::Factory::CreateMinimizer(
"Minuit2");
872 min->SetFunction( f );
873 min->SetMaxFunctionCalls(10000);
874 min->SetTolerance(0.05);
878 const unsigned int n_el = 100;
879 const double d_el = (max_el - min_el) /
double(n_el - 1);
882 const double max_thetal =
kPi / 4.0;
883 const unsigned int n_thetal = 10;
884 const double d_thetal = (max_thetal - min_thetal) /
double(n_thetal - 1);
887 const double max_thetapi =
kPi / 2.0;
888 const unsigned int n_thetapi = 10;
889 const double d_thetapi = (max_thetapi - min_thetapi) /
double(n_thetapi - 1);
895 const unsigned int n_phipi = 10;
896 const double d_phipi = (max_phipi - min_phipi) /
double(n_phipi - 1);
898 min->SetLimitedVariable ( 0 ,
"E_lep" , max_el -kASmallNum , d_el , min_el , max_el );
899 min->SetLimitedVariable ( 1 ,
"theta_l" , min_thetal +kASmallNum , d_thetal , min_thetal , max_thetal );
900 min->SetLimitedVariable ( 2 ,
"theta_pi" , min_thetapi+kASmallNum , d_thetapi , min_thetapi, max_thetapi );
901 min->SetLimitedVariable ( 3 ,
"phi_pi" , min_phipi +kASmallNum , d_phipi , min_phipi , max_phipi );
904 max_xsec = -min->MinValue();
910 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 912 SLOG(
"COHKinematics",
pDEBUG) <<
"Max xsec in phase space = " << max_xsec;
946 <<
"*** Could not select valid kinematics after " 947 << iters <<
" iterations";
950 exception.
SetReason(
"Couldn't select kinematics");
996 gROOT->GetListOfFunctions()->Remove(
fEnvelope);
double COHImportanceSamplingEnvelope(double *x, double *par)
virtual double MaxXSec(GHepRecord *evrec) const
const KPhaseSpace & PhaseSpace(void) const
virtual void SetWeight(double wght)
bool IsWeakCC(void) const
bool fGenerateUniformly
uniform over allowed phase space + event weight?
double MaxXSec_ReinSehgal(const Interaction *in) const
THE MAIN GENIE PROJECT NAMESPACE
double fRo
nuclear scale parameter
static const double kNucleonMass
double Q2(const Interaction *const i)
virtual Interaction * Summary(void) const
static RandomGen * Instance()
Access instance.
void SetQ2(double Q2, bool selected=false)
Kinematics * KinePtr(void) const
double fSafetyFactor
maxxsec -> maxxsec * safety_factor
void CalculateKin_BergerSehgalFM(GHepRecord *event_rec) const
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
A simple [min,max] interval for doubles.
virtual void AssertXSecLimits(const Interaction *in, double xsec, double xsec_max) const
static const double kPi0Mass
double fMaxXSecDiffTolerance
max{100*(xsec-maxxsec)/.5*(xsec+maxxsec)} if xsec>maxxsec
Defines the EventGeneratorI interface.
void CalculateKin_AlvarezRuso(GHepRecord *event_rec) const
Generated/set kinematical variables for an event.
void Configure(const Registry &config)
virtual double Weight(void) const
double x(bool selected=false) const
Range1D_t YLim(void) const
y limits
void SetKinematics(const double E_l, const double theta_l, const double phi_l, const double theta_pi, const double phi_pi, const Interaction *interaction, Kinematics *kinematics) const
Range1D_t Q2Lim(void) const
Q2 limits.
Abstract class. Provides a data caching mechanism for for concrete implementations of the EventRecord...
A singleton holding random number generator classes. All random number generation in GENIE should tak...
double MaxXSec_AlvarezRuso(const Interaction *in) const
string AsString(void) const
const XSecAlgorithmI * fXSecModel
void CalculateKin_ReinSehgal(GHepRecord *event_rec) const
void ProcessEventRecord(GHepRecord *event_rec) const
Summary information for an interaction.
double fTMax
upper bound for t = (q - p_pi)^2
An exception thrown by EventRecordVisitorI when the normal processing sequence has to be disrupted (f...
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
void SetFSLeptonP4(const TLorentzVector &p4)
static constexpr double cm2
double fQ2Min
lower bound of integration for Q^2 in Berger-Sehgal Model
virtual void Configure(const Registry &config)
void Sett(double t, bool selected=false)
static const double kASmallNum
void CalculateKin_BergerSehgal(GHepRecord *event_rec) const
static int max(int a, int b)
TRandom3 & RndKine(void) const
rnd number generator used by kinematics generators
Misc GENIE control constants.
const UInt_t kISkipKinematicChk
if set, skip kinematic validity checks
virtual const XSecAlgorithmI * CrossSectionAlg(void) const =0
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
void SetFactor(double factor)
void Setx(double x, bool selected=false)
void UpdateWQ2FromXY(const Interaction *in)
static RunningThreadInfo * Instance(void)
double ComputeMaxXSec(const Interaction *in) const
virtual const AlgId & Id(void) const
Get algorithm ID.
void SetW(double W, bool selected=false)
void SwitchOnFastForward(void)
double Energy(const Interaction *in) const
void SetReason(string reason)
virtual TBits * EventFlags(void) const
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
static const double kPionMass
A registry. Provides the container for algorithm configuration parameters.
void Sety(double y, bool selected=false)
~COHKinematicsGenerator()
double MaxXSec_BergerSehgal(const Interaction *in) const
void UpdateXFromQ2Y(const Interaction *in)
bool CheckKinematics(const double E_l, const double theta_l, const double phi_l, const double theta_pi, const double phi_pi, const Interaction *interaction) const
double fEMin
min E for which maxxsec is cached - forcing explicit calc.
double pionMass(const Interaction *in) const
void SetHadSystP4(const TLorentzVector &p4)
virtual double XSec(void) const
static constexpr double fermi
InitialState * InitStatePtr(void) const
TF2 * fEnvelope
2-D envelope used for importance sampling
const InitialState & InitState(void) const
const ProcessInfo & ProcInfo(void) const
void ClearRunningValues(void)
bool GetParamDef(const RgKey &name, T &p, const T &def) const
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
const Target & Tgt(void) const
static const unsigned int kRjMaxIterations
const EventGeneratorI * RunningThread(void)
#define SLOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a short string (using the FUNCTION and...
double ProbeE(RefFrame_t rf) const
GENIE's GHEP MC event record.
double MaxXSec_BergerSehgalFM(const Interaction *in) const
Keep info on the event generation thread currently on charge. This is used so that event generation m...
Root of GENIE utility namespaces.
cet::coded_exception< error, detail::translate > exception
double fQ2Max
upper bound of integration for Q^2 in Berger-Sehgal Model
TLorentzVector * GetProbeP4(RefFrame_t rf=kRfHitNucRest) const
const UInt_t kISkipProcessChk
if set, skip process validity checks
virtual void SetDiffXSec(double xsec, KinePhaseSpace_t ps)
void throwOnTooManyIterations(unsigned int iters, GHepRecord *evrec) const
Initial State information.
static const double kPionMass2