22 #include "Framework/Conventions/GBuild.h" 32 using namespace genie;
71 double dW = (W.
max-W.
min)/(kNW-1);
73 for(
int iw=0; iw<kNW; iw++) {
76 double dQ2 = (Q2.
max-Q2.
min);
97 const double kdV = kdx*kdy;
99 double cW=-1, cQ2 = -1;
103 for(
int ix=0; ix<kNx; ix++) {
104 double x = kminx+ix*kdx;
105 for(
int iy=0; iy<kNy; iy++) {
106 double y = kminy+iy*kdy;
123 <<
"Couldn't compute phase space volume for " 153 <<
"Computing Jacobian for transformation: " 162 if ( fromps == tops )
191 J = 1. / (kine.
x() * kine.
y());
200 J = TMath::Log(10.)*kine.
x() * TMath::Log(10.)*kine.
Q2();
209 J = 1. / (kine.
Q2() * kine.
y());
252 J = TMath::Power(2*M*Ev,2) *
y;
266 J = 2*TMath::Power(M*Ev,2) * y/
W;
291 double mv = init_state.
Probe()->Mass();
293 double pv = std::sqrt(
std::max(0., Ev*Ev - mv*mv) );
296 const TLorentzVector& hit_nuc_P4 = init_state.
Tgt().
HitNucP4();
297 double M = hit_nuc_P4.M();
304 double El = Ev - ( (W*W + Q2 - M*M) / (2.*M) );
305 double pl = std::sqrt(
std::max(0., El*El - ml*ml) );
309 J = W / ( 2. * pv * pl * M );
313 std::ostringstream
msg;
314 msg <<
"Can not compute Jacobian for transforming: " 317 SLOG(
"KineLimits",
pFATAL) <<
"*** " << msg.str();
323 if(!forward) J = 1./
J;
334 bool matched = ( (a==inpa&&b==inpb) || (a==inpb&&b==inpa) );
337 if(a==inpa&&b==inpb) fwd =
true;
349 double M2 = TMath::Power(M,2);
350 double s = M2 + 2*M*Ev;
355 W.
max = TMath::Sqrt(s) - ml;
367 double Ev,
double M,
double ml,
double W,
double Q2min_cut)
375 double M2 = TMath::Power(M, 2.);
376 double ml2 = TMath::Power(ml, 2.);
377 double W2 = TMath::Power(W, 2.);
378 double s = M2 + 2*M*Ev;
384 double auxC = 0.5*(s-M2)/s;
385 double aux1 = s + ml2 - W2;
386 double aux2 = aux1*aux1 - 4*s*ml2;
388 (aux2 < 0) ? ( aux2 = 0 ) : ( aux2 = TMath::Sqrt(aux2) );
390 Q2.
max = -ml2 + auxC * (aux1 + aux2);
391 Q2.
min = -ml2 + auxC * (aux1 - aux2);
394 Q2.
max = TMath::Max(0., Q2.
max);
395 Q2.
min = TMath::Max(0., Q2.
min);
398 if(Q2.
min < Q2min_cut) {Q2.
min = Q2min_cut; }
405 double Ev,
double M,
double ml,
double W,
double q2min_cut)
417 double Ev,
double M,
double ml,
double Q2min_cut)
426 if(W.
min<0)
return Q2;
433 double Ev,
double M,
double ml,
double q2min_cut)
449 double M2 = TMath::Power(M, 2.);
450 double ml2 = TMath::Power(ml,2.);
451 double s = M2 + 2*M*Ev;
473 assert(xl.
min>0 && xl.
max>0);
475 const unsigned int N=100;
476 const double logxmin = TMath::Log10(xl.
min);
477 const double logxmax = TMath::Log10(xl.
max);
478 const double dlogx = (logxmax-logxmin) / (
double)(N-1);
480 for(
unsigned int i=0; i<
N; i++) {
481 double x = TMath::Power(10, logxmin + i*dlogx);
488 if(y.
max >= 0 && y.
max <= 1 && y.
min >= 0 && y.
min <= 1) {
500 double Ev,
double M,
double ml,
double x)
509 double Ev2 = TMath::Power(Ev,2);
510 double ml2 = TMath::Power(ml,2);
518 double a = 0.5 * ml2/(M*Ev*
x);
520 double c = 1 + 0.5*x*M/Ev;
521 double d = TMath::Max(0., TMath::Power(1-a,2.) - b);
523 double A = 0.5 * (1-a-0.5*
b)/c;
524 double B = 0.5 * TMath::Sqrt(d)/
c;
538 double M2 = TMath::Power(M,2);
539 double ml2 = TMath::Power(ml,2);
540 double s = M2 + 2*M*El + ml2;
545 W.
max = TMath::Sqrt(s) - ml;
557 double El,
double ml,
double M,
double W)
565 double M2 = TMath::Power(M, 2.);
566 double ml2 = TMath::Power(ml, 2.);
567 double W2 = TMath::Power(W, 2.);
568 double s = M2 + 2*M*El + ml2;
574 double auxC = 0.5*(s - M2 - ml2)/s;
575 double aux1 = s + ml2 - W2;
576 double aux2 = aux1*aux1 - 4*s*ml2;
578 (aux2 < 0) ? ( aux2 = 0 ) : ( aux2 = TMath::Sqrt(aux2) );
580 Q2.
max = -ml2 + auxC * (aux1 + aux2);
581 Q2.
min = -ml2 + auxC * (aux1 - aux2);
584 Q2.
max = TMath::Max(0., Q2.
max);
585 Q2.
min = TMath::Max(0., Q2.
min);
595 double El,
double ml,
double M,
double W)
607 double El,
double ml,
double M)
616 if(W.
min<0)
return Q2;
623 double El,
double ml,
double M)
639 double M2 = TMath::Power(M, 2.);
640 double ml2 = TMath::Power(ml,2.);
641 double s = M2 + 2*M*El + ml2;
645 assert (s > M2 + ml2);
663 assert(xl.
min>0 && xl.
max>0);
665 const unsigned int N=100;
666 const double logxmin = TMath::Log10(xl.
min);
667 const double logxmax = TMath::Log10(xl.
max);
668 const double dlogx = (logxmax-logxmin) / (
double)(N-1);
670 for(
unsigned int i=0; i<
N; i++) {
671 double x = TMath::Power(10, logxmin + i*dlogx);
678 if(y.
max >= 0 && y.
max <= 1 && y.
min >= 0 && y.
min <= 1) {
690 double El,
double ml,
double M,
double x)
699 double El2 = TMath::Power(El,2);
700 double ml2 = TMath::Power(ml,2);
708 double a = 0.5 * ml2/(M*El*
x);
710 double c = 1 + 0.5*x*M/El;
711 double d = TMath::Max(0., TMath::Power(1-a,2.) - b);
713 double A = 0.5 * (1-a-0.5*
b)/c;
714 double B = 0.5 * TMath::Sqrt(d)/
c;
743 double Mn2 = Mn * Mn;
744 double mlep2 = mlep * mlep;
745 double s = Mn2 + 2.0 * Mn * Ev;
746 double W2min =
CohW2Min(Mn, m_produced);
750 double b = mlep2 /
s;
751 double c = W2min /
s;
752 double lambda = a * a + b * b + c * c - 2.0 * a * b - 2.0 * a * c - 2.0 * b *
c;
754 double A = (s - Mn * Mn) / 2.0;
755 double B = 1 - TMath::Sqrt(lambda);
756 double C = 0.5 * (W2min + mlep2 - Mn2 * (W2min - mlep2) / s );
759 <<
"Q2 kinematic limits calculation failed for CohQ2Lim. " 760 <<
"Assuming Q2min = 0.0";
762 Q2.
min = TMath::Max(0., A * B - C);
765 <<
"Q2 kinematic limits calculation failed for CohQ2Lim. " 766 <<
"Assuming Q2min = 0.0";
782 double Ev,
double Q2)
795 double s = Mn * Mn + 2.0 * Mn * Ev;
796 double Mnterm = 1 - Mn * Mn /
s;
797 double Mlterm = 1 - mlep * mlep /
s;
799 double T1 = 0.25 * s * s * Mnterm * Mnterm * Mlterm;
800 double T2 = Q2 - (0.5 * s * Mnterm) + (0.5 * mlep * mlep * Mnterm);
803 W2l.
max = (T1 - T2 * T2 ) *
805 (1.0 / (Q2 + mlep * mlep));
811 double Q2,
double Mn,
double xsi)
817 double nu_min = (W2min + Q2 - Mn * Mn) / (2.0 * Mn);
818 double nu_max = (W2max + Q2 - Mn * Mn) / (2.0 * Mn);
819 double xsiQ = xsi * TMath::Sqrt(Q2);
821 nul.
min = (xsiQ > nu_min) ? xsiQ : nu_min;
828 double Ev,
double Q2,
double xsi)
835 if (W2lim.
min > W2lim.
max) {
837 <<
"Kinematically forbidden region in CohYLim. W2min = " << W2lim.
min 838 <<
"; W2max =" << W2lim.
max;
840 <<
" Mn = " << Mn <<
"; m_had_sys = " << m_produced <<
"; mlep = " 841 << mlep <<
"; Ev = " << Ev <<
"; Q2 = " <<
Q2;
846 ylim.
min = nulim.
min / Ev;
847 ylim.
max = nulim.
max / Ev;
870 return (Mn + m_produced) * (Mn + m_produced);
883 double M2 = TMath::Power(M,2);
884 double ml2 = TMath::Power(ml,2);
885 double s = M2 + 2*M*Ev + ml2;
891 W.
max = TMath::Sqrt(s) - ml;
903 double Ev,
double M,
double ml,
double W,
double Q2min_cut)
911 double M2 = TMath::Power(M, 2.);
912 double ml2 = TMath::Power(ml, 2.);
913 double W2 = TMath::Power(W, 2.);
914 double s = M2 + 2*M*Ev + ml2;
916 double sqs = TMath::Sqrt(s);
917 double E1CM = (s + ml2 - M2) / (2.*sqs);
918 double p1CM = TMath::Max(0., E1CM*E1CM - ml2);
919 p1CM = TMath::Sqrt(p1CM);
920 double E3CM = (s + ml2 - W2) / (2.*sqs);
921 double p3CM = TMath::Max(0., E3CM*E3CM - ml2);
922 p3CM = TMath::Sqrt(p3CM);
928 SLOG(
"KineLimits",
pDEBUG) <<
"E1_CM = " << E1CM;
929 SLOG(
"KineLimits",
pDEBUG) <<
"p1_CM = " << p1CM;
930 SLOG(
"KineLimits",
pDEBUG) <<
"E3_CM = " << E3CM;
931 SLOG(
"KineLimits",
pDEBUG) <<
"p3_CM = " << p3CM;
933 Q2.
min = TMath::Power(p3CM - p1CM,2) - TMath::Power((W2 - M2) / (2.*sqs),2);
934 Q2.
max = TMath::Power(p3CM + p1CM,2) - TMath::Power((W2 - M2) / (2.*sqs),2);
936 SLOG(
"KineLimits",
pDEBUG) <<
"Nominal Q^2 limits: " << Q2.
min <<
" , " << Q2.
max;
938 Q2.
max = TMath::Max(0., Q2.
max);
939 Q2.
min = TMath::Max(0., Q2.
min);
942 if(Q2.
min < Q2min_cut) {Q2.
min = Q2min_cut; }
949 double Ev,
double M,
double ml,
double W,
double q2min_cut)
961 double Ev,
double M,
double ml,
double Q2min_cut)
970 if(W.
min<0)
return Q2;
977 double Ev,
double M,
double ml,
double q2min_cut)
994 double Wmin = Wl.
min;
995 double W2min = Wmin*Wmin;
996 SLOG(
"KineLimits",
pDEBUG) <<
"W^2_min = " << W2min;
998 SLOG(
"KineLimits",
pDEBUG) <<
"Q^2 range : " << Q2l.
min <<
" , " << Q2l.
max;
1001 x.
min = Q2l.
min / (Q2l.
min + W2min - M2);
1002 x.
max = Q2l.
max / (Q2l.
max + W2min - M2);
1012 double Wmin = Wl.
min;
1013 double W2min = Wmin*Wmin;
1017 y.
min = (Q2l.
min + W2min - M2) / (2*Ev*M);
1018 y.
max = (Q2l.
max + W2min - M2) / (2*Ev*M);
1025 double Ev,
double M,
double ml,
double x)
1036 double Wmin = Wl.
min;
1037 double W2min = Wmin*Wmin;
1040 y.
min = (W2min - M2) / (1.-x) / (2*Ev*M);
1041 y.
max = 2.* M * x *(Ev*Ev - ml2) / Ev / (2. * M * Ev * x + M2 * x * x + ml2);
1072 double Q2 = kinematics.
Q2();
1078 double x = kinematics.
x();
1079 double y = kinematics.
y();
1081 double Q2 = 2*Mn*Ev*x*
y;
1101 double W = kinematics.
W();
1109 double x = kinematics.
x();
1110 double y = kinematics.
y();
1111 double W2 = TMath::Max(0., M2 + 2*Ev*M*y*(1-x));
1112 double W = TMath::Sqrt(W2);
1120 double Ev,
double M,
double W,
double Q2,
double &
x,
double &
y)
1127 double M2 = TMath::Power(M,2);
1128 double W2 = TMath::Power(W,2);
1130 x = Q2 / (W2-M2+
Q2);
1131 y = (W2-M2+
Q2) / (2*M*Ev);
1133 x = TMath::Min(1.,x);
1134 y = TMath::Min(1.,y);
1135 x = TMath::Max(0.,x);
1136 y = TMath::Max(0.,y);
1139 <<
"(W=" << W <<
",Q2=" << Q2 <<
") => (x="<< x <<
", y=" << y<<
")";
1143 double Ev,
double M,
double &
W,
double &
Q2,
double x,
double y)
1150 double M2 = TMath::Power(M,2);
1151 double W2 = M2 + 2*Ev*M*y*(1-
x);
1153 W = TMath::Sqrt(TMath::Max(0., W2));
1157 <<
"(x=" << x <<
",y=" << y <<
" => (W=" << W <<
",Q2=" << Q2 <<
")";
1161 double Ev,
double M,
double &
W,
double Q2,
double x,
double &
y)
1168 double M2 = TMath::Power(M,2);
1169 y = Q2 / (2 * x * M * Ev);
1170 y = TMath::Min(1.,y);
1172 double W2 = M2 + 2*Ev*M*y*(1-
x);
1173 W = TMath::Sqrt(TMath::Max(0., W2));
1176 <<
"(x=" << x <<
",Q2=" << Q2 <<
" => (W=" << W <<
",Y=" << y <<
")";
1180 double Ev,
double M,
double x,
double y)
1186 double M2 = TMath::Power(M,2);
1187 double W2 = M2 + 2*Ev*M*y*(1-
x);
1188 double W = TMath::Sqrt(TMath::Max(0., W2));
1190 LOG(
"KineLimits",
pDEBUG) <<
"(x=" << x <<
",y=" << y <<
") => W=" <<
W;
1196 double Ev,
double M,
double x,
double y)
1202 double Q2 = 2*x*y*M*Ev;
1204 LOG(
"KineLimits",
pDEBUG) <<
"(x=" << x <<
",y=" << y <<
") => Q2=" <<
Q2;
1210 double Ev,
double M,
double Q2,
double y)
1215 assert(Ev > 0. && M > 0. && Q2 > 0. && y > 0.);
1217 double x = Q2 / (2. * y * M * Ev);
1219 LOG(
"KineLimits",
pDEBUG) <<
"(Ev=" << Ev <<
",Q2=" << Q2
1220 <<
",y=" << y <<
",M=" << M <<
") => x=" <<
x;
1226 double x,
double Q2,
double M,
double mc)
1233 double M2 = TMath::Power(M,2);
1234 double v = 0.5*Q2/(M*
x);
1235 double W2 = TMath::Max(0., M2+2*M*v-Q2);
1236 double W = TMath::Sqrt(W2);
1240 if(xc>=1 || W<=Wmin)
return false;
1245 double x,
double Q2,
double M,
double mc)
1252 double mc2 = TMath::Power(mc,2);
1253 double v = 0.5*Q2/(M*
x);
1254 double xc = x + 0.5*mc2/(M*v);
1259 Range1D_t & range,
double min_cut,
double max_cut)
1270 if (min_cut > range.
max || max_cut < range.
min) {
1285 double x = kine->
x();
1286 double y = kine->
y();
1303 double W = kine->
W();
1304 double Q2 = kine->
Q2();
1321 double x = kine->
x();
1322 double Q2 = kine->
Q2();
1351 double y = kine->
y();
1352 double Q2 = kine->
Q2();
1360 double *
x,
double * par)
1369 double xsmax = par[2];
1370 double Wmax = par[3];
1385 double hwfe = mD+2*gD;
1386 double lwfe = mD-gD/2;
1390 func = xsmax / (1 + 5* TMath::Power((W-lwfe)/gD,2));
1391 }
else if (W > hwfe) {
1393 func = xsmax / (1 + TMath::Power((W-hwfe)/gD,2));
1403 double plateau_edge = Wmax-0.1;
1405 if (W > plateau_edge) {
1410 func = xsmax / (1 + TMath::Power((W-plateau_edge)/gD,2));
1415 func *= (1 - (0.5-QD2));
1422 double *
x,
double * par)
1429 double xpeak = par[0];
1430 double xmin = par[1];
1432 double xsmax = par[2];
1436 if(xb < xpeak/20.) {
1438 double plateau_edge = xpeak/20.;
1439 double slope = xsmax/(xmin - plateau_edge);
1440 func = xsmax - slope * (xb - plateau_edge);
1441 }
else if (xb > 2*xpeak) {
1443 double plateau_edge = 2*xpeak;
1444 double slope = xsmax/(xmax - plateau_edge);
1445 func = xsmax - slope * (xb - plateau_edge);
1454 double *
x,
double * par)
1461 double xsmax = 3*par[0];
1464 if(yb<0.|| yb>1.)
return 0.;
1465 if(xb<0.|| xb>1.)
return 0.;
1467 if(Ev<1)
return xsmax;
1468 if(xb/Ev<1E-4 && yb>0.95)
return 5*xsmax;
1472 double yp = (Ev>2.5) ? 2.5/Ev : 1;
1479 xs0 = xsmax - (yb-yp)*xsmax;
1481 double d = TMath::Power( (xb-xp)/0.075, 2);
1485 func = xsmax - (yb-yp)*xsmax;
double COHImportanceSamplingEnvelope(double *x, double *par)
static const double kMaxX
double W(bool selected=false) const
const KPhaseSpace & PhaseSpace(void) const
bool TransformMatched(KinePhaseSpace_t ia, KinePhaseSpace_t ib, KinePhaseSpace_t a, KinePhaseSpace_t b, bool &fwd)
double J(double q0, double q3, double Enu, double ml)
THE MAIN GENIE PROJECT NAMESPACE
Range1D_t InelXLim(double El, double ml, double M)
Range1D_t InelWLim(double El, double ml, double M)
void XQ2toWY(double Ev, double M, double &W, double Q2, double x, double &y)
void msg(const char *fmt,...)
void UpdateXYFromWQ2(const Interaction *in)
double Q2(const Interaction *const i)
Range1D_t Inelq2Lim_W(double El, double ml, double M, double W)
double DISImportanceSamplingEnvelope(double *x, double *par)
void SetQ2(double Q2, bool selected=false)
Kinematics * KinePtr(void) const
Range1D_t CohW2Lim(double Mn, double m_produced, double mlep, double Ev, double Q2)
Range1D_t InelWLim(double Ev, double M, double ml)
int RecoilNucleonPdg(void) const
recoil nucleon pdg
A simple [min,max] interval for doubles.
Range1D_t Darkq2Lim(double Ev, double M, double ml, double q2min_cut=-1 *controls::kMinQ2Limit)
bool IsQuasiElastic(void) const
static const double kMaxY
Range1D_t DarkQ2Lim_W(double Ev, double M, double ml, double W, double Q2min_cut=controls::kMinQ2Limit)
bool IsWithinLimits(double x, Range1D_t range)
Generated/set kinematical variables for an event.
static const double kPhotontest
TParticlePDG * Probe(void) const
double x(bool selected=false) const
static const double kLightestChmHad
static const double kMinQ2Limit
static const double kMinY
Range1D_t CEvNSQ2Lim(double Ev)
bool IsCoherentProduction(void) const
Range1D_t Q2Lim_W(void) const
Q2 limits @ fixed W.
Range1D_t InelQ2Lim(double Ev, double M, double ml, double Q2min_cut=controls::kMinQ2Limit)
enum genie::EKinePhaseSpace KinePhaseSpace_t
Range1D_t Limits(KineVar_t kvar) const
Return the kinematical variable limits.
void UpdateWYFromXQ2(const Interaction *in)
double XYtoW(double Ev, double M, double x, double y)
double PhaseSpaceVolume(const Interaction *const i, KinePhaseSpace_t ps)
Range1D_t InelQ2Lim_W(double Ev, double M, double ml, double W, double Q2min_cut=controls::kMinQ2Limit)
double y(bool selected=false) const
double W(const Interaction *const i)
Summary information for an interaction.
const TLorentzVector & HitNucP4(void) const
double SlowRescalingVar(double x, double Q2, double M, double mc)
Range1D_t InelXLim(double Ev, double M, double ml)
Range1D_t DarkQ2Lim(double Ev, double M, double ml, double Q2min_cut=controls::kMinQ2Limit)
Range1D_t Inelq2Lim_W(double Ev, double M, double ml, double W, double q2min_cut=-1 *controls::kMinQ2Limit)
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Range1D_t CohNuLim(double W2min, double W2max, double Q2, double Mn, double xsi)
bool KVSet(KineVar_t kv) const
double QD2toQ2(double QD2)
Exception used inside Interaction classes.
Range1D_t Darkq2Lim_W(double Ev, double M, double ml, double W, double q2min_cut=-1 *controls::kMinQ2Limit)
static string AsString(KinePhaseSpace_t kps)
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
double XYtoQ2(double Ev, double M, double x, double y)
void XYtoWQ2(double Ev, double M, double &W, double &Q2, double x, double y)
void WQ2toXY(double Ev, double M, double W, double Q2, double &x, double &y)
Range1D_t Cohq2Lim(double Mn, double m_produced, double mlep, double Ev)
const Kinematics & Kine(void) const
static const double kNeutronMass
static const double kMinX
static const double kASmallNum
double GetKV(KineVar_t kv) const
bool IsAboveCharmThreshold(double x, double Q2, double M, double mc)
static int max(int a, int b)
static constexpr double ps
double Q2toQD2(double Q2)
Range1D_t InelQ2Lim(double El, double ml, double M)
Range1D_t Inelq2Lim(double El, double ml, double M)
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
void Setx(double x, bool selected=false)
static const double kMQD2
void UpdateWQ2FromXY(const Interaction *in)
Range1D_t InelYLim(double El, double ml, double M)
void SetW(double W, bool selected=false)
TLorentzVector * HitNucP4Ptr(void) const
Range1D_t InelQ2Lim_W(double El, double ml, double M, double W)
static PDGLibrary * Instance(void)
double CohW2Min(double Mn, double m_produced)
static const double kPionMass
Range1D_t InelYLim_X(double Ev, double M, double ml, double x)
double Q2YtoX(double Ev, double M, double Q2, double y)
Range1D_t CohYLim(double Mn, double m_produced, double mlep, double Ev, double Q2, double xsi)
Range1D_t DarkYLim(double Ev, double M, double ml)
void Sety(double y, bool selected=false)
void UpdateXFromQ2Y(const Interaction *in)
static const double kAVerySmallNum
Range1D_t InelYLim_X(double El, double ml, double M, double x)
Range1D_t CohQ2Lim(double Mn, double m_produced, double mlep, double Ev)
double Jacobian(const Interaction *const i, KinePhaseSpace_t f, KinePhaseSpace_t t)
const InitialState & InitState(void) const
const ProcessInfo & ProcInfo(void) const
Range1D_t DarkYLim_X(double Ev, double M, double ml, double x)
Range1D_t Inelq2Lim(double Ev, double M, double ml, double q2min_cut=-1 *controls::kMinQ2Limit)
TParticlePDG * Find(int pdgc, bool must_exist=true)
double Q2(bool selected=false) const
const Target & Tgt(void) const
#define SLOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a short string (using the FUNCTION and...
Range1D_t InelYLim(double Ev, double M, double ml)
double ProbeE(RefFrame_t rf) const
Range1D_t DarkXLim(double Ev, double M, double ml)
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils...
void ApplyCutsToKineLimits(Range1D_t &r, double min, double max)
Range1D_t DarkWLim(double Ev, double M, double ml)
Initial State information.
double RESImportanceSamplingEnvelope(double *x, double *par)