GENIEGenerator
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
KineUtils.cxx
Go to the documentation of this file.
1 //____________________________________________________________________________
2 /*
3  Copyright (c) 2003-2024, The GENIE Collaboration
4  For the full text of the license visit http://copyright.genie-mc.org
5 
6  Costas Andreopoulos <c.andreopoulos \at cern.ch>
7  University of Liverpool
8 
9  Afroditi Papadopoulou <apapadop \at mit.edu>
10  Massachusetts Institute of Technology
11 
12  Changes required to implement the GENIE Boosted Dark Matter module
13  were installed by Josh Berger (Univ. of Wisconsin)
14 */
15 //____________________________________________________________________________
16 
17 #include <cstdlib>
18 #include <limits>
19 #include <TMath.h>
20 
22 #include "Framework/Conventions/GBuild.h"
29 
30 #include <sstream>
31 
32 using namespace genie;
33 using namespace genie::constants;
34 
35 //____________________________________________________________________________
37  const Interaction * const in, KinePhaseSpace_t ps)
38 {
39  double vol = 0;
40 
41  const KPhaseSpace & phase_space = in->PhaseSpace();
42 
43  switch(ps) {
44 
45  case(kPSQ2fE):
46  {
47  Range1D_t Q2 = phase_space.Limits(kKVQ2);
48  vol = Q2.max - Q2.min;
49  return vol;
50  break;
51  }
52  case(kPSq2fE):
53  {
54  Range1D_t q2 = phase_space.Limits(kKVq2);
55  vol = q2.min - q2.max;
56  return vol;
57  break;
58  }
59  case(kPSWfE):
60  {
61  Range1D_t W = phase_space.Limits(kKVW);
62  vol = W.max - W.min;
63  return vol;
64  break;
65  }
66  case(kPSWQ2fE):
67  {
68  Range1D_t W = phase_space.Limits(kKVW);
69  if(W.max<0) return 0;
70  const int kNW = 100;
71  double dW = (W.max-W.min)/(kNW-1);
72  Interaction interaction(*in);
73  for(int iw=0; iw<kNW; iw++) {
74  interaction.KinePtr()->SetW(W.min + iw*dW);
75  Range1D_t Q2 = interaction.PhaseSpace().Q2Lim_W();
76  double dQ2 = (Q2.max-Q2.min);
77  vol += (dW*dQ2);
78  }
79  return vol;
80  break;
81  }
82  case(kPSxyfE):
83  {
84  Range1D_t W = phase_space.Limits(kKVW);
85  if(W.max<0) return 0;
86 
87  const InitialState & init_state = in->InitState();
88  double Ev = init_state.ProbeE(kRfHitNucRest);
89  double M = init_state.Tgt().HitNucP4Ptr()->M();
90 
91  const int kNx = 100;
92  const int kNy = 100;
93  const double kminx = controls::kMinX;
94  const double kminy = controls::kMinY;
95  const double kdx = (controls::kMaxX - kminx) / (kNx-1);
96  const double kdy = (controls::kMaxY - kminy) / (kNy-1);
97  const double kdV = kdx*kdy;
98 
99  double cW=-1, cQ2 = -1;
100 
101  Interaction interaction(*in);
102 
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;
107 
108  XYtoWQ2(Ev, M, cW, cQ2, x, y);
109  if(!math::IsWithinLimits(cW, W)) continue;
110 
111  interaction.KinePtr()->SetW(cW);
112  Range1D_t Q2 = interaction.PhaseSpace().Q2Lim_W();
113  if(!math::IsWithinLimits(cQ2, Q2)) continue;
114 
115  vol += kdV;
116  }
117  }
118  return vol;
119  break;
120  }
121  default:
122  SLOG("KineLimits", pERROR)
123  << "Couldn't compute phase space volume for "
124  << KinePhaseSpace::AsString(ps) << "using \n" << *in;
125  return 0;
126  }
127  return 0;
128 }
129 //____________________________________________________________________________
131  const Interaction * const i, KinePhaseSpace_t fromps, KinePhaseSpace_t tops)
132 {
133 // Returns the Jacobian for a kinematical transformation:
134 // from_ps (x,y,z,...) -> to_ps (u,v,w,...).
135 //
136 // Note on the convention used here:
137 // The differential cross-section d^n sigma / dxdydz..., is transformed to
138 // d^n sigma / dudvdw... using the Jacobian, computed in this function, as:
139 // as d^n sigma / dudvdw... = J * d^n sigma / dxdydz...
140 // where:
141 //
142 // dxdxdz... = J * dudvdw...
143 //
144 // | dx/du dx/dv dx/dw ... |
145 // | |
146 // d {x,y,z,...} | dy/du dy/dv dy/dw ... |
147 // J = ------------ = | |
148 // d {u,v,w,...} | dz/du dz/dv dz/dw ... |
149 // | |
150 // | ... ... ... ... |
151 //
152  SLOG("KineLimits", pDEBUG)
153  << "Computing Jacobian for transformation: "
154  << KinePhaseSpace::AsString(fromps) << " --> "
155  << KinePhaseSpace::AsString(tops);
156 
157  double J=0;
158  bool forward;
159  const Kinematics & kine = i->Kine();
160 
161  // cover the simple case
162  if ( fromps == tops )
163  {
164  forward = true;
165  J = 1.;
166  }
167  //
168  // transformation: {Q2}|E -> {lnQ2}|E
169  //
170  else
171  if ( TransformMatched(fromps,tops,kPSQ2fE,kPSlogQ2fE,forward) )
172  {
173  J = 1. / kine.Q2();
174  }
175 
176  //
177  // transformation: {QD2}|E -> {Q2}|E
178  //
179  else
180  if ( TransformMatched(fromps,tops,kPSQD2fE,kPSQ2fE,forward) )
181  {
182  J = TMath::Power(1+kine.Q2()/controls::kMQD2,-2)/controls::kMQD2;
183  }
184 
185  //
186  // transformation: {x,y}|E -> {lnx,lny}|E
187  //
188  else
189  if ( TransformMatched(fromps,tops,kPSxyfE,kPSlogxlogyfE,forward) )
190  {
191  J = 1. / (kine.x() * kine.y());
192  }
193 
194  //
195  // transformation: {log10x,log10Q2}|E -> {x,Q2}|E
196  //
197  else
198  if ( TransformMatched(fromps,tops,kPSxQ2fE,kPSlog10xlog10Q2fE,forward) )
199  {
200  J = TMath::Log(10.)*kine.x() * TMath::Log(10.)*kine.Q2();
201  }
202 
203  //
204  // transformation: {Q2,y}|E -> {lnQ2,lny}|E
205  //
206  else
207  if ( TransformMatched(fromps,tops,kPSQ2yfE,kPSlogQ2logyfE,forward) )
208  {
209  J = 1. / (kine.Q2() * kine.y());
210  }
211 
212  //
213  // transformation: {x,y}|E -> {x,Q2}|E
214  //
215  else
216  if ( TransformMatched(fromps,tops,kPSxQ2fE,kPSxyfE,forward) )
217  {
218  const InitialState & init_state = i->InitState();
219  double Ev = init_state.ProbeE(kRfHitNucRest);
220  double M = init_state.Tgt().HitNucP4Ptr()->M();
221  double x = kine.x();
222  J = 2*x*Ev*M;
223  }
224 
225  //
226  // transformation: {Q2,y}|E -> {x,y}|E
227  //
228  else
229  if ( TransformMatched(fromps,tops,kPSQ2yfE,kPSxyfE,forward) )
230  {
231  const InitialState & init_state = i->InitState();
232  double Ev = init_state.ProbeE(kRfHitNucRest);
233  double M = init_state.Tgt().HitNucP4Ptr()->M();
234  double y = kine.y();
235  J = 2*y*Ev*M;
236  }
237 
238  //
239  // transformation: {W,Q2}|E -> {W,lnQ2}|E
240  //
241  else if ( TransformMatched(fromps,tops,kPSWQ2fE,kPSWlogQ2fE,forward) )
242  {
243  J = 1. / kine.Q2();
244  }
245 
246  //
247  // transformation: {W,QD2}|E -> {W,Q2}|E
248  //
249  else
250  if ( TransformMatched(fromps,tops,kPSWQD2fE,kPSWQ2fE,forward) )
251  {
252  J = TMath::Power(1+kine.Q2()/controls::kMQD2,-2)/controls::kMQD2;
253  }
254 
255  //
256  // transformation: {W2,Q2}|E --> {x,y}|E
257  //
258  else
259  if ( TransformMatched(fromps,tops,kPSW2Q2fE,kPSxyfE,forward) )
260  {
261  const InitialState & init_state = i->InitState();
262  double Ev = init_state.ProbeE(kRfHitNucRest);
263  double M = init_state.Tgt().HitNucP4Ptr()->M();
264  double y = kine.y();
265  J = TMath::Power(2*M*Ev,2) * y;
266  }
267 
268  //
269  // transformation: {W,Q2}|E -> {x,y}|E
270  //
271  else
272  if ( TransformMatched(fromps,tops,kPSWQ2fE,kPSxyfE,forward) )
273  {
274  const InitialState & init_state = i->InitState();
275  double Ev = init_state.ProbeE(kRfHitNucRest);
276  double M = init_state.Tgt().HitNucP4Ptr()->M();
277  double y = kine.y();
278  double W = kine.W();
279  J = 2*TMath::Power(M*Ev,2) * y/W;
280  }
281 
282  // Transformation: {Omegalep,Omegapi}|E -> {Omegalep,Thetapi}|E
283  else if ( TransformMatched(fromps,tops,kPSElOlOpifE,kPSElOlTpifE,forward) ) {
284  // Use symmetry to turn 4d integral into 3d * 2pi
285  J = 2*constants::kPi;
286  }
287 
288  // Transformation: {Tl,ctl} -> {W,Q2}|E
289  // IMPORTANT NOTE: This transformation can't be done exactly in two
290  // dimensions due to Fermi motion. For fixed {W,Q2}, the azimuthal angle phi
291  // is uniformly distributed in the hit nucleon rest frame, while for fixed
292  // {Tl,ctl} it is uniform in the lab frame (i.e., the target nucleus rest
293  // frame). A 3D Jacobian is needed in order to account for the relationship
294  // between these two azimuthal angles. In the implementation below, I choose
295  // to neglect Fermi motion. Under this approximation, the lab and hit nucleon
296  // rest frames are the same. Doing things this way is a stopgap solution for
297  // the new XSecShape_CCMEC reweighting dial developed for MicroBooNE. A full
298  // solution should use the Jacobian that connects the 3D phase spaces which
299  // each include phi. - S. Gardiner, 29 July 2020
300  else if ( TransformMatched(fromps, tops, kPSTlctl, kPSWQ2fE, forward) )
301  {
302  // Probe properties (mass, energy, momentum)
303  const InitialState& init_state = i->InitState();
304  double mv = init_state.Probe()->Mass();
305  double Ev = init_state.ProbeE( kRfLab );
306  double pv = std::sqrt( std::max(0., Ev*Ev - mv*mv) );
307 
308  // Invariant mass of the initial hit nucleon
309  const TLorentzVector& hit_nuc_P4 = init_state.Tgt().HitNucP4();
310  double M = hit_nuc_P4.M();
311 
312  // Outgoing lepton mass
313  double ml = i->FSPrimLepton()->Mass();
314 
315  double W = i->Kine().GetKV( kKVW );
316  double Q2 = i->Kine().GetKV( kKVQ2 );
317  double El = Ev - ( (W*W + Q2 - M*M) / (2.*M) );
318  double pl = std::sqrt( std::max(0., El*El - ml*ml) );
319 
320  // Compute the Jacobian for {Tl, ctl} --> {W, Q2}
321  // (it will be inverted below for the inverse transformation)
322  J = W / ( 2. * pv * pl * M );
323  }
324 
325  else {
326  std::ostringstream msg;
327  msg << "Can not compute Jacobian for transforming: "
328  << KinePhaseSpace::AsString(fromps) << " --> "
329  << KinePhaseSpace::AsString(tops);
330  SLOG("KineLimits", pFATAL) << "*** " << msg.str();
332  //exit(1);
333  }
334 
335  // if any of the above transforms was reverse, invert the jacobian
336  if(!forward) J = 1./J;
337 
338  return J;
339 }
340 //____________________________________________________________________________
343  KinePhaseSpace_t a, KinePhaseSpace_t b, bool & fwd)
344 {
345 
346  // match a->b or b->a
347  bool matched = ( (a==inpa&&b==inpb) || (a==inpb&&b==inpa) );
348 
349  if(matched) {
350  if(a==inpa&&b==inpb) fwd = true; // matched forward transform: a->b
351  else fwd = false; // matched reverse transform: b->a
352  }
353  return matched;
354 }
355 //____________________________________________________________________________
356 // Kinematical Limits:
357 //____________________________________________________________________________
358 Range1D_t genie::utils::kinematics::InelWLim(double Ev, double M, double ml)
359 {
360 // Computes W limits for inelastic v interactions
361 //
362  double M2 = TMath::Power(M,2);
363  double s = M2 + 2*M*Ev;
364  assert (s>0);
365 
366  Range1D_t W;
368  W.max = TMath::Sqrt(s) - ml;
369  if(W.max<=W.min) {
370  W.min = -1;
371  W.max = -1;
372  return W;
373  }
376  return W;
377 }
378 //____________________________________________________________________________
380  double Ev, double M, double ml, double W, double Q2min_cut)
381 {
382 // Computes Q2 limits (>0) @ the input W for inelastic v interactions
383 
384  Range1D_t Q2;
385  Q2.min = -1;
386  Q2.max = -1;
387 
388  double M2 = TMath::Power(M, 2.);
389  double ml2 = TMath::Power(ml, 2.);
390  double W2 = TMath::Power(W, 2.);
391  double s = M2 + 2*M*Ev;
392 
393  SLOG("KineLimits", pDEBUG) << "s = " << s;
394  SLOG("KineLimits", pDEBUG) << "Ev = " << Ev;
395  assert (s>0);
396 
397  double auxC = 0.5*(s-M2)/s;
398  double aux1 = s + ml2 - W2;
399  double aux2 = aux1*aux1 - 4*s*ml2;
400 
401  (aux2 < 0) ? ( aux2 = 0 ) : ( aux2 = TMath::Sqrt(aux2) );
402 
403  Q2.max = -ml2 + auxC * (aux1 + aux2); // => 0
404  Q2.min = -ml2 + auxC * (aux1 - aux2); // => 0
405 
406  // guard against overflows
407  Q2.max = TMath::Max(0., Q2.max);
408  Q2.min = TMath::Max(0., Q2.min);
409 
410  // limit the minimum Q2
411  if(Q2.min < Q2min_cut) {Q2.min = Q2min_cut; }
412  if(Q2.max < Q2.min ) {Q2.min = -1; Q2.max = -1;}
413 
414  return Q2;
415 }
416 //____________________________________________________________________________
418  double Ev, double M, double ml, double W, double q2min_cut)
419 {
420 // Computes q2 (<0) limits @ the input W for inelastic v interactions
421 
422  Range1D_t Q2 = utils::kinematics::InelQ2Lim_W(Ev,M,ml,W,-1.*q2min_cut);
423  Range1D_t q2;
424  q2.min = - Q2.max;
425  q2.max = - Q2.min;
426  return q2;
427 }
428 //____________________________________________________________________________
430  double Ev, double M, double ml, double Q2min_cut)
431 {
432 // Computes Q2 (>0) limits irrespective of W for inelastic v interactions
433 
434  Range1D_t Q2;
435  Q2.min = -1;
436  Q2.max = -1;
437 
439  if(W.min<0) return Q2;
440 
441  Q2 = utils::kinematics::InelQ2Lim_W(Ev,M,ml,W.min,Q2min_cut);
442  return Q2;
443 }
444 //____________________________________________________________________________
446  double Ev, double M, double ml, double q2min_cut)
447 {
448 // Computes Q2 (>0) limits irrespective of W for inelastic v interactions
449 
450  Range1D_t Q2 = utils::kinematics::InelQ2Lim(Ev,M,ml,-1.*q2min_cut);
451  Range1D_t q2;
452  q2.min = - Q2.max;
453  q2.max = - Q2.min;
454  return q2;
455 }
456 //____________________________________________________________________________
457 Range1D_t genie::utils::kinematics::InelXLim(double Ev, double M, double ml)
458 
459 {
460 // Computes Bjorken x limits for inelastic v interactions
461 
462  double M2 = TMath::Power(M, 2.);
463  double ml2 = TMath::Power(ml,2.);
464  double s = M2 + 2*M*Ev;
465 
466  SLOG("KineLimits", pDEBUG) << "s = " << s;
467  SLOG("KineLimits", pDEBUG) << "Ev = " << Ev;
468  assert (s>M2);
469 
470  Range1D_t x;
471  x.min = ml2/(s-M2) + controls::kAVerySmallNum;
472  x.max = 1. - controls::kAVerySmallNum;
473 
474  return x;
475 }
476 //____________________________________________________________________________
477 Range1D_t genie::utils::kinematics::InelYLim(double Ev, double M, double ml)
478 {
479 // Computes y limits for inelastic v interactions
480 
481  Range1D_t y;
482  y.min = 999;
483  y.max = -999;
484 
485  Range1D_t xl = kinematics::InelXLim(Ev,M,ml);
486  assert(xl.min>0 && xl.max>0);
487 
488  const unsigned int N=100;
489  const double logxmin = TMath::Log10(xl.min);
490  const double logxmax = TMath::Log10(xl.max);
491  const double dlogx = (logxmax-logxmin) / (double)(N-1);
492 
493  for(unsigned int i=0; i<N; i++) {
494  double x = TMath::Power(10, logxmin + i*dlogx);
495 
496  Range1D_t y_x = kinematics::InelYLim_X(Ev,M,ml,x);
497  if(y_x.min>=0 && y_x.min<=1) y.min = TMath::Min(y.min, y_x.min);
498  if(y_x.max>=0 && y_x.max<=1) y.max = TMath::Max(y.max, y_x.max);
499  }
500 
501  if(y.max >= 0 && y.max <= 1 && y.min >= 0 && y.min <= 1) {
502  y.min = TMath::Max(y.min, controls::kAVerySmallNum);
503  y.max = TMath::Min(y.max, 1 - controls::kAVerySmallNum);
504  } else {
505  y.min = -1;
506  y.max = -1;
507  }
508  SLOG("KineLimits", pDEBUG) << "y = [" << y.min << ", " << y.max << "]";
509  return y;
510 }
511 //____________________________________________________________________________
513  double Ev, double M, double ml, double x)
514 
515 {
516 // Computes y limits @ the input x for inelastic v interactions
517 
518  Range1D_t y;
519  y.min = -1;
520  y.max = -1;
521 
522  double Ev2 = TMath::Power(Ev,2);
523  double ml2 = TMath::Power(ml,2);
524 
525  SLOG("KineLimits", pDEBUG) << "x = " << x;
526  SLOG("KineLimits", pDEBUG) << "Ev = " << Ev;
527 
528  assert (Ev>0);
529  assert (x>0&&x<1);
530 
531  double a = 0.5 * ml2/(M*Ev*x);
532  double b = ml2/Ev2;
533  double c = 1 + 0.5*x*M/Ev;
534  double d = TMath::Max(0., TMath::Power(1-a,2.) - b);
535 
536  double A = 0.5 * (1-a-0.5*b)/c;
537  double B = 0.5 * TMath::Sqrt(d)/c;
538 
539  y.min = TMath::Max(0., A-B) + controls::kAVerySmallNum;
540  y.max = TMath::Min(1., A+B) - controls::kAVerySmallNum;
541 
542  return y;
543 }
544 //____________________________________________________________________________
545 // Kinematical Limits for em interactions:
546 //____________________________________________________________________________
548 {
549 // Computes W limits for inelastic em interactions
550 //
551  double M2 = TMath::Power(M,2);
552  double ml2 = TMath::Power(ml,2); // added lepton mass squared to be used in s calculation
553  double s = M2 + 2*M*El + ml2; // non-negligible mass for em interactions
554  assert (s>0);
555 
556  Range1D_t W;
558  W.max = TMath::Sqrt(s) - ml;
559  if(W.max<=W.min) {
560  W.min = -1;
561  W.max = -1;
562  return W;
563  }
566  return W;
567 }
568 //____________________________________________________________________________
570  double El, double ml, double M, double W)
571 {
572 // Computes Q2 limits (>0) @ the input W for inelastic em interactions
573 
574  Range1D_t Q2;
575  Q2.min = -1;
576  Q2.max = -1;
577 
578  double M2 = TMath::Power(M, 2.);
579  double ml2 = TMath::Power(ml, 2.);
580  double W2 = TMath::Power(W, 2.);
581  double s = M2 + 2*M*El + ml2; // added lepton mass squared to be used in s calculation (non-negligible mass for em interactions)
582 
583  SLOG("KineLimits", pDEBUG) << "s = " << s;
584  SLOG("KineLimits", pDEBUG) << "El = " << El;
585  assert (s>0);
586 
587  double auxC = 0.5*(s - M2 - ml2)/s; // subtract ml2 to account for the non-negligible mass of the incoming lepton
588  double aux1 = s + ml2 - W2;
589  double aux2 = aux1*aux1 - 4*s*ml2;
590 
591  (aux2 < 0) ? ( aux2 = 0 ) : ( aux2 = TMath::Sqrt(aux2) );
592 
593  Q2.max = -ml2 + auxC * (aux1 + aux2); // => 0
594  Q2.min = -ml2 + auxC * (aux1 - aux2); // => 0
595 
596  // guard against overflows
597  Q2.max = TMath::Max(0., Q2.max);
598  Q2.min = TMath::Max(0., Q2.min);
599 
600  // limit the minimum Q2
601  if(Q2.min < utils::kinematics::electromagnetic::kMinQ2Limit) {Q2.min = utils::kinematics::electromagnetic::kMinQ2Limit; } // use the relevant threshold for em scattering
602  if(Q2.max < Q2.min ) {Q2.min = -1; Q2.max = -1;}
603 
604  return Q2;
605 }
606 //____________________________________________________________________________
608  double El, double ml, double M, double W)
609 {
610 // Computes q2 (<0) limits @ the input W for inelastic em interactions
611 
613  Range1D_t q2;
614  q2.min = - Q2.max;
615  q2.max = - Q2.min;
616  return q2;
617 }
618 //____________________________________________________________________________
620  double El, double ml, double M)
621 {
622 // Computes Q2 (>0) limits irrespective of W for inelastic em interactions
623 
624  Range1D_t Q2;
625  Q2.min = -1;
626  Q2.max = -1;
627 
629  if(W.min<0) return Q2;
630 
632  return Q2;
633 }
634 //____________________________________________________________________________
636  double El, double ml, double M)
637 {
638 // Computes Q2 (>0) limits irrespective of W for inelastic em interactions
639 
641  Range1D_t q2;
642  q2.min = - Q2.max;
643  q2.max = - Q2.min;
644  return q2;
645 }
646 //____________________________________________________________________________
648 
649 {
650 // Computes Bjorken x limits for inelastic em interactions
651 
652  double M2 = TMath::Power(M, 2.);
653  double ml2 = TMath::Power(ml,2.);
654  double s = M2 + 2*M*El + ml2; // added lepton mass squared to be used in s calculation (non-negligible mass for em interactions)
655 
656  SLOG("KineLimits", pDEBUG) << "s = " << s;
657  SLOG("KineLimits", pDEBUG) << "El = " << El;
658  assert (s > M2 + ml2); // added lepton mass squared to be used in s calculation (non-negligible mass for em interactions)
659 
660  Range1D_t x;
661  x.min = ml2/(s - M2 - ml2) + controls::kASmallNum; // subtracted lepton mass squared to be used in s calculation (non-negligible mass for em interactions)
662  x.max = 1. - controls::kASmallNum;
663 
664  return x;
665 }
666 //____________________________________________________________________________
668 {
669 // Computes y limits for inelastic v interactions
670 
671  Range1D_t y;
672  y.min = 999;
673  y.max = -999;
674 
676  assert(xl.min>0 && xl.max>0);
677 
678  const unsigned int N=100;
679  const double logxmin = TMath::Log10(xl.min);
680  const double logxmax = TMath::Log10(xl.max);
681  const double dlogx = (logxmax-logxmin) / (double)(N-1);
682 
683  for(unsigned int i=0; i<N; i++) {
684  double x = TMath::Power(10, logxmin + i*dlogx);
685 
687  if(y_x.min>=0 && y_x.min<=1) y.min = TMath::Min(y.min, y_x.min);
688  if(y_x.max>=0 && y_x.max<=1) y.max = TMath::Max(y.max, y_x.max);
689  }
690 
691  if(y.max >= 0 && y.max <= 1 && y.min >= 0 && y.min <= 1) {
692  y.min = TMath::Max(y.min, controls::kASmallNum);
693  y.max = TMath::Min(y.max, 1 - controls::kASmallNum);
694  } else {
695  y.min = -1;
696  y.max = -1;
697  }
698  SLOG("KineLimits", pDEBUG) << "y = [" << y.min << ", " << y.max << "]";
699  return y;
700 }
701 //____________________________________________________________________________
703  double El, double ml, double M, double x)
704 
705 {
706 // Computes y limits @ the input x for inelastic em interactions
707 
708  Range1D_t y;
709  y.min = -1;
710  y.max = -1;
711 
712  double El2 = TMath::Power(El,2);
713  double ml2 = TMath::Power(ml,2);
714 
715  SLOG("KineLimits", pDEBUG) << "x = " << x;
716  SLOG("KineLimits", pDEBUG) << "El = " << El;
717 
718  assert (El>0);
719  assert (x>0&&x<1);
720 
721  double a = 0.5 * ml2/(M*El*x);
722  double b = ml2/El2;
723  double c = 1 + 0.5*x*M/El;
724  double d = TMath::Max(0., TMath::Power(1-a,2.) - b);
725 
726  double A = 0.5 * (1-a-0.5*b)/c;
727  double B = 0.5 * TMath::Sqrt(d)/c;
728 
729  y.min = TMath::Max(0., A-B) + controls::kASmallNum;
730  y.max = TMath::Min(1., A+B) - controls::kASmallNum;
731 
732  return y;
733 }
734 //____________________________________________________________________________
736 {
737 // Computes x limits for coherent v interactions
738 
740  return x;
741 }
742 //____________________________________________________________________________
743 Range1D_t genie::utils::kinematics::CohQ2Lim(double Mn, double m_produced, double mlep, double Ev)
744 {
745  // The expressions for Q^2 min appears in PRD 74, 054007 (2006) by
746  // Kartavtsev, Paschos, and Gounaris
747 
748  // That expression is specified for the case of the pion.
749  // GENIE has also Coherent interaction with Single gamma production and Rho
750  // so that formula will be adapted for a generic m_produced
751 
752  Range1D_t Q2;
753  Q2.min = 0.0;
754  Q2.max = std::numeric_limits<double>::max(); // Value must be overriden in user options
755 
756  double Mn2 = Mn * Mn;
757  double mlep2 = mlep * mlep;
758  double s = Mn2 + 2.0 * Mn * Ev;
759  double W2min = CohW2Min(Mn, m_produced);
760 
761  // Looks like Q2min = A * B - C, where A, B, and C are complicated
762  double a = 1.0;
763  double b = mlep2 / s;
764  double c = W2min / s;
765  double lambda = a * a + b * b + c * c - 2.0 * a * b - 2.0 * a * c - 2.0 * b * c;
766  if (lambda > 0) {
767  double A = (s - Mn * Mn) / 2.0;
768  double B = 1 - TMath::Sqrt(lambda);
769  double C = 0.5 * (W2min + mlep2 - Mn2 * (W2min - mlep2) / s );
770  if (A * B - C < 0) {
771  SLOG("KineLimits", pERROR)
772  << "Q2 kinematic limits calculation failed for CohQ2Lim. "
773  << "Assuming Q2min = 0.0";
774  }
775  Q2.min = TMath::Max(0., A * B - C);
776  } else {
777  SLOG("KineLimits", pERROR)
778  << "Q2 kinematic limits calculation failed for CohQ2Lim. "
779  << "Assuming Q2min = 0.0";
780  }
781 
782  return Q2;
783 }
784 //____________________________________________________________________________
785 Range1D_t genie::utils::kinematics::Cohq2Lim(double Mn, double m_produced, double mlep, double Ev)
786 {
787  Range1D_t Q2 = utils::kinematics::CohQ2Lim(Mn, m_produced, mlep, Ev);
788  Range1D_t q2;
789  q2.min = - Q2.max;
790  q2.max = - Q2.min;
791  return q2;
792 }
793 //____________________________________________________________________________
794 Range1D_t genie::utils::kinematics::CohW2Lim(double Mn, double m_produced, double mlep,
795  double Ev, double Q2)
796 {
797  // These expressions for W^2 min and max appear in PRD 74, 054007 (2006) by
798  // Kartavtsev, Paschos, and Gounaris
799 
800  // That expression is specified for the case of the pion.
801  // GENIE has also Coherent interaction with Single gamma production and Rho
802  // so that formula will be adapted for a generic m_produced
803 
804  Range1D_t W2l;
805  W2l.min = -1;
806  W2l.max = -1;
807 
808  double s = Mn * Mn + 2.0 * Mn * Ev;
809  double Mnterm = 1 - Mn * Mn / s;
810  double Mlterm = 1 - mlep * mlep / s;
811  // Here T1, T2 are generically "term 1" and "term 2" in a long expression
812  double T1 = 0.25 * s * s * Mnterm * Mnterm * Mlterm;
813  double T2 = Q2 - (0.5 * s * Mnterm) + (0.5 * mlep * mlep * Mnterm);
814 
815  W2l.min = CohW2Min(Mn, m_produced);
816  W2l.max = (T1 - T2 * T2 ) *
817  (1.0 / Mnterm) *
818  (1.0 / (Q2 + mlep * mlep));
819 
820  return W2l;
821 }
822 //____________________________________________________________________________
824  double Q2, double Mn, double xsi)
825 {
826  Range1D_t nul;
827  nul.min = -1;
828  nul.max = -1;
829 
830  double nu_min = (W2min + Q2 - Mn * Mn) / (2.0 * Mn);
831  double nu_max = (W2max + Q2 - Mn * Mn) / (2.0 * Mn);
832  double xsiQ = xsi * TMath::Sqrt(Q2);
833 
834  nul.min = (xsiQ > nu_min) ? xsiQ : nu_min;
835  nul.max = nu_max;
836 
837  return nul;
838 }
839 //____________________________________________________________________________
840 Range1D_t genie::utils::kinematics::CohYLim(double Mn, double m_produced, double mlep,
841  double Ev, double Q2, double xsi)
842 {
843  Range1D_t ylim;
844  ylim.min = -1;
845  ylim.max = -1;
846 
847  Range1D_t W2lim = genie::utils::kinematics::CohW2Lim(Mn, m_produced, mlep, Ev, Q2);
848  if (W2lim.min > W2lim.max) {
849  LOG("KineLimits", pDEBUG)
850  << "Kinematically forbidden region in CohYLim. W2min = " << W2lim.min
851  << "; W2max =" << W2lim.max;
852  LOG("KineLimits", pDEBUG)
853  << " Mn = " << Mn << "; m_had_sys = " << m_produced << "; mlep = "
854  << mlep << "; Ev = " << Ev << "; Q2 = " << Q2;
855  return ylim;
856  }
858  Q2, Mn, xsi);
859  ylim.min = nulim.min / Ev;
860  ylim.max = nulim.max / Ev;
861 
862  return ylim;
863 }
864 //____________________________________________________________________________
866 {
867 // Computes y limits for coherent v interactions
868 
870  1.-ml/EvL - controls::kASmallNum);
871  return y;
872 }
873 //____________________________________________________________________________
874 double genie::utils::kinematics::CohW2Min(double Mn, double m_produced)
875 {
876  // These expressions for W^2 min and max appear in PRD 74, 054007 (2006) by
877  // Kartavtsev, Paschos, and Gounaris
878 
879  // That expression is specified for the case of the pion.
880  // GENIE has also Coherent interaction with Single gamma production and Rho
881  // so that formula will be adapted for a generic m_produced
882 
883  return (Mn + m_produced) * (Mn + m_produced);
884 }
885 //____________________________________________________________________________
887 {
889  return Q2;
890 }
891 //____________________________________________________________________________
892 Range1D_t genie::utils::kinematics::DarkWLim(double Ev, double M, double ml)
893 {
894 // Computes W limits for inelastic v interactions
895 //
896  double M2 = TMath::Power(M,2);
897  double ml2 = TMath::Power(ml,2);
898  double s = M2 + 2*M*Ev + ml2;
899  assert (s>0);
900 
901  Range1D_t W;
903 // W.min = kNeutronMass + kPionMass;
904  W.max = TMath::Sqrt(s) - ml;
905  if(W.max<=W.min) {
906  W.min = -1;
907  W.max = -1;
908  return W;
909  }
912  return W;
913 }
914 //____________________________________________________________________________
916  double Ev, double M, double ml, double W, double Q2min_cut)
917 {
918 // Computes Q2 limits (>0) @ the input W for inelastic v interactions
919 
920  Range1D_t Q2;
921  Q2.min = -1;
922  Q2.max = -1;
923 
924  double M2 = TMath::Power(M, 2.);
925  double ml2 = TMath::Power(ml, 2.);
926  double W2 = TMath::Power(W, 2.);
927  double s = M2 + 2*M*Ev + ml2;
928  assert(s > 0);
929  double sqs = TMath::Sqrt(s);
930  double E1CM = (s + ml2 - M2) / (2.*sqs);
931  double p1CM = TMath::Max(0., E1CM*E1CM - ml2);
932  p1CM = TMath::Sqrt(p1CM);
933  double E3CM = (s + ml2 - W2) / (2.*sqs);
934  double p3CM = TMath::Max(0., E3CM*E3CM - ml2);
935  p3CM = TMath::Sqrt(p3CM);
936 
937  SLOG("KineLimits", pDEBUG) << "s = " << s;
938  SLOG("KineLimits", pDEBUG) << "Ev = " << Ev;
939  SLOG("KineLimits", pDEBUG) << "M = " << M;
940  SLOG("KineLimits", pDEBUG) << "W = " << W;
941  SLOG("KineLimits", pDEBUG) << "E1_CM = " << E1CM;
942  SLOG("KineLimits", pDEBUG) << "p1_CM = " << p1CM;
943  SLOG("KineLimits", pDEBUG) << "E3_CM = " << E3CM;
944  SLOG("KineLimits", pDEBUG) << "p3_CM = " << p3CM;
945 
946  Q2.min = TMath::Power(p3CM - p1CM,2) - TMath::Power((W2 - M2) / (2.*sqs),2);
947  Q2.max = TMath::Power(p3CM + p1CM,2) - TMath::Power((W2 - M2) / (2.*sqs),2);
948 
949  SLOG("KineLimits", pDEBUG) << "Nominal Q^2 limits: " << Q2.min << " , " << Q2.max;
950  // guard against overflows
951  Q2.max = TMath::Max(0., Q2.max);
952  Q2.min = TMath::Max(0., Q2.min);
953 
954  // limit the minimum Q2
955  if(Q2.min < Q2min_cut) {Q2.min = Q2min_cut; }
956  if(Q2.max < Q2.min ) {Q2.min = -1; Q2.max = -1;}
957 
958  return Q2;
959 }
960 //____________________________________________________________________________
962  double Ev, double M, double ml, double W, double q2min_cut)
963 {
964 // Computes q2 (<0) limits @ the input W for inelastic v interactions
965 
966  Range1D_t Q2 = utils::kinematics::DarkQ2Lim_W(Ev,M,ml,W,-1.*q2min_cut);
967  Range1D_t q2;
968  q2.min = - Q2.max;
969  q2.max = - Q2.min;
970  return q2;
971 }
972 //____________________________________________________________________________
974  double Ev, double M, double ml, double Q2min_cut)
975 {
976 // Computes Q2 (>0) limits irrespective of W for inelastic v interactions
977 
978  Range1D_t Q2;
979  Q2.min = -1;
980  Q2.max = -1;
981 
983  if(W.min<0) return Q2;
984 
985  Q2 = utils::kinematics::DarkQ2Lim_W(Ev,M,ml,W.min,Q2min_cut);
986  return Q2;
987 }
988 //____________________________________________________________________________
990  double Ev, double M, double ml, double q2min_cut)
991 {
992 // Computes Q2 (>0) limits irrespective of W for inelastic v interactions
993 
994  Range1D_t Q2 = utils::kinematics::DarkQ2Lim(Ev,M,ml,-1.*q2min_cut);
995  Range1D_t q2;
996  q2.min = - Q2.max;
997  q2.max = - Q2.min;
998  return q2;
999 }
1000 //____________________________________________________________________________
1001 Range1D_t genie::utils::kinematics::DarkXLim(double Ev, double M, double ml)
1002 
1003 {
1004 // Computes Bjorken x limits for inelastic interactions
1005 // For the dark matter case, it is relatively straightforward
1006  Range1D_t Wl = utils::kinematics::DarkWLim(Ev, M, ml);
1007  double Wmin = Wl.min;
1008  double W2min = Wmin*Wmin;
1009  SLOG("KineLimits", pDEBUG) << "W^2_min = " << W2min;
1010  Range1D_t Q2l = utils::kinematics::DarkQ2Lim_W(Ev, M, ml, Wmin);
1011  SLOG("KineLimits", pDEBUG) << "Q^2 range : " << Q2l.min << " , " << Q2l.max;
1012  double M2 = M*M;
1013  Range1D_t x;
1014  x.min = Q2l.min / (Q2l.min + W2min - M2);
1015  x.max = Q2l.max / (Q2l.max + W2min - M2);
1016 
1017  SLOG("KineLimits", pDEBUG) << "x = [" << x.min << ", " << x.max << "]";
1018  return x;
1019 }
1020 //____________________________________________________________________________
1021 Range1D_t genie::utils::kinematics::DarkYLim(double Ev, double M, double ml)
1022 {
1023  // For dark inelastic scattering, can compute exactly and is much simpler
1024  Range1D_t Wl = utils::kinematics::DarkWLim(Ev, M, ml);
1025  double Wmin = Wl.min;
1026  double W2min = Wmin*Wmin;
1027  Range1D_t Q2l = utils::kinematics::DarkQ2Lim_W(Ev, M, ml, Wmin);
1028  double M2 = M*M;
1029  Range1D_t y;
1030  y.min = (Q2l.min + W2min - M2) / (2*Ev*M);
1031  y.max = (Q2l.max + W2min - M2) / (2*Ev*M);
1032 
1033  SLOG("KineLimits", pDEBUG) << "y = [" << y.min << ", " << y.max << "]";
1034  return y;
1035 }
1036 //____________________________________________________________________________
1038  double Ev, double M, double ml, double x)
1039 
1040 {
1041  // Computes y limits @ the input x for inelastic interactions
1042  // We hit y_min when W = W_min and y_max when Q^2 = Q_min^2 or Q^2 = Q_max^2
1043 
1044  Range1D_t y;
1045  y.min = -1;
1046  y.max = -1;
1047 
1048  Range1D_t Wl = utils::kinematics::DarkWLim(Ev, M, ml);
1049  double Wmin = Wl.min;
1050  double W2min = Wmin*Wmin;
1051  double M2 = M*M;
1052  double ml2 = ml*ml;
1053  y.min = (W2min - M2) / (1.-x) / (2*Ev*M);
1054  y.max = 2.* M * x *(Ev*Ev - ml2) / Ev / (2. * M * Ev * x + M2 * x * x + ml2);
1055 
1056  return y;
1057 }
1058 //____________________________________________________________________________
1059 // Kinematical Transformations:
1060 //____________________________________________________________________________
1062 {
1063 // Q2 -> QD2 transformation. QD2 takes out the dipole form of the form factors
1064 // making the differential cross section to be flatter and speeding up the
1065 // kinematical selection.
1066 
1067  assert(Q2>0);
1068  return TMath::Power(1+Q2/controls::kMQD2, -1);
1069 }
1070 //____________________________________________________________________________
1072 {
1073  assert(QD2>0);
1074  return controls::kMQD2*(1/QD2-1);
1075 }
1076 //____________________________________________________________________________
1077 double genie::utils::kinematics::Q2(const Interaction * const interaction)
1078 {
1079 // Get Q^2 from kinematics object
1080 // If Q^2 is not set but x,y are, then compute Q2 from x,y
1081 
1082  const Kinematics & kinematics = interaction->Kine();
1083 
1084  if (kinematics.KVSet(kKVQ2) || kinematics.KVSet(kKVq2)) {
1085  double Q2 = kinematics.Q2();
1086  return Q2;
1087  }
1088  if (kinematics.KVSet(kKVy)) {
1089  const InitialState & init_state = interaction->InitState();
1090  double Mn = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
1091  double x = kinematics.x();
1092  double y = kinematics.y();
1093  double Ev = init_state.ProbeE(kRfHitNucRest);
1094  double Q2 = 2*Mn*Ev*x*y;
1095  return Q2;
1096  }
1097  SLOG("KineLimits", pERROR) << "Couldn't compute Q^2 for \n"<< *interaction;
1098  return 0;
1099 }
1100 //____________________________________________________________________________
1101 double genie::utils::kinematics::W(const Interaction * const interaction)
1102 {
1103  const ProcessInfo & process_info = interaction->ProcInfo();
1104 
1105  if(process_info.IsQuasiElastic()) {
1106  // hadronic inv. mass is equal to the recoil nucleon on-shell mass
1107  int rpdgc = interaction->RecoilNucleonPdg();
1108  double M = PDGLibrary::Instance()->Find(rpdgc)->Mass();
1109  return M;
1110  }
1111 
1112  const Kinematics & kinematics = interaction->Kine();
1113  if(kinematics.KVSet(kKVW)) {
1114  double W = kinematics.W();
1115  return W;
1116  }
1117  if(kinematics.KVSet(kKVx) && kinematics.KVSet(kKVy)) {
1118  const InitialState & init_state = interaction->InitState();
1119  double Ev = init_state.ProbeE(kRfHitNucRest);
1120  double M = init_state.Tgt().HitNucP4Ptr()->M();
1121  double M2 = M*M;
1122  double x = kinematics.x();
1123  double y = kinematics.y();
1124  double W2 = TMath::Max(0., M2 + 2*Ev*M*y*(1-x));
1125  double W = TMath::Sqrt(W2);
1126  return W;
1127  }
1128  SLOG("KineLimits", pERROR) << "Couldn't compute W for \n"<< *interaction;
1129  return 0;
1130 }
1131 //___________________________________________________________________________
1133  double Ev, double M, double W, double Q2, double & x, double & y)
1134 {
1135 // Converts (W,Q2) => (x,y)
1136 // Uses the system: a) W^2 - M^2 = 2*Ev*M*y*(1-x) and b) Q^2 = 2*x*y*M*Ev
1137 // Ev is the neutrino energy at the struck nucleon rest frame
1138 // M is the nucleon mass - it does not need to be on the mass shell
1139 
1140  double M2 = TMath::Power(M,2);
1141  double W2 = TMath::Power(W,2);
1142 
1143  x = Q2 / (W2-M2+Q2);
1144  y = (W2-M2+Q2) / (2*M*Ev);
1145 
1146  x = TMath::Min(1.,x);
1147  y = TMath::Min(1.,y);
1148  x = TMath::Max(0.,x);
1149  y = TMath::Max(0.,y);
1150 
1151  LOG("KineLimits", pDEBUG)
1152  << "(W=" << W << ",Q2=" << Q2 << ") => (x="<< x << ", y=" << y<< ")";
1153 }
1154 //___________________________________________________________________________
1156  double Ev, double M, double & W, double & Q2, double x, double y)
1157 {
1158 // Converts (x,y) => (W,Q2)
1159 // Uses the system: a) W^2 - M^2 = 2*Ev*M*y*(1-x) and b) Q^2 = 2*x*y*M*Ev
1160 // Ev is the neutrino energy at the struck nucleon rest frame
1161 // M is the nucleon mass - it does not need to be on the mass shell
1162 
1163  double M2 = TMath::Power(M,2);
1164  double W2 = M2 + 2*Ev*M*y*(1-x);
1165 
1166  W = TMath::Sqrt(TMath::Max(0., W2));
1167  Q2 = 2*x*y*M*Ev;
1168 
1169  LOG("KineLimits", pDEBUG)
1170  << "(x=" << x << ",y=" << y << " => (W=" << W << ",Q2=" << Q2 << ")";
1171 }
1172 //___________________________________________________________________________
1174  double Ev, double M, double & W, double Q2, double x, double & y)
1175 {
1176 // Converts (x,Q2) => (W,Y)
1177 // Uses the system: a) W^2 - M^2 = 2*Ev*M*y*(1-x) and b) Q^2 = 2*x*y*M*Ev
1178 // Ev is the neutrino energy at the struck nucleon rest frame
1179 // M is the nucleon mass - it does not need to be on the mass shell
1180 
1181  double M2 = TMath::Power(M,2);
1182  y = Q2 / (2 * x * M * Ev);
1183  y = TMath::Min(1.,y);
1184 
1185  double W2 = M2 + 2*Ev*M*y*(1-x);
1186  W = TMath::Sqrt(TMath::Max(0., W2));
1187 
1188  LOG("KineLimits", pDEBUG)
1189  << "(x=" << x << ",Q2=" << Q2 << " => (W=" << W << ",Y=" << y << ")";
1190 }
1191 //___________________________________________________________________________
1193  double Ev, double M, double x, double y)
1194 {
1195 // Converts (x,y) => W
1196 // Ev is the neutrino energy at the struck nucleon rest frame
1197 // M is the nucleon mass - it does not need to be on the mass shell
1198 
1199  double M2 = TMath::Power(M,2);
1200  double W2 = M2 + 2*Ev*M*y*(1-x);
1201  double W = TMath::Sqrt(TMath::Max(0., W2));
1202 
1203  LOG("KineLimits", pDEBUG) << "(x=" << x << ",y=" << y << ") => W=" << W;
1204 
1205  return W;
1206 }
1207 //___________________________________________________________________________
1209  double Ev, double M, double x, double y)
1210 {
1211 // Converts (x,y) => Q2
1212 // Ev is the neutrino energy at the struck nucleon rest frame
1213 // M is the nucleon mass - it does not need to be on the mass shell
1214 
1215  double Q2 = 2*x*y*M*Ev;
1216 
1217  LOG("KineLimits", pDEBUG) << "(x=" << x << ",y=" << y << ") => Q2=" << Q2;
1218 
1219  return Q2;
1220 }
1221 //___________________________________________________________________________
1223  double Ev, double M, double Q2, double y)
1224 {
1225 // Converts (Q^2,y) => x
1226 // Ev is the neutrino energy at the struck nucleon rest frame
1227 // M is the nucleon mass - it does not need to be on the mass shell
1228  assert(Ev > 0. && M > 0. && Q2 > 0. && y > 0.);
1229 
1230  double x = Q2 / (2. * y * M * Ev);
1231 
1232  LOG("KineLimits", pDEBUG) << "(Ev=" << Ev << ",Q2=" << Q2
1233  << ",y=" << y << ",M=" << M << ") => x=" << x;
1234 
1235  return x;
1236 }
1237 //___________________________________________________________________________
1239  double x, double Q2, double M, double mc)
1240 {
1241 // x : scaling variable (plain or modified)
1242 // Q2: momentum transfer
1243 // M : hit nucleon "mass" (nucleon can be off the mass shell)
1244 // mc: charm mass
1245 //
1246  double M2 = TMath::Power(M,2);
1247  double v = 0.5*Q2/(M*x);
1248  double W2 = TMath::Max(0., M2+2*M*v-Q2);
1249  double W = TMath::Sqrt(W2);
1250  double Wmin = M + kLightestChmHad;
1251  double xc = utils::kinematics::SlowRescalingVar(x,Q2,M,mc);
1252 
1253  if(xc>=1 || W<=Wmin) return false;
1254  else return true;
1255 }
1256 //____________________________________________________________________________
1258  double x, double Q2, double M, double mc)
1259 {
1260 // x : scaling variable (plain or modified)
1261 // Q2: momentum transfer
1262 // M : hit nucleon "mass" (nucleon can be off the mass shell)
1263 // mc: charm mass
1264 
1265  double mc2 = TMath::Power(mc,2);
1266  double v = 0.5*Q2/(M*x);
1267  double xc = x + 0.5*mc2/(M*v);
1268  return xc;
1269 }
1270 //____________________________________________________________________________
1272  Range1D_t & range, double min_cut, double max_cut)
1273 {
1274  // if the min,max are within the existing limits, the cut can be applied
1275  // by narrowing down the xisting limits
1276  if ( utils::math::IsWithinLimits(min_cut, range ) ) range.min = min_cut;
1277  if ( utils::math::IsWithinLimits(max_cut, range ) ) range.max = max_cut;
1278 
1279  // if the min-cut is above the existing max-limit or
1280  // if the max-cut is below the existing min-limit then
1281  // the range should be invalidated
1282 
1283  if (min_cut > range.max || max_cut < range.min) {
1284 
1285  range.min = 0;
1286  range.max = 0;
1287  }
1288 }
1289 //___________________________________________________________________________
1291 {
1292  Kinematics * kine = in->KinePtr();
1293 
1294  if(kine->KVSet(kKVx) && kine->KVSet(kKVy)) {
1295  const InitialState & init_state = in->InitState();
1296  double Ev = init_state.ProbeE(kRfHitNucRest);
1297  double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off mass shell
1298  double x = kine->x();
1299  double y = kine->y();
1300 
1301  double Q2=-1,W=-1;
1302  kinematics::XYtoWQ2(Ev,M,W,Q2,x,y);
1303  kine->SetQ2(Q2);
1304  kine->SetW(W);
1305  }
1306 }
1307 //___________________________________________________________________________
1309 {
1310  Kinematics * kine = in->KinePtr();
1311 
1312  if(kine->KVSet(kKVW) && kine->KVSet(kKVQ2)) {
1313  const InitialState & init_state = in->InitState();
1314  double Ev = init_state.ProbeE(kRfHitNucRest);
1315  double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off mass shell
1316  double W = kine->W();
1317  double Q2 = kine->Q2();
1318 
1319  double x=-1,y=-1;
1320  kinematics::WQ2toXY(Ev,M,W,Q2,x,y);
1321  kine->Setx(x);
1322  kine->Sety(y);
1323  }
1324 }
1325 //___________________________________________________________________________
1327 {
1328  Kinematics * kine = in->KinePtr();
1329 
1330  if(kine->KVSet(kKVx) && kine->KVSet(kKVQ2)) {
1331  const InitialState & init_state = in->InitState();
1332  double Ev = init_state.ProbeE(kRfHitNucRest);
1333  double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off mass shell
1334  double x = kine->x();
1335  double Q2 = kine->Q2();
1336 
1337  double W=-1,y=-1;
1338  kinematics::XQ2toWY(Ev,M,W,Q2,x,y);
1339  kine->SetW(W);
1340  kine->Sety(y);
1341  }
1342 }
1343 //___________________________________________________________________________
1345 {
1346  Kinematics * kine = in->KinePtr();
1347 
1348  if(kine->KVSet(kKVy) && kine->KVSet(kKVQ2)) {
1349 
1350  const ProcessInfo & pi = in->ProcInfo();
1351  const InitialState & init_state = in->InitState();
1352  double M = 0.0;
1353  double Ev = 0.0;
1354 
1355  if (pi.IsCoherentProduction()) {
1356  M = in->InitState().Tgt().Mass(); // nucleus mass
1357  Ev = init_state.ProbeE(kRfLab);
1358  }
1359  else {
1360  M = in->InitState().Tgt().HitNucP4Ptr()->M(); //can be off m/shell
1361  Ev = init_state.ProbeE(kRfHitNucRest);
1362  }
1363 
1364  double y = kine->y();
1365  double Q2 = kine->Q2();
1366 
1367  double x = kinematics::Q2YtoX(Ev,M,Q2,y);
1368  kine->Setx(x);
1369  }
1370 }
1371 //____________________________________________________________________________
1373  double * x, double * par)
1374 {
1375  //-- inputs
1376  double Q2 = x[0]; // Q2
1377  double W = x[1]; // invariant mass
1378 
1379  //-- parameters
1380  double mD = par[0]; // resonance mass
1381  double gD = par[1]; // resonance width
1382  double xsmax = par[2]; // safety factor * max cross section in (W,Q2)
1383  double Wmax = par[3]; // kinematically allowed Wmax
1384  //double E = par[4]; // neutrino energy
1385 
1386 // return 3*xsmax; ////////////////NOTE
1387 
1388  xsmax*=5;
1389 
1390  double func = 0;
1391 
1392  if(Wmax > mD) {
1393  // -- if the resonance mass is within the kinematical range then
1394  // -- the envelope is defined as a plateau above the resonance peak,
1395  // -- a steeply falling leading edge (low-W side) and a more slowly
1396  // -- falling trailing edge (high-W side)
1397 
1398  double hwfe = mD+2*gD; // high W falling edge
1399  double lwfe = mD-gD/2; // low W falling edge
1400 
1401  if(W < lwfe) {
1402  //low-W falling edge
1403  func = xsmax / (1 + 5* TMath::Power((W-lwfe)/gD,2));
1404  } else if (W > hwfe) {
1405  //high-W falling edge
1406  func = xsmax / (1 + TMath::Power((W-hwfe)/gD,2));
1407  } else {
1408  // plateau
1409  func = xsmax;
1410  }
1411  } else {
1412  // -- if the resonance mass is above the kinematical range then the
1413  // -- envelope is a small plateau just bellow Wmax and a falling edge
1414  // -- at lower W
1415 
1416  double plateau_edge = Wmax-0.1;
1417 
1418  if (W > plateau_edge) {
1419  // plateau
1420  func = xsmax;
1421  } else {
1422  //low-W falling edge
1423  func = xsmax / (1 + TMath::Power((W-plateau_edge)/gD,2));
1424  }
1425  }
1426 
1427  if(Q2>controls::kMQD2) {
1428  func *= (0.5 + 1./(1 + Q2/controls::kMQD2));
1429  }
1430 
1431  return func;
1432 }
1433 //___________________________________________________________________________
1435  double * x, double * par)
1436 {
1437  //-- inputs
1438  double xb = x[0]; // scaling variable x brorken
1439  //double y = x[1]; // inelasticity y
1440 
1441  //-- parameters
1442  double xpeak = par[0]; // x peak
1443  double xmin = par[1]; // min x
1444  double xmax = 1.; // max x
1445  double xsmax = par[2]; // safety factor * max cross section in (x,y)
1446 
1447  double func = 0;
1448 
1449  if(xb < xpeak/20.) {
1450  //low-x falling edge
1451  double plateau_edge = xpeak/20.;
1452  double slope = xsmax/(xmin - plateau_edge);
1453  func = xsmax - slope * (xb - plateau_edge);
1454  } else if (xb > 2*xpeak) {
1455  //high-x falling edge
1456  double plateau_edge = 2*xpeak;
1457  double slope = xsmax/(xmax - plateau_edge);
1458  func = xsmax - slope * (xb - plateau_edge);
1459  } else {
1460  // plateau
1461  func = xsmax;
1462  }
1463  return func;
1464 }
1465 //___________________________________________________________________________
1467  double * x, double * par)
1468 {
1469  //-- inputs
1470  double xb = x[0]; // x
1471  double yb = x[1]; // y
1472 
1473  //-- parameters
1474  double xsmax = 3*par[0]; // safety factor * max cross section in (x,y)
1475  double Ev = par[1]; // neutrino energy;
1476 
1477  if(yb<0.|| yb>1.) return 0.;
1478  if(xb<0.|| xb>1.) return 0.;
1479 
1480  if(Ev<1) return xsmax;
1481  if(xb/Ev<1E-4 && yb>0.95) return 5*xsmax;
1482 
1483  double func = 0;
1484  double xp = 0.1;
1485  double yp = (Ev>2.5) ? 2.5/Ev : 1;
1486 
1487  if(xb>xp) {
1488  double xs0=0;
1489  if(yb<yp) {
1490  xs0 = xsmax;
1491  } else {
1492  xs0 = xsmax - (yb-yp)*xsmax;
1493  }
1494  double d = TMath::Power( (xb-xp)/0.075, 2);
1495  func = xs0/(1 + d);
1496  } else {
1497  if(yb>yp) {
1498  func = xsmax - (yb-yp)*xsmax;
1499  } else {
1500  func = xsmax;
1501  }
1502  }
1503  return func;
1504 }
1505 //___________________________________________________________________________
double COHImportanceSamplingEnvelope(double *x, double *par)
Definition: KineUtils.cxx:1466
static const double kMaxX
Definition: Controls.h:44
double W(bool selected=false) const
Definition: Kinematics.cxx:157
const KPhaseSpace & PhaseSpace(void) const
Definition: Interaction.h:73
bool TransformMatched(KinePhaseSpace_t ia, KinePhaseSpace_t ib, KinePhaseSpace_t a, KinePhaseSpace_t b, bool &fwd)
Definition: KineUtils.cxx:341
double J(double q0, double q3, double Enu, double ml)
Definition: MECUtils.cxx:147
Range1D_t InelXLim(double El, double ml, double M)
Definition: KineUtils.cxx:647
Range1D_t InelWLim(double El, double ml, double M)
Definition: KineUtils.cxx:547
#define pERROR
Definition: Messenger.h:59
void XQ2toWY(double Ev, double M, double &W, double Q2, double x, double &y)
Definition: KineUtils.cxx:1173
void UpdateXYFromWQ2(const Interaction *in)
Definition: KineUtils.cxx:1308
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1077
Range1D_t Inelq2Lim_W(double El, double ml, double M, double W)
Definition: KineUtils.cxx:607
double DISImportanceSamplingEnvelope(double *x, double *par)
Definition: KineUtils.cxx:1434
void SetQ2(double Q2, bool selected=false)
Definition: Kinematics.cxx:255
Kinematics * KinePtr(void) const
Definition: Interaction.h:76
Range1D_t CohW2Lim(double Mn, double m_produced, double mlep, double Ev, double Q2)
Definition: KineUtils.cxx:794
Range1D_t InelWLim(double Ev, double M, double ml)
Definition: KineUtils.cxx:358
int RecoilNucleonPdg(void) const
recoil nucleon pdg
A simple [min,max] interval for doubles.
Definition: Range1.h:42
Range1D_t Darkq2Lim(double Ev, double M, double ml, double q2min_cut=-1 *controls::kMinQ2Limit)
Definition: KineUtils.cxx:989
bool IsQuasiElastic(void) const
Definition: ProcessInfo.cxx:69
static const double kMaxY
Definition: Controls.h:46
#define pFATAL
Definition: Messenger.h:56
Range1D_t DarkQ2Lim_W(double Ev, double M, double ml, double W, double Q2min_cut=controls::kMinQ2Limit)
Definition: KineUtils.cxx:915
bool IsWithinLimits(double x, Range1D_t range)
Definition: MathUtils.cxx:258
static constexpr double s
Definition: Units.h:95
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
TParticlePDG * Probe(void) const
double x(bool selected=false) const
Definition: Kinematics.cxx:99
static const double kLightestChmHad
static const double kMinY
Definition: Controls.h:45
Range1D_t CEvNSQ2Lim(double Ev)
Definition: KineUtils.cxx:886
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)
Definition: KineUtils.cxx:429
enum genie::EKinePhaseSpace KinePhaseSpace_t
double Mass(void) const
Definition: Target.cxx:224
Range1D_t Limits(KineVar_t kvar) const
Return the kinematical variable limits.
void UpdateWYFromXQ2(const Interaction *in)
Definition: KineUtils.cxx:1326
double XYtoW(double Ev, double M, double x, double y)
Definition: KineUtils.cxx:1192
double PhaseSpaceVolume(const Interaction *const i, KinePhaseSpace_t ps)
Definition: KineUtils.cxx:36
Range1D_t InelQ2Lim_W(double Ev, double M, double ml, double W, double Q2min_cut=controls::kMinQ2Limit)
Definition: KineUtils.cxx:379
double y(bool selected=false) const
Definition: Kinematics.cxx:112
static constexpr double b
Definition: Units.h:78
double W(const Interaction *const i)
Definition: KineUtils.cxx:1101
Summary information for an interaction.
Definition: Interaction.h:56
const TLorentzVector & HitNucP4(void) const
Definition: Target.h:91
double SlowRescalingVar(double x, double Q2, double M, double mc)
Definition: KineUtils.cxx:1257
Range1D_t InelXLim(double Ev, double M, double ml)
Definition: KineUtils.cxx:457
Range1D_t DarkQ2Lim(double Ev, double M, double ml, double Q2min_cut=controls::kMinQ2Limit)
Definition: KineUtils.cxx:973
Range1D_t Inelq2Lim_W(double Ev, double M, double ml, double W, double q2min_cut=-1 *controls::kMinQ2Limit)
Definition: KineUtils.cxx:417
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
Range1D_t CohNuLim(double W2min, double W2max, double Q2, double Mn, double xsi)
Definition: KineUtils.cxx:823
bool KVSet(KineVar_t kv) const
Definition: Kinematics.cxx:317
double QD2toQ2(double QD2)
Definition: KineUtils.cxx:1071
Exception used inside Interaction classes.
Range1D_t Darkq2Lim_W(double Ev, double M, double ml, double W, double q2min_cut=-1 *controls::kMinQ2Limit)
Definition: KineUtils.cxx:961
static constexpr double A
Definition: Units.h:74
static string AsString(KinePhaseSpace_t kps)
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
const double a
double XYtoQ2(double Ev, double M, double x, double y)
Definition: KineUtils.cxx:1208
void XYtoWQ2(double Ev, double M, double &W, double &Q2, double x, double y)
Definition: KineUtils.cxx:1155
void WQ2toXY(double Ev, double M, double W, double Q2, double &x, double &y)
Definition: KineUtils.cxx:1132
Range1D_t Cohq2Lim(double Mn, double m_produced, double mlep, double Ev)
Definition: KineUtils.cxx:785
const Kinematics & Kine(void) const
Definition: Interaction.h:71
Kinematical phase space.
Definition: KPhaseSpace.h:33
static const double kNeutronMass
static const double kMinX
Definition: Controls.h:43
double func(double x, double y)
Range1D_t CohXLim(void)
Definition: KineUtils.cxx:735
static const double kASmallNum
Definition: Controls.h:40
double GetKV(KineVar_t kv) const
Definition: Kinematics.cxx:323
bool IsAboveCharmThreshold(double x, double Q2, double M, double mc)
Definition: KineUtils.cxx:1238
static constexpr double ps
Definition: Units.h:99
double Q2toQD2(double Q2)
Definition: KineUtils.cxx:1061
Range1D_t InelQ2Lim(double El, double ml, double M)
Definition: KineUtils.cxx:619
Range1D_t Inelq2Lim(double El, double ml, double M)
Definition: KineUtils.cxx:635
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
void Setx(double x, bool selected=false)
Definition: Kinematics.cxx:231
static const double kMQD2
Definition: Controls.h:65
void UpdateWQ2FromXY(const Interaction *in)
Definition: KineUtils.cxx:1290
Range1D_t InelYLim(double El, double ml, double M)
Definition: KineUtils.cxx:667
double max
Definition: Range1.h:53
void SetW(double W, bool selected=false)
Definition: Kinematics.cxx:279
TLorentzVector * HitNucP4Ptr(void) const
Definition: Target.cxx:247
Range1D_t InelQ2Lim_W(double El, double ml, double M, double W)
Definition: KineUtils.cxx:569
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:68
double CohW2Min(double Mn, double m_produced)
Definition: KineUtils.cxx:874
Range1D_t InelYLim_X(double Ev, double M, double ml, double x)
Definition: KineUtils.cxx:512
double Q2YtoX(double Ev, double M, double Q2, double y)
Definition: KineUtils.cxx:1222
Range1D_t CohYLim(double Mn, double m_produced, double mlep, double Ev, double Q2, double xsi)
Definition: KineUtils.cxx:840
Range1D_t DarkYLim(double Ev, double M, double ml)
Definition: KineUtils.cxx:1021
void Sety(double y, bool selected=false)
Definition: Kinematics.cxx:243
void UpdateXFromQ2Y(const Interaction *in)
Definition: KineUtils.cxx:1344
static const double kAVerySmallNum
Definition: Controls.h:39
Range1D_t InelYLim_X(double El, double ml, double M, double x)
Definition: KineUtils.cxx:702
Range1D_t CohQ2Lim(double Mn, double m_produced, double mlep, double Ev)
Definition: KineUtils.cxx:743
double Jacobian(const Interaction *const i, KinePhaseSpace_t f, KinePhaseSpace_t t)
Definition: KineUtils.cxx:130
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
double min
Definition: Range1.h:52
Range1D_t DarkYLim_X(double Ev, double M, double ml, double x)
Definition: KineUtils.cxx:1037
Range1D_t Inelq2Lim(double Ev, double M, double ml, double q2min_cut=-1 *controls::kMinQ2Limit)
Definition: KineUtils.cxx:445
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:86
double Q2(bool selected=false) const
Definition: Kinematics.cxx:125
const Target & Tgt(void) const
Definition: InitialState.h:66
#define SLOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a short string (using the FUNCTION and...
Definition: Messenger.h:84
Range1D_t InelYLim(double Ev, double M, double ml)
Definition: KineUtils.cxx:477
double ProbeE(RefFrame_t rf) const
Range1D_t DarkXLim(double Ev, double M, double ml)
Definition: KineUtils.cxx:1001
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)
Definition: KineUtils.cxx:1271
Range1D_t DarkWLim(double Ev, double M, double ml)
Definition: KineUtils.cxx:892
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63
double RESImportanceSamplingEnvelope(double *x, double *par)
Definition: KineUtils.cxx:1372