HNIntranuke2014.cxx
Go to the documentation of this file.
1 
2 //____________________________________________________________________________
3 /*
4  Copyright (c) 2003-2017, GENIE Neutrino MC Generator Collaboration
5  For the full text of the license visit http://copyright.genie-mc.org
6  or see $GENIE/LICENSE
7 
8  Author: Steve Dytman <dytman+@pitt.edu>, Pittsburgh Univ.
9  Aaron Meyer <asm58@pitt.edu>, Pittsburgh Univ.
10  Alex Bell, Pittsburgh Univ.
11  Hugh Gallagher <gallag@minos.phy.tufts.edu>, Tufts Univ.
12  Costas Andreopoulos <costas.andreopoulos \at stfc.ac.uk>, Rutherford Lab.
13  September 20, 2005
14 
15  For the class documentation see the corresponding header file.
16 
17  Important revisions after version 2.0.0 :
18  @ Nov 30, 2007 - SD
19  Changed the hadron tracking algorithm to take into account the radial
20  nuclear density dependence. Using the somewhat empirical approach of
21  increasing the nuclear radius by a const (tunable) number times the tracked
22  particle's de Broglie wavelength as this helps getting the hadron+nucleus
23  cross sections right.
24  @ Mar 08, 2008 - CA
25  Fixed code retrieving the remnant nucleus which stopped working as soon as
26  simulation of nuclear de-excitation started pushing photons in the target
27  nucleus daughter list.
28  @ Jun 20, 2008 - CA
29  Fix a mem leak: The (clone of the) GHepParticle being re-scattered was not
30  deleted after it was added at the GHEP event record.
31  @ Jul 15, 2010 - AM
32  The hN mode is now implemented in Intranuke. Similar to hA mode, but particles
33  produced by reactions are stepped through the nucleus like probe particles.
34  Particles react with nucleons instead of the entire nucleus, and final states
35  are determined after reactions are finished, not before.
36  @ Dec 15, 2014 - SD, NG
37  Update fates to include Compound Nucleus final state correctly.
38  @ Jan 9, 2015 - SD, NG, TG
39  Added 2014 version of INTRANUKE codes (new class) for independent development.
40 */
41 //____________________________________________________________________________
42 
43 #include <cstdlib>
44 #include <sstream>
45 
46 #include <TMath.h>
47 
49 #include "Algorithm/AlgFactory.h"
50 #include "Conventions/GBuild.h"
51 #include "Conventions/Constants.h"
52 #include "Conventions/Controls.h"
54 #include "GHEP/GHepFlags.h"
55 #include "GHEP/GHepStatus.h"
56 #include "GHEP/GHepRecord.h"
57 #include "GHEP/GHepParticle.h"
64 #include "Messenger/Messenger.h"
65 #include "Numerical/RandomGen.h"
66 #include "Numerical/Spline.h"
67 #include "PDG/PDGLibrary.h"
68 #include "PDG/PDGCodes.h"
69 #include "PDG/PDGCodeList.h"
70 #include "PDG/PDGUtils.h"
71 #include "Utils/PrintUtils.h"
72 #include "Utils/NuclearUtils.h"
73 
74 using std::ostringstream;
75 
76 using namespace genie;
77 using namespace genie::utils;
78 using namespace genie::utils::intranuke2014;
79 using namespace genie::constants;
80 using namespace genie::controls;
81 
82 //___________________________________________________________________________
83 //___________________________________________________________________________
84 // Methods specific to INTRANUKE's HN-mode
85 //___________________________________________________________________________
86 //___________________________________________________________________________
88 Intranuke2014("genie::HNIntranuke2014")
89 {
90 
91 }
92 //___________________________________________________________________________
94 Intranuke2014("genie::HNIntranuke2014",config)
95 {
96 
97 }
98 //___________________________________________________________________________
100 {
101 
102 }
103 //___________________________________________________________________________
105 {
106  LOG("HNIntranuke2014", pNOTICE)
107  << "************ Running hN2014 MODE INTRANUKE ************";
108 
109  LOG("HNIntranuke2014", pWARN)
111  "Experimental code (INTRANUKE/hN model) - Run at your own risk");
112 
114 
115  LOG("HNIntranuke2014", pINFO) << "Done with this event";
116 }
117 //___________________________________________________________________________
119 {
120 // Simulate a hadron interaction for the input particle p in HN mode
121 //
122  if(!p || !ev)
123  {
124  LOG("HNIntranuke2014", pERROR) << "** Null input!";
125  return;
126  }
127 
128  // check particle id
129  int pdgc = p->Pdg();
130  bool is_pion = (pdgc==kPdgPiP || pdgc==kPdgPiM || pdgc==kPdgPi0);
131  // bool is_kaon = (pdgc==kPdgKP); // UNUSED - comment to quiet compiler warnings
132  bool is_baryon = (pdgc==kPdgProton || pdgc==kPdgNeutron);
133  bool is_gamma = (pdgc==kPdgGamma);
134  if(!(is_pion || is_baryon || is_gamma))
135  {
136  LOG("HNIntranuke2014", pERROR) << "** Cannot handle particle: " << p->Name();
137  return;
138  }
139  try
140  {
141  // select a fate for the input particle
142  INukeFateHN_t fate = this->HadronFateHN(p);
143 
144  // store the fate
145  ev->Particle(p->FirstMother())->SetRescatterCode((int)fate);
146 
147  if(fate == kIHNFtUndefined)
148  {
149  LOG("HNIntranuke2014", pERROR) << "** Couldn't select a fate";
150  LOG("HNIntranuke2014", pERROR) << "** Num Protons: " << fRemnZ
151  << ", Num Neutrons: "<<(fRemnA-fRemnZ);
152  LOG("HNIntranuke2014", pERROR) << "** Particle: " << "\n" << (*p);
153  //LOG("HNIntranuke2014", pERROR) << "** Event Record: " << "\n" << (*ev);
154  //p->SetStatus(kIStUndefined);
155  p->SetStatus(kIStStableFinalState);
156  ev->AddParticle(*p);
157  return;
158  }
159 
160  LOG("HNIntranuke2014", pNOTICE)
161  << "Selected " << p->Name() << " fate: " << INukeHadroFates::AsString(fate);
162 
163  // handle the reaction
164  if(fate == kIHNFtCEx || fate == kIHNFtElas)
165  {
166  this->ElasHN(ev,p,fate);
167  }
168  else if(fate == kIHNFtAbs) {this-> AbsorbHN(ev,p,fate);}
169  else if(fate == kIHNFtInelas && pdgc != kPdgGamma)
170  {
171 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
172  LOG("HNIntranuke2014", pDEBUG)
173  << "Invoking InelasticHN() for a : " << p->Name()
174  << " whose fate is : " << INukeHadroFates::AsString(fate);
175 #endif
176 
177  this-> InelasticHN(ev,p);
178  }
179  else if(fate == kIHNFtInelas && pdgc == kPdgGamma) {this-> GammaInelasticHN(ev,p,fate);}
180  else if(fate == kIHNFtNoInteraction)
181  {
183  ev->AddParticle(*p);
184  return;
185  }
187 
188 
189  }
191  {
192  this->SimulateHadronicFinalState(ev,p);
193  LOG("HNIntranuke2014", pNOTICE)
194  << "retry call to SimulateHadronicFinalState ";
195  LOG("HNIntranuke2014", pNOTICE) << exception;
196 
197  }
198 }
199 //___________________________________________________________________________
201 {
202 // Select a hadron fate in HN mode
203 //
204  RandomGen * rnd = RandomGen::Instance();
205 
206  // get pdgc code & kinetic energy in MeV
207  int pdgc = p->Pdg();
208  double ke = p->KinE() / units::MeV;
209 
210  LOG("HNIntranuke2014", pINFO)
211  << "Selecting hN fate for " << p->Name() << " with KE = " << ke << " MeV";
212 
213  // try to generate a hadron fate
214  unsigned int iter = 0;
215  while(iter++ < kRjMaxIterations) {
216 
217  // handle pions
218  //
219  if (pdgc==kPdgPiP || pdgc==kPdgPiM || pdgc==kPdgPi0) {
220 
221  double frac_cex = this->FateWeight(pdgc, kIHNFtCEx)
222  * fHadroData2014->Frac(pdgc, kIHNFtCEx, ke, fRemnA, fRemnZ);
223  double frac_elas = this->FateWeight(pdgc, kIHNFtElas)
224  * fHadroData2014->Frac(pdgc, kIHNFtElas, ke, fRemnA, fRemnZ);
225  double frac_inel = this->FateWeight(pdgc, kIHNFtInelas)
226  * fHadroData2014->Frac(pdgc, kIHNFtInelas, ke, fRemnA, fRemnZ);
227  double frac_abs = this->FateWeight(pdgc, kIHNFtAbs)
228  * fHadroData2014->Frac(pdgc, kIHNFtAbs, ke, fRemnA, fRemnZ);
229  frac_cex *= fNucCEXFac; // scaling factors
230  frac_abs *= fNucAbsFac;
231  frac_elas *= fNucQEFac;
232  if(pdgc==kPdgPi0) frac_abs*= 0.665; //isospin factor
233 
234  LOG("HNIntranuke2014", pINFO)
235  << "\n frac{" << INukeHadroFates::AsString(kIHNFtCEx) << "} = " << frac_cex
236  << "\n frac{" << INukeHadroFates::AsString(kIHNFtElas) << "} = " << frac_elas
237  << "\n frac{" << INukeHadroFates::AsString(kIHNFtInelas) << "} = " << frac_inel
238  << "\n frac{" << INukeHadroFates::AsString(kIHNFtAbs) << "} = " << frac_abs;
239 
240  // compute total fraction (can be <1 if fates have been switched off)
241  double tf = frac_cex +
242  frac_elas +
243  frac_inel +
244  frac_abs;
245 
246  double r = tf * rnd->RndFsi().Rndm();
247 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
248  LOG("HNIntranuke2014", pDEBUG) << "r = " << r << " (max = " << tf << ")";
249 #endif
250 
251  double cf=0; // current fraction
252 
253  if(r < (cf += frac_cex )) return kIHNFtCEx; //cex
254  if(r < (cf += frac_elas )) return kIHNFtElas; //elas
255  if(r < (cf += frac_inel )) return kIHNFtInelas; //inelas
256  if(r < (cf += frac_abs )) return kIHNFtAbs; //abs
257 
258  LOG("HNIntranuke2014", pWARN)
259  << "No selection after going through all fates! "
260  << "Total fraction = " << tf << " (r = " << r << ")";
261  ////////////////////////////
262  return kIHNFtUndefined;
263  }
264 
265  // handle nucleons
266  else if (pdgc==kPdgProton || pdgc==kPdgNeutron) {
267 
268  double frac_elas = this->FateWeight(pdgc, kIHNFtElas)
269  * fHadroData2014->Frac(pdgc, kIHNFtElas, ke, fRemnA, fRemnZ);
270  double frac_inel = this->FateWeight(pdgc, kIHNFtInelas)
271  * fHadroData2014->Frac(pdgc, kIHNFtInelas, ke, fRemnA, fRemnZ);
272  double frac_cmp = this->FateWeight(pdgc, kIHNFtCmp)
273  * fHadroData2014->Frac(pdgc, kIHNFtCmp, ke, fRemnA , fRemnZ);
274 
275  LOG("HNIntranuke2014", pINFO)
276  << "\n frac{" << INukeHadroFates::AsString(kIHNFtElas) << "} = " << frac_elas
277  << "\n frac{" << INukeHadroFates::AsString(kIHNFtInelas) << "} = " << frac_inel;
278 
279  // compute total fraction (can be <1 if fates have been switched off)
280  double tf = frac_elas + frac_inel + frac_cmp;
281 
282  double r = tf * rnd->RndFsi().Rndm();
283 
284 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
285  LOG("HNIntranuke2014", pDEBUG) << "r = " << r << " (max = " << tf << ")";
286 #endif
287 
288  double cf=0; // current fraction
289  if(r < (cf += frac_elas )) return kIHNFtElas; // elas
290  if(r < (cf += frac_inel )) return kIHNFtInelas; // inelas
291  if(r < (cf += frac_cmp )) return kIHNFtCmp; // cmp
292 
293  LOG("HNIntranuke2014", pWARN)
294  << "No selection after going through all fates! "
295  << "Total fraction = " << tf << " (r = " << r << ")";
296  //////////////////////////
297  return kIHNFtUndefined;
298  }
299 
300  // handle gamma -- does not currently consider the elastic case
301  else if (pdgc==kPdgGamma) return kIHNFtInelas;
302 
303  // handle kaon -- currently only elastic case
304  else if (pdgc==kPdgKP) return kIHNFtElas;
305 
306  }//iterations
307 
308  return kIHNFtUndefined;
309 }
310 //___________________________________________________________________________
311 double HNIntranuke2014::FateWeight(int pdgc, INukeFateHN_t fate) const
312 {
313  // turn fates off if the remnant nucleus does not have the number of p,n
314  // required
315 
316  int np = fRemnZ;
317  int nn = fRemnA - fRemnZ;
318 
319  if (np < 1 && nn < 1)
320  {
321  LOG("HNIntranuke2014", pERROR) << "** Nothing left in nucleus!!! **";
322  return 0;
323  }
324  else
325  {
326  if (fate == kIHNFtCEx && pdgc==kPdgPiP ) { return (nn>=1) ? 1. : 0.; }
327  if (fate == kIHNFtCEx && pdgc==kPdgPiM ) { return (np>=1) ? 1. : 0.; }
328  if (fate == kIHNFtAbs) { return ((nn>=1) && (np>=1)) ? 1. : 0.; }
329  if (fate == kIHNFtCmp ) { return ((pdgc==kPdgProton||pdgc==kPdgNeutron)&&fDoCompoundNucleus&&fRemnA>5) ? 1. : 0.; }
330 
331  }
332  return 1.;
333 }
334 //___________________________________________________________________________
336  GHepRecord * ev, GHepParticle * p, INukeFateHN_t fate) const
337 {
338  // handles pi+d->2p, pi-d->nn, pi0 d->pn absorbtion, all using pi+d values
339 
340  int pdgc = p->Pdg();
341 
342 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
343  LOG("HNIntranuke2014", pDEBUG)
344  << "AbsorbHN() is invoked for a : " << p->Name()
345  << " whose fate is : " << INukeHadroFates::AsString(fate);
346 #endif
347 
348  // check fate
349  if(fate!=kIHNFtAbs)
350  {
351  LOG("HNIntranuke2014", pWARN)
352  << "AbsorbHN() cannot handle fate: " << INukeHadroFates::AsString(fate);
353  return;
354  }
355 
356  // random number generator
357  RandomGen * rnd = RandomGen::Instance();
358 
359  // Notes on the kinematics
360  // -- Simple variables are used for efficiency
361  // -- Variables are numbered according to particle
362  // -- -- #1 -> incoming particle
363  // -- -- #2 -> target (here, 2_1 and 2_2 for individual particles)
364  // -- -- #3 -> scattered incoming (Particle tracked in hA mode)
365  // -- -- #4 -> other scattered particle
366  // -- Suffix "L" is for lab frame, suffix "CM" is for center of mass frame
367  // -- Subscript "z" is for parallel component, "t" is for transverse
368 
369  int pcode, t1code, t2code, scode, s2code; // particles
370  double M1, M2_1, M2_2, M3, M4; // rest energies, in GeV
371  double E1L, P1L, E2L, P2L, E3L, P3L, E4L, P4L;
372  double P1zL, P2zL;
373  double beta, gm; // speed and gamma for CM frame in lab
374  double Et, E2CM;
375  double C3CM, S3CM; // cos and sin of scattering angle
376  double Theta1, Theta2, theta5;
377  double PHI3; // transverse scattering angle
378  double E1CM, E3CM, E4CM, P3CM;
379  double P3zL, P3tL, P4zL, P4tL;
380  double E2_1L, E2_2L;
381  TVector3 tP2_1L, tP2_2L, tP1L, tP2L, tPtot, P1zCM, P2zCM;
382  TVector3 tP3L, tP4L;
383  TVector3 bDir, tTrans, tbeta, tVect;
384 
385  // Library instance for reference
386  PDGLibrary * pLib = PDGLibrary::Instance();
387 
388  // Handle fermi target
389  Target target(ev->TargetNucleus()->Pdg());
390 
391  // Target should be a deuteron, but for now
392  // handling it as seperate nucleons
393  if(pdgc==211) // pi-plus
394  {
395  pcode = 211;
396  t1code = 2212; // proton
397  t2code = 2112; // neutron
398  scode = 2212;
399  s2code = 2212;
400  }
401  else if(pdgc==-211) // pi-minus
402  {
403  pcode = -211;
404  t1code = 2212;
405  t2code = 2112;
406  scode = 2112;
407  s2code = 2112;
408  }
409  else if(pdgc==111) // pi-zero
410  {
411  pcode = 111;
412  t1code = 2212;
413  t2code = 2112;
414  scode = 2212;
415  s2code = 2112;
416  }
417  else
418  {
419  LOG("HNIntranuke2014", pWARN)
420  << "AbsorbHN() cannot handle probe: " << pdgc;
421  return;
422  }
423 
424  // assign proper masses
425  M1 = pLib->Find(pcode) ->Mass();
426  M2_1 = pLib->Find(t1code)->Mass();
427  M2_2 = pLib->Find(t2code)->Mass();
428  M3 = pLib->Find(scode) ->Mass();
429  M4 = pLib->Find(s2code)->Mass();
430 
431  // handle fermi momentum
432  if(fDoFermi)
433  {
434  target.SetHitNucPdg(t1code);
436  tP2_1L=fFermiFac * fNuclmodel->Momentum3();
437  E2_1L = TMath::Sqrt(tP2_1L.Mag2() + M2_1*M2_1);
438 
439  target.SetHitNucPdg(t2code);
441  tP2_2L=fFermiFac * fNuclmodel->Momentum3();
442  E2_2L = TMath::Sqrt(tP2_2L.Mag2() + M2_2*M2_2);
443  }
444  else
445  {
446  tP2_1L.SetXYZ(0.0, 0.0, 0.0);
447  E2_1L = M2_1;
448 
449  tP2_2L.SetXYZ(0.0, 0.0, 0.0);
450  E2_2L = M2_2;
451  }
452 
453  E2L = E2_1L + E2_2L;
454 
455  // adjust p to reflect scattering
456  // get random scattering angle
457  C3CM = fHadroData2014->IntBounce(p,t1code,scode,fate);
458  if (C3CM<-1.)
459  {
461  ev->AddParticle(*p);
462  return;
463  }
464  S3CM = TMath::Sqrt(1.0 - C3CM*C3CM);
465 
466  // Get lab energy and momenta
467  E1L = p->E();
468  if(E1L<0.001) E1L=0.001;
469  P1L = TMath::Sqrt(E1L*E1L - M1*M1);
470  tP1L = p->P4()->Vect();
471  tP2L = tP2_1L + tP2_2L;
472  P2L = tP2L.Mag();
473  tPtot = tP1L + tP2L;
474 
475  // get unit vectors and angles needed for later
476  bDir = tPtot.Unit();
477  Theta1 = tP1L.Angle(bDir);
478  Theta2 = tP2L.Angle(bDir);
479 
480  // get parallel and transverse components
481  P1zL = P1L*TMath::Cos(Theta1);
482  P2zL = P2L*TMath::Cos(Theta2);
483  tVect.SetXYZ(1,0,0);
484  if(TMath::Abs((tVect - bDir).Mag())<.01) tVect.SetXYZ(0,1,0);
485  theta5 = tVect.Angle(bDir);
486  tTrans = (tVect - TMath::Cos(theta5)*bDir).Unit();
487 
488  // calculate beta and gamma
489  tbeta = tPtot * (1.0 / (E1L + E2L));
490  beta = tbeta.Mag();
491  gm = 1.0 / TMath::Sqrt(1.0 - beta*beta);
492 
493  // boost to CM frame to get scattered particle momenta
494  E1CM = gm*E1L - gm*beta*P1zL;
495  P1zCM = gm*P1zL*bDir - gm*tbeta*E1L;
496  E2CM = gm*E2L - gm*beta*P2zL;
497  P2zCM = gm*P2zL*bDir - gm*tbeta*E2L;
498  Et = E1CM + E2CM;
499  E3CM = (Et*Et + (M3*M3) - (M4*M4)) / (2.0*Et);
500  E4CM = Et - E3CM;
501  P3CM = TMath::Sqrt(E3CM*E3CM - M3*M3);
502 
503  // boost back to lab
504  P3zL = gm*beta*E3CM + gm*P3CM*C3CM;
505  P3tL = P3CM*S3CM;
506  P4zL = gm*beta*E4CM + gm*P3CM*(-C3CM);
507  P4tL = P3CM*(-S3CM);
508 
509  P3L = TMath::Sqrt(P3zL*P3zL + P3tL*P3tL);
510  P4L = TMath::Sqrt(P4zL*P4zL + P4tL*P4tL);
511 
512  // check for too small values
513  // may introduce error, so warn if it occurs
514  if(!(TMath::Finite(P3L))||P3L<.001)
515  {
516  LOG("HNIntranuke2014",pINFO)
517  << "Particle 3 " << M3 << " momentum small or non-finite: " << P3L
518  << "\n" << "--> Assigning .001 as new momentum";
519  P3tL = 0;
520  P3zL = .001;
521  P3L = .001;
522  E3L = TMath::Sqrt(P3L*P3L + M3*M3);
523  }
524 
525  if(!(TMath::Finite(P4L))||P4L<.001)
526  {
527  LOG("HNIntranuke2014",pINFO)
528  << "Particle 4 " << M4 << " momentum small or non-finite: " << P4L
529  << "\n" << "--> Assigning .001 as new momentum";
530  P4tL = 0;
531  P4zL = .001;
532  P4L = .001;
533  E4L = TMath::Sqrt(P4L*P4L + M4*M4);
534  }
535 
536  // pauli blocking
537  if(P3L < fFermiMomentum || P4L < fFermiMomentum)
538  {
539 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
540  LOG("HNIntranuke2014",pINFO) << "AbsorbHN failed: Pauli blocking";
541 #endif
543 
544  //utils::intranuke2014::StepParticle(p,fFreeStep,fTrackingRadius); //not needed
545 
546  ev->AddParticle(*p);
547  return;
548  }
549 
550  // handle remnant nucleus updates
551  fRemnZ--;
552  fRemnA -=2;
553  fRemnP4 -= TLorentzVector(tP2_1L,E2_1L);
554  fRemnP4 -= TLorentzVector(tP2_2L,E2_2L);
555 
556  // get random phi angle, distributed uniformally in 360 deg
557  PHI3 = 2 * kPi * rnd->RndFsi().Rndm();
558 
559  tP3L = P3zL*bDir + P3tL*tTrans;
560  tP4L = P4zL*bDir + P4tL*tTrans;
561 
562  tP3L.Rotate(PHI3,bDir); // randomize transverse components
563  tP4L.Rotate(PHI3,bDir);
564 
565  E3L = TMath::Sqrt(P3L*P3L + M3*M3);
566  E4L = TMath::Sqrt(P4L*P4L + M4*M4);
567 
568  // create t particle w/ appropriate momenta, code, and status
569  // set target's mom to be the mom of the hadron that was cloned
570  GHepParticle * t = new GHepParticle(*p);
571  t->SetFirstMother(p->FirstMother());
572  t->SetLastMother(p->LastMother());
573 
574  TLorentzVector t4P4L(tP4L,E4L);
575  t->SetPdgCode(s2code);
576  t->SetMomentum(t4P4L);
578 
579  // adjust p to reflect scattering
580  TLorentzVector t4P3L(tP3L,E3L);
581  p->SetPdgCode(scode);
582  p->SetMomentum(t4P3L);
584 
585 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
586  LOG("HNIntranuke2014", pDEBUG)
587  << "|p3| = " << (P3L) << ", E3 = " << (E3L);
588  LOG("HNIntranuke2014", pDEBUG)
589  << "|p4| = " << (P4L) << ", E4 = " << (E4L);
590 #endif
591 
592  ev->AddParticle(*p);
593  ev->AddParticle(*t);
594 
595  delete t; // delete particle clone
596 }
597 //___________________________________________________________________________
599  GHepRecord * ev, GHepParticle * p, INukeFateHN_t fate) const
600 {
601  // scatters particle within nucleus
602 
603 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
604  LOG("HNIntranuke2014", pDEBUG)
605  << "ElasHN() is invoked for a : " << p->Name()
606  << " whose fate is : " << INukeHadroFates::AsString(fate);
607 #endif
608 
609  if(fate!=kIHNFtCEx && fate!=kIHNFtElas)
610  {
611  LOG("HNIntranuke2014", pWARN)
612  << "ElasHN() cannot handle fate: " << INukeHadroFates::AsString(fate);
613  return;
614  }
615 
616  // Random number generator
617  RandomGen * rnd = RandomGen::Instance();
618 
619  // vars for incoming particle, target, and scattered pdg codes
620  int pcode = p->Pdg();
621  int tcode, scode, s2code;
622  double ppcnt = (double) fRemnZ / (double) fRemnA; // % of protons
623 
624  // Select a target randomly, weighted to #
625  // -- Unless, of course, the fate is CEx,
626  // -- in which case the target may be deterministic
627  // Also assign scattered particle code
628  if(fate==kIHNFtCEx)
629  {
630  if(pcode==kPdgPiP) {tcode = kPdgNeutron; scode = kPdgPi0; s2code = kPdgProton;}
631  else if(pcode==kPdgPiM) {tcode = kPdgProton; scode = kPdgPi0; s2code = kPdgNeutron;}
632  else
633  {
634  // for pi0
635  tcode = (rnd->RndFsi().Rndm()<=ppcnt)?(kPdgProton) :(kPdgNeutron);
636  scode = (tcode == kPdgProton) ?(kPdgPiP) :(kPdgPiM);
637  s2code = (tcode == kPdgProton) ?(kPdgNeutron):(kPdgProton);
638  }
639  }
640  else
641  {
642  tcode = (rnd->RndFsi().Rndm()<=ppcnt)?(kPdgProton):(kPdgNeutron);
643  scode = pcode;
644  s2code = tcode;
645  }
646 
647  // get random scattering angle
648  double C3CM = fHadroData2014->IntBounce(p,tcode,scode,fate);
649  if (C3CM<-1.)
650  {
652  ev->AddParticle(*p);
653  return;
654  }
655 
656  // create scattered particle
657  GHepParticle * t = new GHepParticle(*p);
658  t->SetPdgCode(tcode);
659  double Mt = t->Mass();
660  //t->SetMomentum(TLorentzVector(0,0,0,Mt));
661 
662  // handle fermi momentum
663  if(fDoFermi)
664  {
665  // Handle fermi target
666  Target target(ev->TargetNucleus()->Pdg());
667 
668  target.SetHitNucPdg(tcode);
670  TVector3 tP3L = fFermiFac * fNuclmodel->Momentum3();
671  double tE = TMath::Sqrt(tP3L.Mag2() + Mt*Mt);
672  t->SetMomentum(TLorentzVector(tP3L,tE));
673  }
674  else
675  {
676  t->SetMomentum(TLorentzVector(0,0,0,Mt));
677  }
678 
679  bool pass = utils::intranuke2014::TwoBodyCollision(ev,pcode,tcode,scode,s2code,C3CM,
681 
682 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
683  LOG("HNIntranuke2014",pDEBUG)
684  << "|p3| = " << P3L << ", E3 = " << E3L;
685  LOG("HNIntranuke2014",pDEBUG)
686  << "|p4| = " << P4L << ", E4 = " << E4L;
687 #endif
688 
689  if (pass==true)
690  {
691 
692  //utils::intranuke2014::StepParticle(p,fFreeStep,fTrackingRadius); //not needed
693  //utils::intranuke2014::StepParticle(t,fFreeStep,fTrackingRadius);
694 
695  ev->AddParticle(*p);
696  ev->AddParticle(*t);
697  } else
698  {
699 
700  delete t; //fixes a troublesome memory leak
701 
702  LOG("HNIntranuke2014", pINFO) << "Elastic in hN failed calling TwoBodyCollision";
704  exception.SetReason("hN scattering kinematics through TwoBodyCollision failed");
705  throw exception;
706  }
707 
708  delete t;
709 
710 }
711 //___________________________________________________________________________
713 {
714  // Aaron Meyer (Jan 2010)
715  // Updated version of InelasticHN
716 
717  GHepParticle * s1 = new GHepParticle(*p);
718  GHepParticle * s2 = new GHepParticle(*p);
719  GHepParticle * s3 = new GHepParticle(*p);
720 
722  {
723  // set status of particles and return
724 
728 
729  ev->AddParticle(*s1);
730  ev->AddParticle(*s2);
731  ev->AddParticle(*s3);
732  }
733  else
734  {
735 
736  delete s1; //fixes potential memory leak
737  delete s2;
738  delete s3;
739 
740  LOG("HNIntranuke2014", pNOTICE) << "Error: could not create pion production final state";
742  exception.SetReason("PionProduction in hN failed");
743  throw exception;
744  }
745 
746  delete s1;
747  delete s2;
748  delete s3;
749  return;
750 
751 }
752 //___________________________________________________________________________
754 {
755  // This function handles pion photoproduction reactions
756 
757 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
758  LOG("HNIntranuke2014", pDEBUG)
759  << "GammaInelasticHN() is invoked for a : " << p->Name()
760  << " whose fate is : " << INukeHadroFates::AsString(fate);
761 #endif
762 
763  if(fate!=kIHNFtInelas && p->Pdg()!=kPdgGamma)
764  {
765  LOG("HNIntranuke2014", pWARN)
766  << "GammaInelasticHN() cannot handle fate: " << INukeHadroFates::AsString(fate);
767  return;
768  }
769 
770  // random number generator
771  RandomGen * rnd = RandomGen::Instance();
772 
773  // vars for incoming particle, target, and scattered reaction products
774  double ppcnt = (double) fRemnZ / (double) fRemnA; // % of protons
775  int pcode = p->Pdg();
776  int tcode = (rnd->RndFsi().Rndm()<=ppcnt)?(kPdgProton):(kPdgNeutron);
777  int scode, s2code;
778  double ke = p->KinE() / units::MeV;
779 
780  LOG("HNIntranuke2014", pNOTICE)
781  << "Particle code: " << pcode << ", target: " << tcode;
782 
783 
784  if (rnd->RndFsi().Rndm() * (fHadroData2014 -> XSecGamp_fs() -> Evaluate(ke) +
785  fHadroData2014 -> XSecGamn_fs() -> Evaluate(ke) )
786  <= fHadroData2014 -> XSecGamp_fs() -> Evaluate(ke) ) { scode = kPdgProton; }
787  else { scode = kPdgNeutron; }
788 
789  //scode=fHadroData2014->AngleAndProduct(p,tcode,C3CM,fate);
790  //double C3CM = 0.0; // cos of scattering angle
791  double C3CM = fHadroData2014->IntBounce(p,tcode,scode,fate);
792 
793  if ((tcode == kPdgProton ) && (scode==kPdgProton )) {s2code=kPdgPi0;}
794  else if ((tcode == kPdgProton ) && (scode==kPdgNeutron)) {s2code=kPdgPiP;}
795  else if ((tcode == kPdgNeutron) && (scode==kPdgProton )) {s2code=kPdgPiM;}
796  else if ((tcode == kPdgNeutron) && (scode==kPdgNeutron)) {s2code=kPdgPi0;}
797  else {
798  LOG("HNIntranuke2014", pERROR)
799  << "Error: could not determine particle final states";
800  ev->AddParticle(*p);
801  return;
802  }
803 
804  LOG("HNIntranuke2014", pNOTICE)
805  << "GammaInelastic fate: " << INukeHadroFates::AsString(fate);
806  LOG("HNIntranuke2014", pNOTICE)
807  << " final state: " << scode << " and " << s2code;
808  LOG("HNIntranuke2014", pNOTICE)
809  << " scattering angle: " << C3CM;
810 
811  GHepParticle * t = new GHepParticle(*p);
812  t->SetPdgCode(tcode);
813  double Mt = t->Mass();
814 
815  // handle fermi momentum
816  if(fDoFermi)
817  {
818  // Handle fermi target
819  Target target(ev->TargetNucleus()->Pdg());
820 
821  target.SetHitNucPdg(tcode);
823  TVector3 tP3L = fFermiFac * fNuclmodel->Momentum3();
824  double tE = TMath::Sqrt(tP3L.Mag2() + Mt*Mt);
825  t->SetMomentum(TLorentzVector(tP3L,tE));
826  }
827  else
828  {
829  t->SetMomentum(TLorentzVector(0,0,0,Mt));
830  }
831 
832  bool pass = utils::intranuke2014::TwoBodyCollision(ev,pcode,tcode,scode,s2code,C3CM,
834 
835 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
836  LOG("HNIntranuke2014",pDEBUG)
837  << "|p3| = " << P3L << ", E3 = " << E3L;
838  LOG("HNIntranuke2014",pDEBUG)
839  << "|p4| = " << P4L << ", E4 = " << E4L;
840 #endif
841 
842  if (pass==true)
843  {
844  //p->SetStatus(kIStStableFinalState);
845  //t->SetStatus(kIStStableFinalState);
846  ev->AddParticle(*p);
847  ev->AddParticle(*t);
848  } else
849  {
850  ev->AddParticle(*p);
851  }
852 
853  delete t;
854 
855 }
856 //___________________________________________________________________________
858 {
859 
860  // handle compound nucleus option
861  // -- Call the PreEquilibrium function
863  { // random number generator
864  RandomGen * rnd = RandomGen::Instance();
865 
866  double rpreeq = rnd->RndFsi().Rndm(); // sdytman test
867 
868  if((p->KinE() < fEPreEq) )
869  {
870  if(fRemnA>5&&rpreeq<0.12)
871  {
872  GHepParticle * sp = new GHepParticle(*p);
873  sp->SetFirstMother(mom);
876  delete sp;
877  return 2;
878  }
879  else
880  {
881  // nothing left to interact with!
882  LOG("HNIntranuke2014", pNOTICE)
883  << "*** Nothing left to interact with, escaping.";
884  GHepParticle * sp = new GHepParticle(*p);
885  sp->SetFirstMother(mom);
887  ev->AddParticle(*sp);
888  delete sp;
889  return 1;
890  }
891  }
892  }
893  return 0;
894 }
895 //___________________________________________________________________________
897 {
898  return (this->HandleCompoundNucleus(ev,p,p->FirstMother())==2);
899 }
900 //___________________________________________________________________________
902 {
904  const Registry * gc = confp->GlobalParameterList();
905 
906  // load hadronic cross sections
908 
909  // fermi momentum setup
911  fNuclmodel = dynamic_cast<const NuclearModelI *>
912  (fAlgf->GetAlgorithm("genie::FGMBodekRitchie","Default"));
913 
914  // other intranuke config params
915  fR0 = fConfig->GetDoubleDef ("R0", gc->GetDouble("NUCL-R0")); // fm
916  fNR = fConfig->GetDoubleDef ("NR", gc->GetDouble("NUCL-NR"));
917  fNucRmvE = fConfig->GetDoubleDef ("NucRmvE", gc->GetDouble("INUKE-NucRemovalE")); // GeV
918  fDelRPion = fConfig->GetDoubleDef ("DelRPion", gc->GetDouble("HNINUKE-DelRPion"));
919  fDelRNucleon = fConfig->GetDoubleDef ("DelRNucleon", gc->GetDouble("HNINUKE-DelRNucleon"));
920  fHadStep = fConfig->GetDoubleDef ("HadStep", gc->GetDouble("INUKE-HadStep")); // fm
921  fNucAbsFac = fConfig->GetDoubleDef ("NucAbsFac", gc->GetDouble("INUKE-NucAbsFac"));
922  fNucQEFac = fConfig->GetDoubleDef ("NucQEFac", gc->GetDouble("INUKE-NucQEFac"));
923  fNucCEXFac = fConfig->GetDoubleDef ("NucCEXFac", gc->GetDouble("INUKE-NucCEXFac"));
924  fEPreEq = fConfig->GetDoubleDef ("EPreEq", gc->GetDouble("INUKE-Energy_Pre_Eq"));
925  fFermiFac = fConfig->GetDoubleDef ("FermiFac", gc->GetDouble("INUKE-FermiFac"));
926  fFermiMomentum = fConfig->GetDoubleDef ("FermiMomentum",gc->GetDouble("INUKE-FermiMomentum"));
927  fDoFermi = fConfig->GetBoolDef ("DoFermi", gc->GetBool("INUKE-DoFermi"));
928  fFreeStep = fConfig->GetDoubleDef ("FreeStep", gc->GetDouble("INUKE-FreeStep"));
929  fDoCompoundNucleus = fConfig->GetBoolDef ("DoCompoundNucleus", gc->GetBool("INUKE-DoCompoundNucleus"));
930 
931 
932  // report
933  LOG("HNIntranuke2014", pINFO) << "Settings for Intranuke2014 mode: " << INukeMode::AsString(kIMdHN);
934  LOG("HNIntranuke2014", pWARN) << "R0 = " << fR0 << " fermi";
935  LOG("HNIntranuke2014", pWARN) << "NR = " << fNR;
936  LOG("HNIntranuke2014", pWARN) << "DelRPion = " << fDelRPion;
937  LOG("HNIntranuke2014", pWARN) << "DelRNucleon = " << fDelRNucleon;
938  LOG("HNIntranuke2014", pWARN) << "HadStep = " << fHadStep << " fermi";
939  LOG("HNIntranuke2014", pWARN) << "NucAbsFac = " << fNucAbsFac;
940  LOG("HNIntranuke2014", pWARN) << "NucQEFac = " << fNucQEFac;
941  LOG("HNIntranuke2014", pWARN) << "NucCEXFac = " << fNucCEXFac;
942  LOG("HNIntranuke2014", pWARN) << "EPreEq = " << fEPreEq;
943  LOG("HNIntranuke2014", pWARN) << "FermiFac = " << fFermiFac;
944  LOG("HNIntranuke2014", pWARN) << "FreeStep = " << fFreeStep; // free step in fm
945  LOG("HNIntranuke2014", pWARN) << "FermiMomtm = " << fFermiMomentum;
946  LOG("HNIntranuke2014", pWARN) << "DoFermi? = " << ((fDoFermi)?(true):(false));
947  LOG("HNIntranuke2014", pWARN) << "DoCmpndNuc? = " << ((fDoCompoundNucleus)?(true):(false));
948 }
949 //___________________________________________________________________________
static string AsString(INukeMode_t mode)
Definition: INukeMode.h:41
RgDbl GetDoubleDef(RgKey key, RgDbl def_opt, bool set_def=true)
Definition: Registry.cxx:549
An exception thrown by SimulateHadronState for kinematics problems. TwoBodyCollision/Kinematics used ...
void SetFirstMother(int m)
Definition: GHepParticle.h:131
Basic constants.
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:148
bool HandleCompoundNucleusHN(GHepRecord *ev, GHepParticle *p) const
TRandom3 & RndFsi(void) const
rnd number generator used by intranuclear cascade monte carlos
Definition: RandomGen.h:60
void SimulateHadronicFinalState(GHepRecord *ev, GHepParticle *p) const
#include "Numerical/GSFunc.h"
Definition: AlgCmp.h:26
AlgFactory * fAlgf
algorithm factory instance
Definition: Intranuke2014.h:94
#define pERROR
Definition: Messenger.h:50
double fNR
param multiplying the nuclear radius, determining how far to track hadrons beyond the "nuclear bounda...
bool TwoBodyCollision(GHepRecord *ev, int pcode, int tcode, int scode, int s2code, double C3CM, GHepParticle *p, GHepParticle *t, int &RemnA, int &RemnZ, TLorentzVector &RemnP4, EINukeMode mode=kIMdHA)
Intranuke utility functions.
double E(void) const
Get energy.
Definition: GHepParticle.h:90
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:79
A singleton class holding all configuration registries built while parsing all loaded XML configurati...
Definition: AlgConfigPool.h:41
double fEPreEq
threshold for pre-equilibrium reaction
static const double MeV
Definition: Units.h:130
double fR0
effective nuclear size param
int fRemnZ
remnant nucleus Z
Definition: Intranuke2014.h:97
bool IsInNucleus(const GHepParticle *p) const
Pure abstract base class. Defines the NuclearModelI interface to be implemented by any physics model ...
Definition: NuclearModelI.h:42
INukeHadroData2014 * fHadroData2014
a collection of h+N,h+A data & calculations
Definition: Intranuke2014.h:93
void SetMomentum(const TLorentzVector &p4)
double Mass(void) const
Mass that corresponds to the PDG code.
RgDbl GetDouble(RgKey key) const
Definition: Registry.cxx:488
double fFermiMomentum
whether or not particle collision is pauli blocked
Definition: tf_graph.h:23
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:30
virtual void ProcessEventRecord(GHepRecord *event_rec) const
int LastMother(void) const
Definition: GHepParticle.h:68
int Pdg(void) const
Definition: GHepParticle.h:64
int FirstMother(void) const
Definition: GHepParticle.h:67
string Name(void) const
Name that corresponds to the PDG code.
int fRemnA
remnant nucleus A
Definition: Intranuke2014.h:96
void InelasticHN(GHepRecord *ev, GHepParticle *p) const
double fDelRPion
factor by which Pion Compton wavelength gets multiplied to become nuclear size enhancement ...
double fFreeStep
produced particle free stem, in fm
const HLTPathStatus pass
const NuclearModelI * fNuclmodel
nuclear model used to generate fermi momentum
Definition: Intranuke2014.h:95
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:87
double FateWeight(int pdgc, INukeFateHN_t fate) const
void PreEquilibrium(GHepRecord *ev, GHepParticle *p, int &RemnA, int &RemnZ, TLorentzVector &RemnP4, bool DoFermi, double FermiFac, const NuclearModelI *Nuclmodel, double NucRmvE, EINukeMode mode=kIMdHN)
const int kPdgGamma
Definition: PDGCodes.h:163
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:41
double Frac(int hpdgc, INukeFateHN_t fate, double ke, int targA=0, int targZ=0) const
const Algorithm * GetAlgorithm(const AlgId &algid)
Definition: AlgFactory.cxx:86
const int kPdgKP
Definition: PDGCodes.h:146
const int kPdgPiP
Definition: PDGCodes.h:132
const int kPdgPi0
Definition: PDGCodes.h:134
virtual GHepParticle * TargetNucleus(void) const
Definition: GHepRecord.cxx:323
void ProcessEventRecord(GHepRecord *event_rec) const
#define pINFO
Definition: Messenger.h:53
bool PionProduction(GHepRecord *ev, GHepParticle *p, GHepParticle *s1, GHepParticle *s2, GHepParticle *s3, int &RemnA, int &RemnZ, TLorentzVector &RemnP4, bool DoFermi, double FermiFac, double FermiMomentum, const NuclearModelI *Nuclmodel)
Misc GENIE control constants.
enum genie::EINukeFateHN_t INukeFateHN_t
Registry * GlobalParameterList(void) const
#define pWARN
Definition: Messenger.h:51
double fFermiFac
testing parameter to modify fermi momentum
p
Definition: test.py:228
double fNucCEXFac
charge exchange xsec correction factor (hN Mode)
double fNucRmvE
binding energy to subtract from cascade nucleons
bool fDoCompoundNucleus
whether or not to do compound nucleus considerations
double KinE(bool mass_from_pdg=false) const
Get kinetic energy.
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:45
void GammaInelasticHN(GHepRecord *ev, GHepParticle *p, INukeFateHN_t fate) const
static AlgFactory * Instance()
Definition: AlgFactory.cxx:75
void SetLastMother(int m)
Definition: GHepParticle.h:132
Singleton class to load & serve a TDatabasePDG.
Definition: PDGLibrary.h:27
void SetStatus(GHepStatus_t s)
Definition: GHepParticle.h:125
virtual TVector3 Momentum3(void) const
Definition: NuclearModelI.h:59
static string AsString(INukeFateHN_t fate)
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:66
void AbsorbHN(GHepRecord *ev, GHepParticle *p, INukeFateHN_t fate) const
bool IsNeutronOrProton(int pdgc)
Definition: PDGUtils.cxx:321
TParticlePDG * Find(int pdgc)
Definition: PDGLibrary.cxx:63
RgBool GetBool(RgKey key) const
Definition: Registry.cxx:474
const int kPdgPiM
Definition: PDGCodes.h:133
string PrintFramedMesg(string mesg, unsigned int nl=1, const char f='*')
Definition: PrintUtils.cxx:171
virtual bool GenerateNucleon(const Target &) const =0
Registry * fConfig
config. (either owned or pointing to config pool)
Definition: Algorithm.h:122
virtual void AddParticle(const GHepParticle &p)
Definition: GHepRecord.cxx:526
int HandleCompoundNucleus(GHepRecord *ev, GHepParticle *p, int mom) const
double fDelRNucleon
factor by which Nucleon Compton wavelength gets multiplied to become nuclear size enhancement ...
const int kPdgProton
Definition: PDGCodes.h:62
#define pNOTICE
Definition: Messenger.h:52
static const unsigned int kRjMaxIterations
Definition: Controls.h:27
void ElasHN(GHepRecord *ev, GHepParticle *p, INukeFateHN_t fate) const
TLorentzVector * P4(void) const
Definition: GHepParticle.h:79
const int kPdgNeutron
Definition: PDGCodes.h:64
GENIE&#39;s GHEP MC event record.
Definition: GHepRecord.h:46
static const double kPi
Definition: Constants.h:38
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils...
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:40
RgBool GetBoolDef(RgKey key, RgBool def_opt, bool set_def=true)
Definition: Registry.cxx:539
double fNucAbsFac
absorption xsec correction factor (hN Mode)
double fHadStep
step size for intranuclear hadron transport
INukeFateHN_t HadronFateHN(const GHepParticle *p) const
TLorentzVector fRemnP4
P4 of remnant system.
Definition: Intranuke2014.h:98
Root of GENIE utility namespaces.
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
static AlgConfigPool * Instance()
static INukeHadroData2014 * Instance(void)
double IntBounce(const GHepParticle *p, int target, int s1, INukeFateHN_t fate)
#define pDEBUG
Definition: Messenger.h:54
void SetPdgCode(int c)
bool fDoFermi
whether or not to do fermi mom.