COHKinematicsGenerator.cxx
Go to the documentation of this file.
1 //____________________________________________________________________________
2 /*
3  Copyright (c) 2003-2019, The GENIE Collaboration
4  For the full text of the license visit http://copyright.genie-mc.org
5  or see $GENIE/LICENSE
6 
7  Author: Costas Andreopoulos <costas.andreopoulos \at stfc.ac.uk>
8  University of Liverpool & STFC Rutherford Appleton Lab
9 
10  For the class documentation see the corresponding header file.
11 
12  Important revisions after version 2.0.0 :
13  @ Feb 09, 2009 - CA
14  Moved into the new Coherent package from its previous location (EVGModules
15  package)
16  @ Mar 03, 2009 - CA
17  Renamed COHPiKinematicsGenerator -> COHKinematicsGenerator in
18  anticipation of reusing the code for simulating coherent production of
19  vector mesons.
20  @ May 06, 2009 - CA
21  Fix a problem with the search for the max cross section over the allowed
22  phase space which prevented kinematics to be generated for events near the
23  energy threshold.
24  @ Feb 06, 2013 - CA
25  When the value of the differential cross-section for the selected kinematics
26  is set to the event, set the corresponding KinePhaseSpace_t value too.
27 
28 */
29 //____________________________________________________________________________
30 
31 #include <cstdlib>
32 
33 #include <TROOT.h>
34 #include <TMath.h>
35 #include <TF2.h>
36 #include "Math/Minimizer.h"
37 #include "Math/Factory.h"
38 
55 
56 using namespace genie;
57 using namespace genie::constants;
58 using namespace genie::controls;
59 using namespace genie::utils;
60 
61 //___________________________________________________________________________
63  KineGeneratorWithCache("genie::COHKinematicsGenerator")
64 {
65  fEnvelope = 0;
66 }
67 //___________________________________________________________________________
69  KineGeneratorWithCache("genie::COHKinematicsGenerator", config)
70 {
71  fEnvelope = 0;
72 }
73 //___________________________________________________________________________
75 {
76  if(fEnvelope) delete fEnvelope;
77 }
78 //___________________________________________________________________________
80 {
81  if(fGenerateUniformly) {
82  LOG("COHKinematics", pNOTICE)
83  << "Generating kinematics uniformly over the allowed phase space";
84  }
85 
86  //-- Access cross section algorithm for running thread
88  const EventGeneratorI * evg = rtinfo->RunningThread();
89  fXSecModel = evg->CrossSectionAlg();
90  if (fXSecModel->Id().Name() == "genie::ReinSehgalCOHPiPXSec") {
92  } else if (fXSecModel->Id().Name() == "genie::BergerSehgalCOHPiPXSec2015") {
94  } else if (fXSecModel->Id().Name() == "genie::BergerSehgalFMCOHPiPXSec2015") {
96  } else if ((fXSecModel->Id().Name() == "genie::AlvarezRusoCOHPiPXSec")) {
98  }
99  else {
100  LOG("COHKinematicsGenerator",pFATAL) <<
101  "ProcessEventRecord >> Cannot calculate kinematics for " <<
102  fXSecModel->Id().Name();
103  }
104 }
105 //___________________________________________________________________________
107 {
108  // Get the Primary Interacton object
109  Interaction * interaction = evrec->Summary();
110  interaction->SetBit(kISkipProcessChk);
111  interaction->SetBit(kISkipKinematicChk); // TODO: Why turn this off?
112 
113  // Initialise a random number generator
115 
116  //-- For the subsequent kinematic selection with the rejection method:
117  // Calculate the max differential cross section or retrieve it from the
118  // cache. Throw an exception and quit the evg thread if a non-positive
119  // value is found.
120  //
121  // TODO: We are not offering the "fGenerateUniformly" option here.
122  double xsec_max = this->MaxXSec(evrec);
123 
124  //-- Get the kinematical limits for the generated x,y
125  const KPhaseSpace & kps = interaction->PhaseSpace();
126  Range1D_t y = kps.YLim();
128  assert(y.min>0. && y.max>0. && y.min<1. && y.max<1. && y.min<y.max);
129 
130  const double ymin = y.min + kASmallNum;
131  const double ymax = y.max - kASmallNum;
132  const double dy = ymax - ymin;
133  const double Q2min = Q2.min + kASmallNum;
134  const double Q2max = Q2.max - kASmallNum;
135  const double dQ2 = Q2max - Q2min;
136 
137  unsigned int iter = 0;
138  bool accept=false;
139  double xsec=-1, gy=-1, gQ2=-1;
140 
141  while(1) {
142  iter++;
143  if(iter > kRjMaxIterations) this->throwOnTooManyIterations(iter,evrec);
144 
145  //-- Select unweighted kinematics using importance sampling method.
146  // TODO: The importance sampling envelope is not used. Currently,
147  // we just employ a standard rejection-method approach.
148 
149  gy = ymin + dy * rnd->RndKine().Rndm();
150  gQ2 = Q2min + dQ2 * rnd->RndKine().Rndm();
151 
152  LOG("COHKinematics", pINFO) <<
153  "Trying: Q^2 = " << gQ2 << ", y = " << gy; /* << ", t = " << gt; */
154 
155  interaction->KinePtr()->Sety(gy);
156  interaction->KinePtr()->SetQ2(gQ2);
157  kinematics::UpdateXFromQ2Y(interaction);
158 
159  // computing cross section for the current kinematics
160  xsec = fXSecModel->XSec(interaction, kPSQ2yfE);
161 
162  //-- decide whether to accept the current kinematics
163  accept = (xsec_max * rnd->RndKine().Rndm() < xsec);
164 
165  //-- If the generated kinematics are accepted, finish-up module's job
166  if(accept) {
167  LOG("COHKinematics", pNOTICE)
168  << "Selected: Q^2 = " << gQ2 << ", y = " << gy; /* << ", t = " << gt; */
169 
170  // the Berger-Sehgal COH cross section should be a triple differential cross section
171  // d^2xsec/dQ2dydt where t is the the square of the 4p transfer to the
172  // nucleus. The cross section used for kinematical selection should have
173  // the t-dependence integrated out. The t-dependence is of the form
174  // ~exp(-bt). Now that the x,y kinematical variables have been selected
175  // we can generate a t using the t-dependence as a PDF.
176  const InitialState & init_state = interaction->InitState();
177  double Ev = init_state.ProbeE(kRfLab);
178  double Epi = gy*Ev; // pion energy
179  double Epi2 = TMath::Power(Epi,2);
180  double pme2 = kPionMass2/Epi2;
181  double gx = interaction->KinePtr()->x(); // utils::kinematics::Q2YtoX(Ev,kNucleonMass,gQ2,gy); // gQ2 / ( 2. * gy * kNucleonMass * Ev);
182  double xME = kNucleonMass*gx/Epi;
183  double tA = 1. + xME - 0.5*pme2;
184  /* Range1D_t tl = kps.TLim(); // this becomes the bounds for t */
185  double tB = TMath::Sqrt(1.+ 2*xME) * TMath::Sqrt(1.-pme2);
186  double tmin = 2*Epi2 * (tA-tB);
187  double tmax = 2*Epi2 * (tA+tB);
188  double A = (double) init_state.Tgt().A();
189  double A13 = TMath::Power(A,1./3.);
190  double R = fRo * A13 * units::fermi; // nuclear radius
191  double R2 = TMath::Power(R,2.);
192  double b = 0.33333 * R2;
193  double tsum = (TMath::Exp(-b*tmin) - TMath::Exp(-b*tmax))/b;
194  double rt = tsum * rnd->RndKine().Rndm();
195  double gt = -1.*TMath::Log(-1.*b*rt + TMath::Exp(-1.*b*tmin))/b;
196 
197  // TODO: If we re-install the fGenerateUniformly option, we
198  // would compute the event weight here.
199 
200  // reset bits
201  interaction->ResetBit(kISkipProcessChk);
202  interaction->ResetBit(kISkipKinematicChk);
203 
204  // lock selected kinematics & clear running values
205  interaction->KinePtr()->Setx(gx, true);
206  interaction->KinePtr()->Sety(gy, true);
207  interaction->KinePtr()->Sett(gt, true);
208  interaction->KinePtr()->SetW(this->pionMass(interaction), true);
209  interaction->KinePtr()->SetQ2(gQ2, true);
210  interaction->KinePtr()->ClearRunningValues();
211 
212  // set the cross section for the selected kinematics
213  evrec->SetDiffXSec(xsec * TMath::Exp(-b * gt) / tsum, kPSQ2yfE);
214 
215  return;
216  }
217  }// iterations
218 }
219 //___________________________________________________________________________
221 {
222  // Get the Primary Interacton object
223  Interaction * interaction = evrec->Summary();
224  interaction->SetBit(kISkipProcessChk);
225  interaction->SetBit(kISkipKinematicChk);
226 
227  // Initialise a random number generator
229 
230  //-- For the subsequent kinematic selection with the rejection method:
231  // Calculate the max differential cross section or retrieve it from the
232  // cache. Throw an exception and quit the evg thread if a non-positive
233  // value is found.
234  //
235  // TODO: We are not offering the "fGenerateUniformly" option here.
236  double xsec_max = this->MaxXSec(evrec);
237 
238  //-- Get the kinematical limits for the generated x,y
239  const KPhaseSpace & kps = interaction->PhaseSpace();
240  Range1D_t y = kps.YLim();
242  assert(y.min>0. && y.max>0. && y.min<1. && y.max<1. && y.min<y.max);
243 
244  const double ymin = y.min + kASmallNum;
245  const double ymax = y.max - kASmallNum;
246  const double dy = ymax - ymin;
247  const double Q2min = Q2.min + kASmallNum;
248  const double Q2max = Q2.max - kASmallNum;
249  const double dQ2 = Q2max - Q2min;
250  const double tmin = kASmallNum;
251  const double tmax = fTMax - kASmallNum; // TODO: Choose realistic t bounds
252  const double dt = tmax - tmin;
253 
254  //-- Try to select a valid (Q^2,y,t) triple.
255 
256  unsigned int iter = 0;
257  bool accept=false;
258  double xsec=-1, gy=-1, gt=-1, gQ2=-1;
259 
260  while(1) {
261  iter++;
262  if(iter > kRjMaxIterations) this->throwOnTooManyIterations(iter,evrec);
263 
264  //-- Select unweighted kinematics using importance sampling method.
265  // TODO: The importance sampling envelope is not used. Currently,
266  // we just employ a standard rejection-method approach.
267 
268  gy = ymin + dy * rnd->RndKine().Rndm();
269  gt = tmin + dt * rnd->RndKine().Rndm();
270  gQ2 = Q2min + dQ2 * rnd->RndKine().Rndm();
271 
272  LOG("COHKinematics", pINFO) <<
273  "Trying: Q^2 = " << gQ2 << ", y = " << gy << ", t = " << gt;
274 
275  interaction->KinePtr()->Sety(gy);
276  interaction->KinePtr()->Sett(gt);
277  interaction->KinePtr()->SetQ2(gQ2);
278 
279  // computing cross section for the current kinematics
280  xsec = fXSecModel->XSec(interaction, kPSxyfE);
281 
282  //-- decide whether to accept the current kinematics
283  accept = (xsec_max * rnd->RndKine().Rndm() < xsec);
284 
285  //-- If the generated kinematics are accepted, finish-up module's job
286  if(accept) {
287  LOG("COHKinematics", pNOTICE)
288  << "Selected: Q^2 = " << gQ2 << ", y = " << gy << ", t = " << gt;
289 
290  // TODO: If we re-install the fGenerateUniformly option, we
291  // would compute the event weight here.
292 
293  // reset bits
294  interaction->ResetBit(kISkipProcessChk);
295  interaction->ResetBit(kISkipKinematicChk);
296 
297  // lock selected kinematics & clear running values
298  interaction->KinePtr()->SetQ2(gQ2, true);
299  interaction->KinePtr()->Sety(gy, true);
300  interaction->KinePtr()->Sett(gt, true);
301  interaction->KinePtr()->SetW(this->pionMass(interaction), true);
302  interaction->KinePtr()->ClearRunningValues();
303 
304  // set the cross section for the selected kinematics
305  evrec->SetDiffXSec(xsec, kPSxytfE);
306 
307  return;
308  }
309  }// iterations
310 }
311 //___________________________________________________________________________
313 {
314  // Get the Primary Interacton object
315  Interaction * interaction = evrec->Summary();
316  interaction->SetBit(kISkipProcessChk);
317  interaction->SetBit(kISkipKinematicChk);
318 
319  //-- Get the random number generators
321 
322  //-- For the subsequent kinematic selection with the rejection method:
323  // Calculate the max differential cross section or retrieve it from the
324  // cache. Throw an exception and quit the evg thread if a non-positive
325  // value is found.
326  // If the kinematics are generated uniformly over the allowed phase
327  // space the max xsec is irrelevant
328  double xsec_max = (fGenerateUniformly) ? -1 : this->MaxXSec(evrec);
329 
330  //-- Get the kinematical limits for the generated x,y
331  const KPhaseSpace & kps = interaction->PhaseSpace();
332  Range1D_t y = kps.YLim();
333  assert(y.min>0. && y.max>0. && y.min<1. && y.max<1. && y.min<y.max);
334 
335  const double xmin = kASmallNum;
336  const double xmax = 1.- kASmallNum;
337  const double ymin = y.min + kASmallNum;
338  const double ymax = y.max - kASmallNum;
339  const double dx = xmax - xmin;
340  const double dy = ymax - ymin;
341 
342  //------ Try to select a valid x,y pair
343 
344  unsigned int iter = 0;
345  bool accept=false;
346  double xsec=-1, gx=-1, gy=-1;
347 
348  while(1) {
349  iter++;
350  if(iter > kRjMaxIterations) this->throwOnTooManyIterations(iter,evrec);
351 
352  if(fGenerateUniformly) {
353  //-- Generate a x,y pair uniformly in the kinematically allowed range.
354  gx = xmin + dx * rnd->RndKine().Rndm();
355  gy = ymin + dy * rnd->RndKine().Rndm();
356 
357  } else {
358  //-- Select unweighted kinematics using importance sampling method.
359 
360  if(iter==1) {
361  LOG("COHKinematics", pNOTICE) << "Initializing the sampling envelope";
362  double Ev = interaction->InitState().ProbeE(kRfLab);
363  fEnvelope->SetRange(xmin,ymin,xmax,ymax);
364  fEnvelope->SetParameter(0, xsec_max);
365  fEnvelope->SetParameter(1, Ev);
366  }
367 
368  // Generate W,QD2 using the 2-D envelope as PDF
369  fEnvelope->GetRandom2(gx,gy);
370  }
371 
372  LOG("COHKinematics", pINFO) << "Trying: x = " << gx << ", y = " << gy;
373 
374  interaction->KinePtr()->Setx(gx);
375  interaction->KinePtr()->Sety(gy);
376 
377  // computing cross section for the current kinematics
378  xsec = fXSecModel->XSec(interaction, kPSxyfE);
379 
380  //-- decide whether to accept the current kinematics
381  if(!fGenerateUniformly) {
382  double max = fEnvelope->Eval(gx, gy);
383  double t = max * rnd->RndKine().Rndm();
384 
385  this->AssertXSecLimits(interaction, xsec, max);
386 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
387  LOG("COHKinematics", pDEBUG)
388  << "xsec= " << xsec << ", J= 1, Rnd= " << t;
389 #endif
390  accept = (t<xsec);
391  }
392  else {
393  accept = (xsec>0);
394  }
395 
396  //-- If the generated kinematics are accepted, finish-up module's job
397  if(accept) {
398  LOG("COHKinematics", pNOTICE) << "Selected: x = "<< gx << ", y = "<< gy;
399 
400  // the Rein-Sehgal COH cross section should be a triple differential cross section
401  // d^2xsec/dxdydt where t is the the square of the 4p transfer to the
402  // nucleus. The cross section used for kinematical selection should have
403  // the t-dependence integrated out. The t-dependence is of the form
404  // ~exp(-bt). Now that the x,y kinematical variables have been selected
405  // we can generate a t using the t-dependence as a PDF.
406  const InitialState & init_state = interaction->InitState();
407  double Ev = init_state.ProbeE(kRfLab);
408  double Epi = gy*Ev; // pion energy
409  double Epi2 = TMath::Power(Epi,2);
410  double pme2 = kPionMass2/Epi2;
411  double xME = kNucleonMass*gx/Epi;
412  double tA = 1. + xME - 0.5*pme2;
413  double tB = TMath::Sqrt(1.+ 2*xME) * TMath::Sqrt(1.-pme2);
414  double tmin = 2*Epi2 * (tA-tB);
415  double tmax = 2*Epi2 * (tA+tB);
416  double A = (double) init_state.Tgt().A();
417  double A13 = TMath::Power(A,1./3.);
418  double R = fRo * A13 * units::fermi; // nuclear radius
419  double R2 = TMath::Power(R,2.);
420  double b = 0.33333 * R2;
421  double tsum = (TMath::Exp(-b*tmin) - TMath::Exp(-b*tmax))/b;
422  double rt = tsum * rnd->RndKine().Rndm();
423  double gt = -1.*TMath::Log(-1.*b*rt + TMath::Exp(-1.*b*tmin))/b;
424 
425  LOG("COHKinematics", pNOTICE)
426  << "Selected: t = "<< gt << ", from ["<< tmin << ", "<< tmax << "]";
427 
428  // for uniform kinematics, compute an event weight as
429  // wght = (phase space volume)*(differential xsec)/(event total xsec)
430  if(fGenerateUniformly) {
431  double vol = y.max-y.min; // dx=1, dt: irrelevant
432  double totxsec = evrec->XSec();
433  double wght = (vol/totxsec)*xsec;
434  LOG("COHKinematics", pNOTICE) << "Kinematics wght = "<< wght;
435 
436  // apply computed weight to the current event weight
437  wght *= evrec->Weight();
438  LOG("COHKinematics", pNOTICE) << "Current event wght = " << wght;
439  evrec->SetWeight(wght);
440  }
441 
442  // reset bits
443  interaction->ResetBit(kISkipProcessChk);
444  interaction->ResetBit(kISkipKinematicChk);
445 
446  // lock selected kinematics & clear running values
447  interaction->KinePtr()->Setx(gx, true);
448  interaction->KinePtr()->Sety(gy, true);
449  interaction->KinePtr()->Sett(gt, true);
450  interaction->KinePtr()->SetW(kPionMass, true);
451  interaction->KinePtr()->SetQ2(2*kNucleonMass*gx*gy*Ev, true);
452  interaction->KinePtr()->ClearRunningValues();
453 
454  // set the cross section for the selected kinematics
455  evrec->SetDiffXSec(xsec * TMath::Exp(-b * gt) / tsum, kPSxytfE);
456 
457  return;
458  }
459  }// iterations
460 }
461 //___________________________________________________________________________
463 {
464 
465  LOG("COHKinematics", pNOTICE) << "Using AlvarezRuso Model";
466  // Get the Primary Interacton object
467  Interaction * interaction = evrec->Summary();
468  interaction->SetBit(kISkipProcessChk);
469  interaction->SetBit(kISkipKinematicChk);
470 
471  // Initialise a random number generator
473 
474  //-- For the subsequent kinematic selection with the rejection method:
475  // Calculate the max differential cross section or retrieve it from the
476  // cache. Throw an exception and quit the evg thread if a non-positive
477  // value is found.
478  // If the kinematics are generated uniformly over the allowed phase
479  // space the max xsec is irrelevant
480  double xsec_max = (fGenerateUniformly) ? -1 : this->MaxXSec(evrec);
481 
482  //Set up limits of integration variables
483  // Primary lepton energy
484  const double E_l_min = interaction->FSPrimLepton()->Mass();
485  const double E_l_max = interaction->InitStatePtr()->GetProbeP4(kRfLab)->E() - kPionMass;
486  // Primary lepton angle with respect to the beam axis
487  const double ctheta_l_min = 0.4;
488  const double ctheta_l_max = 1.0 - kASmallNum;
489  // Pion angle with respect to the beam axis
490  const double ctheta_pi_min = 0.4;
491  const double ctheta_pi_max = 1.0 - kASmallNum;
492  // Pion angle transverse to the beam axis
493  const double phi_min = kASmallNum;
494  const double phi_max = (2.0 * kPi) - kASmallNum;
495  //
496  const double d_E_l = E_l_max - E_l_min;
497  const double d_ctheta_l = ctheta_l_max - ctheta_l_min;
498  const double d_ctheta_pi = ctheta_pi_max - ctheta_pi_min;
499  const double d_phi = phi_max - phi_min;
500 
501  //------ Try to select a valid set of kinematics
502  unsigned int iter = 0;
503  bool accept=false;
504  double xsec=-1, g_E_l=-1, g_theta_l=-1, g_phi_l=-1, g_theta_pi=-1, g_phi_pi=-1;
505  double g_ctheta_l, g_ctheta_pi;
506 
507  while(1) {
508  iter++;
509  if(iter > kRjMaxIterations) this->throwOnTooManyIterations(iter,evrec);
510 
511  //Select kinematic point
512  g_E_l = E_l_min + d_E_l * rnd->RndKine().Rndm();
513  g_ctheta_l = ctheta_l_min + d_ctheta_l * rnd->RndKine().Rndm();
514  g_ctheta_pi = ctheta_pi_min + d_ctheta_pi * rnd->RndKine().Rndm();
515  g_phi_l = phi_min + d_phi * rnd->RndKine().Rndm();
516  // random phi is relative to phi_l
517  g_phi_pi = g_phi_l + (phi_min + d_phi * rnd->RndKine().Rndm());
518  g_theta_l = TMath::ACos(g_ctheta_l);
519  g_theta_pi = TMath::ACos(g_ctheta_pi);
520 
521  LOG("COHKinematics", pINFO) << "Trying: Lep(" <<g_E_l << ", " <<
522  g_theta_l << ", " << g_phi_l << ") Pi(" <<
523  g_theta_pi << ", " << g_phi_pi << ")";
524 
525  this->SetKinematics(g_E_l, g_theta_l, g_phi_l, g_theta_pi, g_phi_pi,
526  interaction, interaction->KinePtr());
527 
528  // computing cross section for the current kinematics
529  xsec = fXSecModel->XSec(interaction,kPSElOlOpifE) / (1E-38 * units::cm2);
530 
531  if (!fGenerateUniformly) {
532  //-- decide whether to accept the current kinematics
533  double t = xsec_max * rnd->RndKine().Rndm();
534 
535  LOG("COHKinematics", pINFO) << "Got: xsec = " << xsec << ", t = " <<
536  t << " (max_xsec = " << xsec_max << ")";
537 
538  this->AssertXSecLimits(interaction, xsec, xsec_max);
539 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
540  LOG("COHKinematics", pDEBUG)
541  << "xsec= " << xsec << ", J= 1, Rnd= " << t;
542 #endif
543  accept = (t<xsec);
544  }
545  else {
546  accept = (xsec>0);
547  }
548 
549  //-- If the generated kinematics are accepted, finish-up module's job
550  if(accept) {
551  LOG("COHKinematics", pNOTICE) << "Selected: Lepton(" <<
552  g_E_l << ", " << g_theta_l << ", " <<
553  g_phi_l << ") Pion(" << g_theta_pi << ", " << g_phi_pi << ")";
554 
555  double E_l = g_E_l;
556  double theta_l = g_theta_l;
557  double theta_pi = g_theta_pi;
558  double phi_l = g_phi_l;
559  double phi_pi = g_phi_pi;
560  const TLorentzVector P4_nu = *(interaction->InitStatePtr()->GetProbeP4(kRfLab));
561  double E_nu = P4_nu.E();
562  double E_pi= E_nu-E_l;
563  double m_l = interaction->FSPrimLepton()->Mass();
564  double m_pi = this->pionMass(interaction);
565 
566  double p_l = TMath::Sqrt(E_l*E_l - m_l*m_l);
567  TVector3 lepton_3vector = TVector3(0,0,0);
568  lepton_3vector.SetMagThetaPhi(p_l,theta_l,phi_l);
569  TLorentzVector P4_lep = TLorentzVector(lepton_3vector , E_l );
570 
571  double p_pi = TMath::Sqrt(E_pi*E_pi - m_pi*m_pi);
572  TVector3 pion_3vector = TVector3(0,0,0);
573  pion_3vector.SetMagThetaPhi(p_pi,theta_pi,phi_pi);
574  TLorentzVector P4_pion = TLorentzVector(pion_3vector , E_pi);
575 
576  TLorentzVector q = P4_nu - P4_lep;
577  double Q2 = -q.Mag2();
578  double x = Q2/(2*E_pi*constants::kNucleonMass);
579  double y = E_pi/E_nu;
580 
581  double t = TMath::Abs( (q - P4_pion).Mag2() );
582 
583  // for uniform kinematics, compute an event weight as
584  // wght = (phase space volume)*(differential xsec)/(event total xsec)
585  if(fGenerateUniformly) {
586  // Phase space volume needs checking
587  double vol = d_E_l*d_ctheta_l*d_phi*d_ctheta_pi*d_phi;
588  double totxsec = evrec->XSec();
589  double wght = (vol/totxsec)*xsec;
590  LOG("COHKinematics", pNOTICE) << "Kinematics wght = "<< wght;
591 
592  // apply computed weight to the current event weight
593  wght *= evrec->Weight();
594  LOG("COHKinematics", pNOTICE) << "Current event wght = " << wght;
595  evrec->SetWeight(wght);
596  }
597 
598  // reset bits
599  interaction->ResetBit(kISkipProcessChk);
600  interaction->ResetBit(kISkipKinematicChk);
601  // lock selected kinematics & clear running values
602  interaction->KinePtr()->Setx(x, true);
603  interaction->KinePtr()->Sety(y, true);
604  interaction->KinePtr()->Sett(t, true);
605  interaction->KinePtr()->SetW(kPionMass, true);
606  interaction->KinePtr()->SetQ2(2*kNucleonMass*x*y*E_nu, true);
607  interaction->KinePtr()->ClearRunningValues();
608  // set the cross section for the selected kinematics
609  evrec->SetDiffXSec(xsec,kPSElOlOpifE);
610  return;
611  }
612  }//while
613 }
614 //___________________________________________________________________________
616  const double theta_l,
617  const double phi_l,
618  const double theta_pi,
619  const double phi_pi,
620  const Interaction* interaction,
621  Kinematics* kinematics) const
622 {
623  const TLorentzVector P4_nu = *(interaction->InitStatePtr()->GetProbeP4(kRfLab));
624  double E_nu = P4_nu.E();
625  double E_pi= E_nu-E_l;
626  double m_l = interaction->FSPrimLepton()->Mass();
627  double m_pi;
628  if ( interaction->ProcInfo().IsWeakCC() ) {
629  m_pi = constants::kPionMass;
630  } else {
631  m_pi = constants::kPi0Mass;
632  }
633  double p_l=0.0;
634  if (E_l > m_l) {
635  p_l = TMath::Sqrt(E_l*E_l - m_l*m_l);
636  }
637  TVector3 lepton_3vector = TVector3(0,0,0);
638  lepton_3vector.SetMagThetaPhi(p_l,theta_l,phi_l);
639  TLorentzVector P4_lep = TLorentzVector(lepton_3vector , E_l );
640 
641  double p_pi=0.0;
642  if (E_pi > m_pi) {
643  p_pi = TMath::Sqrt(E_pi*E_pi - m_pi*m_pi);
644  }
645  TVector3 pion_3vector = TVector3(0,0,0);
646  pion_3vector.SetMagThetaPhi(p_pi,theta_pi,phi_pi);
647  TLorentzVector P4_pion = TLorentzVector(pion_3vector , E_pi);
648 
649  double Q2 = -(P4_nu-P4_lep).Mag2();
650  double x = Q2/(2*E_pi*constants::kNucleonMass);
651  double y = E_pi/E_nu;
652 
653  kinematics->Setx(x);
654  kinematics->Sety(y);
655  kinematics::UpdateWQ2FromXY(interaction);
656 
657  kinematics->SetFSLeptonP4(P4_lep );
658  kinematics->SetHadSystP4 (P4_pion); // use Hadronic System variable to store pion momentum
659 }
660 //___________________________________________________________________________
662  const double /* theta_l */ ,
663  const double /* phi_l */ ,
664  const double /* theta_pi */ ,
665  const double /* phi_pi */ ,
666  const Interaction* interaction) const
667 {
668  const TLorentzVector P4_nu = *(interaction->InitStatePtr()->GetProbeP4(kRfLab));
669  double E_nu = P4_nu.E();
670  double E_pi= E_nu-E_l;
671  double m_l = interaction->FSPrimLepton()->Mass();
672  double m_pi;
673  if ( interaction->ProcInfo().IsWeakCC() ) {
674  m_pi = constants::kPionMass;
675  }
676  else {
677  m_pi = constants::kPi0Mass;
678  }
679  if (E_l <= m_l) {
680  return false;
681  }
682  if (E_pi <= m_pi) {
683  return false;
684  }
685  return true;
686 }
687 //___________________________________________________________________________
689 {
690  // Computes the maximum differential cross section in the requested phase
691  // space. This method overloads KineGeneratorWithCache::ComputeMaxXSec
692  // method and the value is cached at a circular cache branch for retrieval
693  // during subsequent event generation.
694 
695 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
696  SLOG("COHKinematics", pDEBUG)
697  << "Scanning the allowed phase space {K} for the max(dxsec/d{K})";
698 #endif
699  double max_xsec = 0.;
700  if (fXSecModel->Id().Name() == "genie::ReinSehgalCOHPiPXSec") {
701  max_xsec = MaxXSec_ReinSehgal(in);
702  } else if ((fXSecModel->Id().Name() == "genie::BergerSehgalCOHPiPXSec2015")) {
703  max_xsec = MaxXSec_BergerSehgal(in);
704  } else if ((fXSecModel->Id().Name() == "genie::BergerSehgalFMCOHPiPXSec2015")) {
705  max_xsec = MaxXSec_BergerSehgalFM(in);
706  } else if ((fXSecModel->Id().Name() == "genie::AlvarezRusoCOHPiPXSec")) {
707  max_xsec = MaxXSec_AlvarezRuso(in);
708  }
709  else {
710  LOG("COHKinematicsGenerator",pFATAL) <<
711  "ComputeMaxXSec >> Cannot calculate max cross-section for " <<
712  fXSecModel->Id().Name();
713  }
714 
715  // Apply safety factor, since value retrieved from the cache might
716  // correspond to a slightly different energy.
717  max_xsec *= fSafetyFactor;
718 
719 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
720  SLOG("COHKinematics", pDEBUG) << in->AsString();
721  SLOG("COHKinematics", pDEBUG) << "Max xsec in phase space = " << max_xsec;
722  SLOG("COHKinematics", pDEBUG) << "Computed using alg = " << fXSecModel->Id();
723 #endif
724 
725  return max_xsec;
726 }
727 //___________________________________________________________________________
729 {
730  double max_xsec = 0;
731  const int NQ2 = 50;
732  const int Ny = 50;
733 
734  const KPhaseSpace & kps = in->PhaseSpace();
735  Range1D_t Q2r = kps.Q2Lim();
736  Q2r.max = fQ2Max;
737 
738  const double logQ2min = TMath::Log10(Q2r.min + kASmallNum);
739  const double logQ2max = TMath::Log10(Q2r.max);
740  const double dlogQ2 = (logQ2max - logQ2min) /(NQ2-1);
741 
742  for(int i=0; i<NQ2; i++) {
743  double Q2 = TMath::Power(10, logQ2min + i * dlogQ2);
744  in->KinePtr()->SetQ2(Q2);
745 
746  Range1D_t yr = kps.YLim();
747  if ((yr.max < 0) || (yr.max < yr.min) ||
748  (yr.max > 1) || (yr.min < 0)) { // forbidden kinematics
749  continue;
750  }
751  const double logymin = TMath::Log10(yr.min);
752  const double logymax = TMath::Log10(yr.max);
753  const double dlogy = (logymax - logymin) /(Ny-1);
754 
755  for(int j=0; j<Ny; j++) {
756  double gy = TMath::Power(10, logymin + j * dlogy);
757  in->KinePtr()->Sety(gy);
758 
759  /* Range1D_t tl = kps.TLim(); // TESTING! - this becomes a loop over t */
761 
762  // Note: We're not stepping through log Q^2, log y - we "unpacked"
763  double xsec = fXSecModel->XSec(in, kPSQ2yfE);
764 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
765  LOG("COHKinematics", pDEBUG)
766  << "xsec(Q2= " << Q2 << ", y= " << gy << ", t = " << gt << ") = " << xsec;
767 #endif
768  max_xsec = TMath::Max(max_xsec, xsec);
769 
770  } // y
771  } // Q2
772  return max_xsec;
773 }
774 //___________________________________________________________________________
776 {
777  double max_xsec = 0;
778  // How many sampling bins in each variable for max xsec calculation?
779  const int NQ2 = 50;
780  const int Ny = 50;
781  const int Nt = 50;
782 
783  const KPhaseSpace & kps = in->PhaseSpace();
784  Range1D_t Q2r = kps.Q2Lim();
785  Q2r.max = fQ2Max;
786 
787  const double logQ2min = TMath::Log10(Q2r.min + kASmallNum);
788  const double logQ2max = TMath::Log10(Q2r.max);
789  const double logtmin = TMath::Log10(kASmallNum);
790  const double logtmax = TMath::Log10(fTMax - kASmallNum);
791  const double dlogQ2 = (logQ2max - logQ2min) /(NQ2-1);
792  const double dlogt = (logtmax - logtmin) /(Nt-1);
793 
794  for(int i=0; i<NQ2; i++) {
795  double Q2 = TMath::Power(10, logQ2min + i * dlogQ2);
796  in->KinePtr()->SetQ2(Q2);
797 
798  Range1D_t yr = kps.YLim();
799  if ((yr.max < 0) || (yr.max < yr.min) ||
800  (yr.max > 1) || (yr.min < 0)) { // forbidden kinematics
801  continue;
802  }
803  const double logymin = TMath::Log10(yr.min);
804  const double logymax = TMath::Log10(yr.max);
805  const double dlogy = (logymax - logymin) /(Ny-1);
806 
807  for(int j=0; j<Ny; j++) {
808  double gy = TMath::Power(10, logymin + j * dlogy);
809 
810  for(int k=0; k<Nt; k++) {
811  double gt = TMath::Power(10, logtmin + k * dlogt);
812 
813  in->KinePtr()->Sety(gy);
814  in->KinePtr()->Sett(gt);
815 
816  double xsec = fXSecModel->XSec(in, kPSxyfE);
817 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
818  LOG("COHKinematics", pDEBUG)
819  << "xsec(Q2= " << Q2 << ", y= " << gy << ", t = " << gt << ") = " << xsec;
820 #endif
821  max_xsec = TMath::Max(max_xsec, xsec);
822 
823  } // t
824  } // y
825  } // Q2
826  return max_xsec;
827 }
828 //___________________________________________________________________________
830 {
831  double max_xsec = 0;
832  double Ev = in->InitState().ProbeE(kRfLab);
833 
834  const int Nx = 50;
835  const int Ny = 50;
836 
837  const KPhaseSpace & kps = in->PhaseSpace();
838  Range1D_t y = kps.YLim();
839 
840  const double logxmin = TMath::Log10(1E-5);
841  const double logxmax = TMath::Log10(1.0);
842  const double logymin = TMath::Log10(y.min);
843  const double logymax = TMath::Log10(y.max);
844 
845  const double dlogx = (logxmax - logxmin) /(Nx-1);
846  const double dlogy = (logymax - logymin) /(Ny-1);
847 
848  for(int i=0; i<Nx; i++) {
849  double gx = TMath::Power(10, logxmin + i * dlogx);
850  for(int j=0; j<Ny; j++) {
851  double gy = TMath::Power(10, logymin + j * dlogy);
852 
853  double Q2 = 2*kNucleonMass*gx*gy*Ev;
854  if(Ev>1.0 && Q2>0.01) continue;
855 
856  in->KinePtr()->Setx(gx);
857  in->KinePtr()->Sety(gy);
858 
859  double xsec = fXSecModel->XSec(in, kPSxyfE);
860 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
861  LOG("COHKinematics", pDEBUG)
862  << "xsec(x= " << gx << ", y= " << gy << ") = " << xsec;
863 #endif
864  max_xsec = TMath::Max(max_xsec, xsec);
865 
866  }//y
867  }//x
868  return max_xsec;
869 }
870 //___________________________________________________________________________
872 {
873  // Computes the maximum differential cross section in the requested phase
874  // space. This method overloads KineGeneratorWithCache::ComputeMaxXSec
875  // method and the value is cached at a circular cache branch for retrieval
876  // during subsequent event generation.
877 
878 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
879  SLOG("COHKinematics", pDEBUG)
880  << "Scanning the allowed phase space {K} for the max(dxsec/d{K})";
881 #endif
882  double max_xsec = 0.;
883  double Ev = in->InitState().ProbeE(kRfLab);
884 
885  const KPhaseSpace & kps = in->PhaseSpace();
886  Range1D_t y = kps.YLim();
887 
888  ROOT::Math::Minimizer * min = ROOT::Math::Factory::CreateMinimizer("Minuit2");
890  f.SetFactor(-1.); // Make it return negative of cross-section so we can minimize
891 
892  min->SetFunction( f );
893  min->SetMaxFunctionCalls(10000);
894  min->SetTolerance(0.05);
895 
896  const double min_el = in->FSPrimLepton()->Mass();
897  const double max_el = Ev - kPionMass;
898  const unsigned int n_el = 100;
899  const double d_el = (max_el - min_el) / double(n_el - 1);
900 
901  const double min_thetal = kASmallNum;
902  const double max_thetal = kPi / 4.0;
903  const unsigned int n_thetal = 10;
904  const double d_thetal = (max_thetal - min_thetal) / double(n_thetal - 1);
905 
906  const double min_thetapi = kASmallNum;
907  const double max_thetapi = kPi / 2.0;
908  const unsigned int n_thetapi = 10;
909  const double d_thetapi = (max_thetapi - min_thetapi) / double(n_thetapi - 1);
910 
911  //~ const double min_phipi = kPi;
912  //~ const double max_phipi = 0.5 * kPi;
913  const double min_phipi = kASmallNum;
914  const double max_phipi = 2*kPi-kASmallNum;
915  const unsigned int n_phipi = 10;
916  const double d_phipi = (max_phipi - min_phipi) / double(n_phipi - 1);
917 
918  min->SetLimitedVariable ( 0 ,"E_lep" , max_el -kASmallNum , d_el , min_el , max_el );
919  min->SetLimitedVariable ( 1 ,"theta_l" , min_thetal +kASmallNum , d_thetal , min_thetal , max_thetal );
920  min->SetLimitedVariable ( 2 ,"theta_pi" , min_thetapi+kASmallNum , d_thetapi , min_thetapi, max_thetapi );
921  min->SetLimitedVariable ( 3 ,"phi_pi" , min_phipi +kASmallNum , d_phipi , min_phipi , max_phipi );
922 
923  min->Minimize();
924  max_xsec = -min->MinValue(); //back to positive xsec
925 
926  // Apply safety factor, since value retrieved from the cache might
927  // correspond to a slightly different energy.
928  max_xsec *= fSafetyFactor;
929 
930 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
931  SLOG("COHKinematics", pDEBUG) << in->AsString();
932  SLOG("COHKinematics", pDEBUG) << "Max xsec in phase space = " << max_xsec;
933  SLOG("COHKinematics", pDEBUG) << "Computed using alg = " << fXSecModel->Id();
934 #endif
935 
936  delete min;
937 
938  return max_xsec;
939 }
940 //___________________________________________________________________________
941 double COHKinematicsGenerator::Energy(const Interaction * interaction) const
942 {
943  // Override the base class Energy() method to cache the max xsec for the
944  // neutrino energy in the LAB rather than in the hit nucleon rest frame.
945 
946  const InitialState & init_state = interaction->InitState();
947  double E = init_state.ProbeE(kRfLab);
948  return E;
949 }
950 //___________________________________________________________________________
952 {
953  double m_pi = 0.0;
954  if ( in->ProcInfo().IsWeakCC() ) {
955  m_pi = constants::kPionMass;
956  } else {
957  m_pi = constants::kPi0Mass;
958  }
959  return m_pi;
960 }
961 //___________________________________________________________________________
963  GHepRecord* evrec) const
964 {
965  LOG("COHKinematics", pWARN)
966  << "*** Could not select valid kinematics after "
967  << iters << " iterations";
968  evrec->EventFlags()->SetBitNumber(kKineGenErr, true);
970  exception.SetReason("Couldn't select kinematics");
971  exception.SwitchOnFastForward();
972  throw exception;
973 }
974 //___________________________________________________________________________
976 {
977  Algorithm::Configure(config);
978  this->LoadConfig();
979 }
980 //____________________________________________________________________________
982 {
983  Algorithm::Configure(config);
984  this->LoadConfig();
985 }
986 //____________________________________________________________________________
988 {
989  //-- COH model parameter Ro
990  GetParam( "COH-Ro", fRo );
991  //-- COH model parameter t_max for t = (q - p_pi)^2
992  GetParam( "COH-t-max", fTMax ) ;
993  //-- COH model bounds of integration for Q^2
994  GetParam( "COH-Q2-min", fQ2Min ) ;
995  GetParam( "COH-Q2-max", fQ2Max ) ;
996 
997  //-- max xsec safety factor (for rejection method) and min cached energy
998  GetParamDef( "MaxXSec-SafetyFactor", fSafetyFactor, 1.6 ) ;
999  GetParamDef( "Cache-MinEnergy", fEMin, -1.0 ) ;
1000 
1001  //-- Generate kinematics uniformly over allowed phase space and compute
1002  // an event weight?
1003  GetParamDef( "UniformOverPhaseSpace", fGenerateUniformly, false ) ;
1004 
1005  //-- Maximum allowed fractional cross section deviation from maxim cross
1006  // section used in rejection method
1007  GetParamDef( "MaxXSec-DiffTolerance", fMaxXSecDiffTolerance, 999999. ) ;
1009 
1010  //-- Envelope employed when importance sampling is used
1011  // (initialize with dummy range)
1012  if(fEnvelope) delete fEnvelope;
1013  fEnvelope = new TF2("CohKinEnvelope",
1015  // stop ROOT from deleting this object of its own volition
1016  gROOT->GetListOfFunctions()->Remove(fEnvelope);
1017 }
1018 //____________________________________________________________________________
1019 
double COHImportanceSamplingEnvelope(double *x, double *par)
Definition: KineUtils.cxx:1343
const double kPi
virtual double MaxXSec(GHepRecord *evrec) const
const KPhaseSpace & PhaseSpace(void) const
Definition: Interaction.h:73
Basic constants.
virtual void SetWeight(double wght)
Definition: GHepRecord.h:131
bool IsWeakCC(void) const
std::map< std::string, double > xmax
bool fGenerateUniformly
uniform over allowed phase space + event weight?
double MaxXSec_ReinSehgal(const Interaction *in) const
THE MAIN GENIE PROJECT NAMESPACE
Definition: GeneratorBase.h:8
double fRo
nuclear scale parameter
static const double kNucleonMass
Definition: Constants.h:78
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:991
virtual Interaction * Summary(void) const
Definition: GHepRecord.cxx:135
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:79
void SetQ2(double Q2, bool selected=false)
Definition: Kinematics.cxx:265
Kinematics * KinePtr(void) const
Definition: Interaction.h:76
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.
static const double fermi
Definition: Units.h:63
int A(void) const
Definition: Target.h:71
A simple [min,max] interval for doubles.
Definition: Range1.h:43
virtual void AssertXSecLimits(const Interaction *in, double xsec, double xsec_max) const
static const double kPi0Mass
Definition: Constants.h:75
#define pFATAL
Definition: Messenger.h:57
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.
Definition: Kinematics.h:40
void Configure(const Registry &config)
virtual double Weight(void) const
Definition: GHepRecord.h:125
double x(bool selected=false) const
Definition: Kinematics.cxx:109
::xsd::cxx::tree::exception< char > exception
Definition: Database.h:225
Range1D_t YLim(void) const
y limits
Definition: config.py:1
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.
static const double cm2
Definition: Units.h:77
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...
Definition: RandomGen.h:30
double MaxXSec_AlvarezRuso(const Interaction *in) const
string AsString(void) const
void CalculateKin_ReinSehgal(GHepRecord *event_rec) const
Double_t ymax
Definition: plot.C:25
void ProcessEventRecord(GHepRecord *event_rec) const
Summary information for an interaction.
Definition: Interaction.h:56
Definition: Cand.cxx:23
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...
double dy[NP][NC]
double dx[NP][NC]
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:97
void SetFSLeptonP4(const TLorentzVector &p4)
Definition: Kinematics.cxx:307
string Name(void) const
Definition: AlgId.h:45
Float_t E
Definition: plot.C:20
double fQ2Min
lower bound of integration for Q^2 in Berger-Sehgal Model
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:70
Kinematical phase space.
Definition: KPhaseSpace.h:34
void Sett(double t, bool selected=false)
Definition: Kinematics.cxx:301
#define R(x)
const double j
Definition: BetheBloch.cxx:29
static const double kASmallNum
Definition: Controls.h:40
void CalculateKin_BergerSehgal(GHepRecord *event_rec) const
#define pINFO
Definition: Messenger.h:63
TRandom3 & RndKine(void) const
rnd number generator used by kinematics generators
Definition: RandomGen.h:51
Misc GENIE control constants.
Double_t xsec[nknots]
Definition: testXsec.C:47
#define pWARN
Definition: Messenger.h:61
const UInt_t kISkipKinematicChk
if set, skip kinematic validity checks
Definition: Interaction.h:48
virtual const XSecAlgorithmI * CrossSectionAlg(void) const =0
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
static float min(const float a, const float b, const float c)
Definition: absgeo.cxx:45
void Setx(double x, bool selected=false)
Definition: Kinematics.cxx:241
float Mag2() const
void UpdateWQ2FromXY(const Interaction *in)
Definition: KineUtils.cxx:1185
static RunningThreadInfo * Instance(void)
double ComputeMaxXSec(const Interaction *in) const
static const double A
Definition: Units.h:82
virtual const AlgId & Id(void) const
Get algorithm ID.
Definition: Algorithm.h:98
double max
Definition: Range1.h:54
void SetW(double W, bool selected=false)
Definition: Kinematics.cxx:289
ifstream in
Definition: comparison.C:7
double Energy(const Interaction *in) const
virtual TBits * EventFlags(void) const
Definition: GHepRecord.h:118
static const double kPionMass
Definition: Constants.h:74
Var Sqrt(const Var &v)
Use to take sqrt of a var.
Definition: Var.cxx:326
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:66
void Sety(double y, bool selected=false)
Definition: Kinematics.cxx:253
double MaxXSec_BergerSehgal(const Interaction *in) const
void UpdateXFromQ2Y(const Interaction *in)
Definition: KineUtils.cxx:1221
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
const hit & b
Definition: hits.cxx:21
assert(nhit_max >=nhit_nbins)
void SetHadSystP4(const TLorentzVector &p4)
Definition: Kinematics.cxx:317
virtual double XSec(void) const
Definition: GHepRecord.h:127
double gQ2
Definition: gtestDISSF.cxx:55
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
TF2 * fEnvelope
2-D envelope used for importance sampling
const InitialState & InitState(void) const
Definition: Interaction.h:69
Double_t ymin
Definition: plot.C:24
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
double min
Definition: Range1.h:53
#define pNOTICE
Definition: Messenger.h:62
void ClearRunningValues(void)
Definition: Kinematics.cxx:357
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
Definition: InitialState.h:67
static const unsigned int kRjMaxIterations
Definition: Controls.h:27
bool gt(unsigned short int a, unsigned short int b)
T max(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:68
const EventGeneratorI * RunningThread(void)
#define SLOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a short string (using the FUNCTION and...
Definition: Messenger.h:85
void kinematics()
Definition: kinematics.C:10
double ProbeE(RefFrame_t rf) const
GENIE&#39;s GHEP MC event record.
Definition: GHepRecord.h:46
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.
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
Definition: Interaction.h:47
virtual void SetDiffXSec(double xsec, KinePhaseSpace_t ps)
Definition: GHepRecord.h:134
void throwOnTooManyIterations(unsigned int iters, GHepRecord *evrec) const
Initial State information.
Definition: InitialState.h:49
#define pDEBUG
Definition: Messenger.h:64
static const double kPionMass2
Definition: Constants.h:87