HG4BertCascIntranuke.cxx
Go to the documentation of this file.
1 #include "Framework/Conventions/GBuild.h"
2 #ifdef __GENIE_GEANT4_INTERFACE_ENABLED__
3 //____________________________________________________________________________
4 /*
5 
6  Author: Dennis Wright <dwright@slac.stanford.edu>
7  31 January 2017
8 
9  See header file for class documentation.
10 
11 */
12 //____________________________________________________________________________
13 
14 #include <cstdlib>
15 #include <sstream>
16 #include <exception>
17 
18 #include "TMath.h"
19 
35 
36 
37 // Geant4 headers
38 #include "Geant4/G4ParticleTypes.hh"
39 #include "Geant4/G4ParticleTable.hh"
40 #include "Geant4/G4IonTable.hh"
41 #include "Geant4/G4LeptonConstructor.hh"
42 #include "Geant4/G4MesonConstructor.hh"
43 #include "Geant4/G4BaryonConstructor.hh"
44 #include "Geant4/G4GenericIon.hh"
45 #include "Geant4/G4ProcessManager.hh"
46 #include "Geant4/G4SystemOfUnits.hh"
47 #include "Geant4/G4ThreeVector.hh"
48 #include "Geant4/G4Fancy3DNucleus.hh"
49 #include "Geant4/G4InuclParticle.hh"
50 #include "Geant4/G4InuclCollider.hh"
51 #include "Geant4/G4InuclElementaryParticle.hh"
52 #include "Geant4/G4InuclNuclei.hh"
53 #include "Geant4/G4KineticTrackVector.hh"
54 
55 
56 using std::ostringstream;
57 
58 using namespace genie;
59 using namespace genie::utils;
60 using namespace genie::utils::intranuke;
61 using namespace genie::constants;
62 using namespace genie::controls;
63 
64 //___________________________________________________________________________
65 HG4BertCascIntranuke::HG4BertCascIntranuke()
66  : EventRecordVisitorI("genie::HG4BertCascIntranuke")
67 {
68  InitG4Particles();
69 }
70 
71 //___________________________________________________________________________
72 HG4BertCascIntranuke::HG4BertCascIntranuke(string config)
73  : EventRecordVisitorI("genie::HG4BertCascIntranuke", config)
74 {
75  InitG4Particles();
76 }
77 //___________________________________________________________________________
78 HG4BertCascIntranuke::~HG4BertCascIntranuke()
79 {
80 }
81 //___________________________________________________________________________
82 int HG4BertCascIntranuke::G4BertCascade(GHepRecord * evrec) const{
83  GHepParticle* probe = evrec->Probe(); // incoming particle
84  GHepParticle* tgtNucl = evrec->TargetNucleus(); // target
85 
86 
87  const G4ParticleDefinition* incidentDef = PDGtoG4Particle(probe->Pdg() );
88 
89  int Zinit = tgtNucl->Z();
90  int Ainit = tgtNucl->A();
91 
92  G4InuclNuclei *theNucleus = new G4InuclNuclei(Ainit,Zinit);
93 
94 
95  const G4ThreeVector tgtmomentum(0.,0.,0.);
96  const G4double tgtEnergy =tgtNucl->Energy()*1000;
97  G4LorentzVector tgt4momentum(tgtmomentum,tgtEnergy);
98 
99  //
100  TLorentzVector *p4Genie = probe->P4();
101  //
102 
103  G4ThreeVector incidentDir(p4Genie->Vect().Unit().Px(),p4Genie->Vect().Unit().Py(),p4Genie->Vect().Unit().Pz());
104 
105  double dynamicMass = std::sqrt(p4Genie->M2() );
106  double incidentKE = p4Genie->E() - dynamicMass;
107  const G4DynamicParticle p(incidentDef, incidentDir, incidentKE/units::MeV, dynamicMass/units::MeV);
108  G4InuclElementaryParticle* incident = new G4InuclElementaryParticle(p,G4InuclParticle::bullet);
109 
110 
111  this->SetTrackingRadius(tgtNucl);
112  this->GenerateVertex(evrec);
113  fRemnA=Ainit;
114  fRemnZ=Zinit;
115 
116  GHepParticle * sp = new GHepParticle(*probe);
117  bool has_interacted = false;
118  while ( this-> IsInNucleus(sp) ) {
119  // advance the hadron by a step
121 
122  // check whether it interacts
123  double d = this->GenerateStep(evrec,sp);
124  has_interacted = (d<fHadStep);
125  if(has_interacted) break;
126  } //stepping
127 
128  if ( has_interacted ) {
129 
130  // Set up output and start the cascade
131  G4CollisionOutput cascadeOutput;
132  G4InuclCollider bertCollider;
133  bertCollider.useCascadeDeexcitation();
134  //collide
135  bertCollider.collide(incident,theNucleus,cascadeOutput);
136 
137  delete incident;
138  delete theNucleus;
139 
140  //
141  // Add Geant4 generated particles to the event record
142  //
143  TLorentzVector remX(0.,0.,0.,0.);
144  int Nfrag = cascadeOutput.numberOfOutgoingNuclei();
145  const std::vector<G4InuclNuclei>& outgoingFragments =
146  cascadeOutput.getOutgoingNuclei();
147  int npdg = 0;
148  fRemnZ = 0;
149  fRemnA = 0;
150 
151 
152  // Now the single hadrons
153  int Nhad = cascadeOutput.numberOfOutgoingParticles();
154  const std::vector<G4InuclElementaryParticle>& outgoingHadrons =
155  cascadeOutput.getOutgoingParticles();
156  for (int l = 0; l < Nhad; l++) {
157  npdg = outgoingHadrons[l].getDefinition()->GetPDGEncoding();
158  G4LorentzVector hadP = outgoingHadrons[l].getMomentum();
159  TLorentzVector tempP(hadP.px(), hadP.py(), hadP.pz(), hadP.e() );
160  GHepParticle new_particle(npdg, kIStStableFinalState,
161  0, 1,-1,-1,tempP, remX);
162  evrec->AddParticle(new_particle);
163  }
164 
165  if (Nfrag > 0) {
166  // Get largest nuclear fragment in output and call it the remnant
167  int maxA = 0;
168  int rem_index = 0;
169  for (int j = 0; j < Nfrag; j++) {
170  if (outgoingFragments[j].getA() > maxA) {
171  maxA = outgoingFragments[j].getA();
172  rem_index = j;
173  }
174  }
175 
176  fRemnZ = outgoingFragments[rem_index].getZ();
177  fRemnA = outgoingFragments[rem_index].getA();
178 
179  // Get remnant momentum from cascade
180  G4LorentzVector nucP = outgoingFragments[rem_index].getMomentum();
181  TLorentzVector remP(nucP.px(), nucP.py(), nucP.pz(), nucP.e() );
182  npdg = outgoingFragments[rem_index].getDefinition()->GetPDGEncoding();
183  GHepParticle largest_Fragment(npdg, kIStFinalStateNuclearRemnant,
184  1,0,-1,-1, remP, remX);
185  evrec->AddParticle(largest_Fragment);
186 
187  // If any nuclear fragments left, add them to the event
188  for (G4int k = 0; k < Nfrag; k++) {
189  if (k != rem_index) {
190  npdg = outgoingFragments[k].getDefinition()->GetPDGEncoding();
191  nucP = outgoingFragments[k].getMomentum(); // need to boost by fRemnP4
192  TLorentzVector tempP(nucP.px(), nucP.py(), nucP.pz(), nucP.e() );
193  GHepParticle nuclear_Fragment(npdg, kIStStableFinalState, 0, 1,-1,-1,
194  tempP, remX);
195  evrec->AddParticle(nuclear_Fragment);
196  }
197  }
198  } // Nfrag > 0
199 
200  } else { // transparent event
201 
202  TLorentzVector p4h (0.,0.,probe->Pz(),probe->E());
203  TLorentzVector x4null(0.,0.,0.,0.);
204  TLorentzVector p4tgt (0.,0.,0.,tgtNucl->Mass());
205  evrec->AddParticle(probe->Pdg(), kIStStableFinalState,
206  0,-1,-1,-1, p4h,x4null);
207  evrec->AddParticle(tgtNucl->Pdg(),kIStFinalStateNuclearRemnant,
208  1,-1,-1,-1,p4tgt,x4null);
209  }
210  delete sp;
211  return 0;
212 }
213 //___________________________________________________________________________
214 double HG4BertCascIntranuke::GenerateStep(GHepRecord* /*evrec*/,
215  GHepParticle* p) const {
216  // Generate a step (in fermis) for particle p in the input event.
217  // Computes the mean free path L and generate an 'interaction' distance d
218  // from an exp(-d/L) distribution
219 
220  RandomGen * rnd = RandomGen::Instance();
221 
222  double LL = utils::intranuke2018::MeanFreePath(p->Pdg(), *p->X4(), *p->P4(),
223  fRemnA, fRemnZ,
226  fXsecNNCorr);
227 
228  double d = -1.*LL * TMath::Log(rnd->RndFsi().Rndm());
229 
230  return d;
231 }
232 //___________________________________________________________________________
233 void HG4BertCascIntranuke::ProcessEventRecord(GHepRecord* evrec) const
234 {
235  // Check the event generation mode: should be lepton-nucleus
236  fGMode = evrec->EventGenerationMode();
237  if ( fGMode== kGMdLeptonNucleus) {
238  LOG("HG4BertCascIntranuke", pINFO)
239  << " Lepton-nucleus event generation mode ";
240  GHepParticle* nucltgt = evrec->TargetNucleus();
241  if (nucltgt) {
242  // Decide tracking radius for the current nucleus (few * R0 * A^1/3)
243  SetTrackingRadius(nucltgt);
244  } else {
245  LOG("HG4BertCascIntranuke", pINFO) << "No nuclear target found - exit ";
246  return;
247  }
248 
249  // Collect hadrons from initial interaction and track them through the
250  // nucleus using Bertini cascade
251  TransportHadrons(evrec);
252  } else if(fGMode == kGMdHadronNucleus ||
254  G4BertCascade(evrec);
255  } else{
256  LOG("HG4BertCascIntranuke", pINFO)
257  << " Inappropriate event generation mode - exit ";
258  return;
259  }
260 
261  LOG("HG4BertCascIntranuke", pINFO) << "Done with this event";
262 }
263 
264 //___________________________________________________________________________
265 void HG4BertCascIntranuke::InitG4Particles() const
266 {
267  G4LeptonConstructor::ConstructParticle();
268  G4MesonConstructor::ConstructParticle();
269  G4BaryonConstructor::ConstructParticle();
270  G4He3::He3();
271  G4Gamma::Gamma();
272  G4ParticleTable* pTable = G4ParticleTable::GetParticleTable();
273  pTable->SetReadiness(true);
274  G4GenericIon* gIon = G4GenericIon::Definition();
275  gIon->SetProcessManager(new G4ProcessManager(gIon) );
276 
277  // testing
278  const G4ParticleDefinition* electron = PDGtoG4Particle(11);
279  const G4ParticleDefinition* proton = PDGtoG4Particle(2212);
280  const G4ParticleDefinition* neutron = PDGtoG4Particle(2112);
281  const G4ParticleDefinition* piplus = PDGtoG4Particle(211);
282  LOG("HG4BertCascIntranuke", pINFO)
283  << "testing initialization of G4 particles \n"
284  << " e 0x" << electron << "\n"
285  << " p 0x" << proton << "\n"
286  << " n 0x" << neutron << "\n"
287  << " pi+ 0x" << piplus << "\n"
288  << "...InitG4Particles complete";
289  if ( electron == 0 || proton == 0 || neutron == 0 || piplus == 0 ) {
290  LOG("HG4BertCascIntranuke", pFATAL)
291  << "something is seriously wrong with g4 particle lookup";
292  exit(42);
293  }
294 }
295 
296 //___________________________________________________________________________
297 void HG4BertCascIntranuke::GenerateVertex(GHepRecord * evrec) const
298 {
299  // Sets a vertex in the nucleus periphery
300  // Called onlt in hadron/photon-nucleus interactions.
301 
302  GHepParticle * nucltgt = evrec->TargetNucleus();
303  assert(nucltgt);
304 
305  RandomGen * rnd = RandomGen::Instance();
306  TVector3 vtx(999999.,999999.,999999.);
307 
308  // *** For h+A events (test mode):
309  // Assume a hadron beam with uniform intensity across an area,
310  // so we need to choose events uniformly within that area.
311  double x=999999., y=999999., epsilon = 0.001;
312  double R2 = TMath::Power(fTrackingRadius,2.);
313  double rp2 = TMath::Power(x,2.) + TMath::Power(y,2.);
314  while(rp2 > R2-epsilon) {
315  x = (fTrackingRadius-epsilon) * rnd->RndFsi().Rndm();
316  y = -fTrackingRadius + 2*fTrackingRadius * rnd->RndFsi().Rndm();
317  y -= ((y>0) ? epsilon : -epsilon);
318  rp2 = TMath::Power(x,2.) + TMath::Power(y,2.);
319  }
320  vtx.SetXYZ(x,y, -1.*TMath::Sqrt(TMath::Max(0.,R2-rp2)) + epsilon);
321 
322  // get the actual unit vector along the incoming hadron direction
323  TVector3 direction = evrec->Probe()->P4()->Vect().Unit();
324 
325  // rotate the vtx position
326  vtx.RotateUz(direction);
327 
328  LOG("HG4BertCascIntranuke", pNOTICE)
329  << "Generated vtx @ R = " << vtx.Mag() << " fm / "
330  << print::Vec3AsString(&vtx);
331 
332  TObjArrayIter piter(evrec);
333  GHepParticle * p = 0;
334  while( (p = (GHepParticle *) piter.Next()) ) {
335  if ( pdg::IsPseudoParticle(p->Pdg()) ) continue;
336  if ( pdg::IsIon(p->Pdg()) ) continue;
337 
338  p->SetPosition(vtx.x(), vtx.y(), vtx.z(), 0.);
339  }
340 }
341 
342 //___________________________________________________________________________
343 void HG4BertCascIntranuke::SetTrackingRadius(const GHepParticle * p) const
344 {
345  assert(p && pdg::IsIon(p->Pdg()));
346  double A = p->A();
347  fTrackingRadius = fR0*TMath::Power(A, 1./3.);
348 
349  // multiply that by some input factor so that hadrons are tracked
350  // beyond the nuclear 'boundary' since the nuclear density distribution
351  // is not zero there
352  fTrackingRadius *= fNR;
353 
354  LOG("HG4BertCascIntranuke", pNOTICE)
355  << "Setting tracking radius to R = " << fTrackingRadius;
356 }
357 
358 //___________________________________________________________________________
359 bool HG4BertCascIntranuke::IsInNucleus(const GHepParticle * p) const
360 {
361  // check whether the input particle is still within the nucleus
362  return (p->X4()->Vect().Mag() < fTrackingRadius + fHadStep);
363 }
364 //___________________________________________________________________________
365 void HG4BertCascIntranuke::TransportHadrons(GHepRecord * evrec) const
366 {
367  // Set up nucleus through which hadrons will be propagated
368  // In frozen nucleus approximation it is the remnant nucleus corrected for
369  // final state lepton charge
370 
371  //rwh-unused//GHepParticle* probe = evrec->Probe(); // incoming neutrino
372  GHepParticle* tgtNucl = evrec->TargetNucleus();
373  GHepParticle* remNucl = evrec->RemnantNucleus();
374  GHepParticle* outLept = evrec->FinalStatePrimaryLepton(); // outgoing lepton
375  GHepParticle* struckNucleon = evrec->HitNucleon();
376 
377  int inucl = evrec->RemnantNucleusPosition();
378 
379  // double sepEnergy = tgtNucl->Mass() - remNucl->Mass() - struckNucleon->Mass();
380  // std::cout << " sepEnergy = " << sepEnergy << std::endl;
381 
382  GHepParticle* p = 0;
383  const G4ParticleDefinition* incidentDef = 0;
384  GHepParticle* incidentBaryon = 0;
385  TObjArrayIter piter(evrec);
386  TObjArrayIter pitter(evrec);
387  int icurr =-1;
388  bool has_incidentBaryon(false), has_secondaries(false);
389  //rwh-unused// bool transparents(false),
390  bool has_remnant(false), has_incidentparticle(false);
391 
392  fRemnA=remNucl->A();
393  fRemnZ=remNucl->Z();
394 
395  while( (p = (GHepParticle*) piter.Next()) ) {
396  icurr++;
397  bool has_interacted(false);
398  if ( ! this->NeedsRescattering(p) ) continue;
399 
400  GHepParticle * sp = new GHepParticle(*p);
401  sp->SetFirstMother(icurr);
402 
403  if ( ! this->CanRescatter(sp) ) {
405  evrec->AddParticle(*sp);
406  evrec->Particle(sp->FirstMother())->SetRescatterCode(1);
407  delete sp;
408  continue;
409  }
410 
411  while ( this-> IsInNucleus(sp) ) {
412  // advance the hadron by a step
414  double d = this->GenerateStep(evrec,sp);
415  has_interacted = (d<fHadStep);
416  if ( has_interacted ) {
417  has_secondaries=true;
418  break;
419  }
420  } //stepping
421  if ( ! has_interacted ) {
423  evrec->AddParticle(*sp);
424  evrec->Particle(sp->FirstMother())->SetRescatterCode(1);
425  //rwh-unused//transparents=true;
426  }
427  if ( ! has_incidentBaryon && sp->Status() == kIStHadronInTheNucleus ) {
428  /*
429  if (sp->Pdg() == kPdgProton ||
430  sp->Pdg() == kPdgNeutron||
431  sp->Pdg() == kPdgLambda ||
432  sp->Pdg() == kPdgSigmaP ||
433  sp->Pdg() == kPdgSigma0 ||
434  sp->Pdg() == kPdgSigmaM ) {
435  */
436  if ( IsBaryon(sp) ) {
437  incidentBaryon = sp;
438  incidentDef = PDGtoG4Particle(sp->Pdg() );
439  has_incidentBaryon=true;
440  } else {
441  if (sp->Pdg() == kPdgProton ||
442  sp->Pdg() == kPdgNeutron||
443  sp->Pdg() == kPdgLambda ||
444  sp->Pdg() == kPdgSigmaP ||
445  sp->Pdg() == kPdgSigma0 ||
446  sp->Pdg() == kPdgSigmaM ) {
447  LOG("G4BertCascInterface::TransportHadrons", pWARN)
448  << "IsBaryon failed to tag " << *sp;
449  }
450  }
451  }
452  delete sp;
453  }
454 
455  if ( has_secondaries ) {
456  if ( ! incidentBaryon ) LOG("G4BertCascInterface::TransportHadrons", pINFO)
457  << "Unrecognized baryon in nucleus";
458 
459  int Zinit = remNucl->Z() - outLept->Charge()/3;
460  if (incidentBaryon) {
461  Zinit += (struckNucleon->Charge() - incidentBaryon->Charge() )/3;
462  }
463 
464  //std::cout << " Zinit = " << Zinit <<" remNucl->Z()= "
465  // <<remNucl->Z()<< std::endl;
466 
467  //rwh-noused//int Ainit = remNucl->A();
468  //std::cout << " Ainit = " << Ainit << std::endl;
469 
470  G4Fancy3DNucleus* g4Nucleus = new G4Fancy3DNucleus();
471 
472  TLorentzVector pIncident;
473 
474  g4Nucleus->Init(remNucl->A(),remNucl->Z());
475  double EE = struckNucleon->E() - tgtNucl->Mass() + g4Nucleus->GetMass()*units::MeV;
476  TLorentzVector struckMomentum(struckNucleon->Px(), struckNucleon->Py(), struckNucleon->Pz(), EE);
477  Double_t PxI(0),PyI(0),PzI(0),EEI(0);
478  int icccur=-1;
479  int pos_in_evrec(0);
480  while( (p = (GHepParticle*) pitter.Next()) ) {
481  icccur++;
482  if (p->Status() == kIStHadronInTheNucleus && this->CanRescatter(p) && p->RescatterCode()!=1) {
483  PxI+=p->P4()->Px();
484  PyI+=p->P4()->Py();
485  PzI+=p->P4()->Pz();
486  EEI+=p->P4()->E();
487  if ( pos_in_evrec==0 ) pos_in_evrec = icccur;
488  if ( ! has_incidentparticle ) { // take the baryon as incident particle
489  /*
490  if (p->Pdg() == kPdgProton ||
491  p->Pdg() == kPdgNeutron ||
492  p->Pdg() == kPdgLambda ||
493  p->Pdg() == kPdgSigmaP ||
494  p->Pdg() == kPdgSigma0 ||
495  p->Pdg() == kPdgSigmaM ) {
496  */
497  if ( IsBaryon(p) ) {
498  incidentDef = PDGtoG4Particle(p->Pdg() );
499  has_incidentparticle=true;
500  }
501  }
502  }
503  }
504  GHepParticle * pinN = evrec->Particle(pos_in_evrec);
505  if ( ! has_incidentparticle) {
506  incidentDef=PDGtoG4Particle(pinN->Pdg()); // if no baryon among the secondaries
507  }
508  pIncident.SetPxPyPzE(PxI,PyI,PzI,EEI);
509 
510 
511  G4ThreeVector incidentDir(pIncident.Vect().Unit().Px(),
512  pIncident.Vect().Unit().Py(),
513  pIncident.Vect().Unit().Pz());
514 
515  double dynamicMass = std::sqrt(pIncident.M2() );
516  double incidentKE = pIncident.E() - dynamicMass;
517  // Create pseudo-particle to supply Bertini collider with bullet
518 
519  G4DynamicParticle dp(incidentDef, incidentDir, incidentKE/units::MeV, dynamicMass/units::MeV);
520 
521  G4InuclElementaryParticle* incident = new G4InuclElementaryParticle(dp,G4InuclParticle::bullet);
522 
523  // Get hadronic secondaries and convert them to G4KineticTracks
524 
525  G4KineticTrackVector* g4secondaries = ConvertGenieSecondariesToG4(evrec);
526 
527  int Nsec = g4secondaries->size();
528  // Set up output and start the cascade
529  G4CollisionOutput cascadeOutput;
530  G4InuclCollider bertCollider;
531  bertCollider.useCascadeDeexcitation();
532  // bertCollider.setVerboseLevel(3);
533  bertCollider.rescatter(incident, g4secondaries, g4Nucleus, cascadeOutput);
534  delete incident;
535  delete g4Nucleus;
536  for (int n = 0; n < Nsec; n++) delete (*g4secondaries)[n];
537  delete g4secondaries;
538 
539  //
540  // Add Geant4 generated particles to the event record
541  //
542  TLorentzVector remX(tgtNucl->Vx(), tgtNucl->Vy(), tgtNucl->Vz(), tgtNucl->Vt() );
543 
544  int rem_nucl = evrec->RemnantNucleusPosition(); // position in array
545  int Nfrag = cascadeOutput.numberOfOutgoingNuclei();
546  const std::vector<G4InuclNuclei>& outgoingFragments = cascadeOutput.getOutgoingNuclei();
547  int npdg = 0;
548  fRemnZ = 0;
549  fRemnA = 0;
550 
551 
552  // Now the single hadrons
553  int Nhad = cascadeOutput.numberOfOutgoingParticles();
554  const std::vector<G4InuclElementaryParticle>& outgoingHadrons = cascadeOutput.getOutgoingParticles();
555  for (int l = 0; l < Nhad; l++) {
556  npdg = outgoingHadrons[l].getDefinition()->GetPDGEncoding();
557 
558  G4LorentzVector hadP = outgoingHadrons[l].getMomentum();
559  TLorentzVector tempP(hadP.px(), hadP.py(), hadP.pz(), hadP.e() );
560 
561  GHepParticle new_particle(npdg, kIStStableFinalState, -1, -1,-1,-1,tempP, remX);
562  evrec->AddParticle(new_particle);
563  }
564 
565  if (Nfrag > 0) {
566  int maxA = 0;
567  int rem_index = 0;
568  for (int j = 0; j < Nfrag; j++) {
569  if (outgoingFragments[j].getA() > maxA) {
570  maxA = outgoingFragments[j].getA();
571  rem_index = j;
572  }
573  }
574 
575  fRemnZ = outgoingFragments[rem_index].getZ();
576  fRemnA = outgoingFragments[rem_index].getA();
577 
578  // Get remnant momentum from cascade
579  G4LorentzVector nucP = outgoingFragments[rem_index].getMomentum();
580  TLorentzVector remP(nucP.px(), nucP.py(), nucP.pz(), nucP.e() );
581 
582  // If any nuclear fragments left, add them to the event
583  for (G4int k = 0; k < Nfrag; k++) {
584  if (k != rem_index) {
585  npdg = outgoingFragments[k].getDefinition()->GetPDGEncoding();
586  nucP = outgoingFragments[k].getMomentum(); // need to boost by fRemnP4
587  TLorentzVector tempP(nucP.px(), nucP.py(), nucP.pz(), nucP.e() );
588 
589  GHepParticle nuclear_Fragment(npdg, kIStStableFinalState, rem_nucl, 0,-1,-1,tempP, remX);
590  evrec->AddParticle(nuclear_Fragment);
591  }
592  }
593 
594  // Get largest nuclear fragment in output and call it the remnant
595  npdg = outgoingFragments[rem_index].getDefinition()->GetPDGEncoding();
596  remP.SetPx(remP.Px()+remNucl->P4()->Px());
597  remP.SetPy(remP.Py()+remNucl->P4()->Py());
598  remP.SetPz(remP.Pz()+remNucl->P4()->Pz());
599 
600  GHepParticle largest_Fragment(npdg, kIStFinalStateNuclearRemnant,rem_nucl,-1,-1,-1, remP, remX);
601  evrec->AddParticle(largest_Fragment);
602  } // Nfrag > 0
603  has_remnant=true;
604  }
605  // Mark the initial remnant nucleus as an intermediate state
606  if(!has_remnant){
607  GHepParticle * sp = new GHepParticle(*evrec->Particle(inucl));
608  sp->SetFirstMother(inucl);
610  evrec->AddParticle(*sp);
611  delete sp;
612  }
613  evrec->Particle(inucl)->SetStatus(kIStIntermediateState);
614 
615  // Geant4 conservation test
616  // this probably isn't configured right ... skip it for now
617  /*
618  if ( Conserve4Momentum(evrec) ) {
619  std::cout << " Energy conservation test " << std::endl;
620  }
621  */
622 
623 }
624 
625 //____________________________________________________________________________
626 bool HG4BertCascIntranuke::NeedsRescattering(const GHepParticle * p) const {
627  // checks whether the particle should be rescattered
628  assert(p);
629  // attempt to rescatter anything marked as 'hadron in the nucleus' (14)
630  return ( p->Status() == kIStHadronInTheNucleus );
631 }
632 
633 //___________________________________________________________________________
634 bool HG4BertCascIntranuke::CanRescatter(const GHepParticle * p) const
635 {
636  assert(p);
637  return ( p->Pdg() == kPdgPiP ||
638  p->Pdg() == kPdgPiM ||
639  p->Pdg() == kPdgPi0 ||
640  p->Pdg() == kPdgProton ||
641  p->Pdg() == kPdgAntiProton ||
642  p->Pdg() == kPdgNeutron ||
643  p->Pdg() == kPdgKP ||
644  p->Pdg() == kPdgKM ||
645  p->Pdg() == kPdgK0 ||
646  p->Pdg() == kPdgK0L ||
647  p->Pdg() == kPdgSigma0 ||
648  p->Pdg() == kPdgSigmaM ||
649  p->Pdg() == kPdgSigmaP ||
650  //p->Pdg() == kPdgSigmaPPc ||
651  p->Pdg() == kPdgXiM ||
652  p->Pdg() == kPdgXi0 ||
653  p->Pdg() == kPdgLambda
654  );
655 }
656 
657 //___________________________________________________________________________
658 bool HG4BertCascIntranuke::IsBaryon(const GHepParticle * p) const
659 {
660  // use PDGLibrary to determine if the particle is baryon
661  assert(p);
662 
663  TParticlePDG * ppdg = PDGLibrary::Instance()->Find(p->Pdg());
664  if ( ! ppdg ) {
665  LOG("G4BertCascInterface", pWARN)
666  << "no entry for PDG " << p->Pdg() << " in PDGLibrary";
667  } else {
668  if ( std::string(ppdg->ParticleClass()) == std::string("Baryon") ) {
669  return true;
670  }
671  }
672  return false;
673 }
674 //___________________________________________________________________________
675 const G4ParticleDefinition* HG4BertCascIntranuke::PDGtoG4Particle(int pdg) const
676 {
677  const G4ParticleDefinition* pDef = 0;
678 
679  if ( abs(pdg) < 1000000000 ) {
680  pDef = G4ParticleTable::GetParticleTable()->FindParticle(pdg);
681  } else if ( pdg < 2000000000 ) {
682  pDef = G4IonTable::GetIonTable()->GetIon(pdg);
683  }
684 
685  if ( ! pDef ) {
686  LOG("HG4BertCascIntranuke", pWARN)
687  << "Unrecognized Bertini particle type: " << pdg;
688  }
689 
690  return pDef;
691 }
692 
693 //___________________________________________________________________________
694 G4KineticTrackVector* HG4BertCascIntranuke::ConvertGenieSecondariesToG4(std::vector<GHepParticle> partList) const
695 {
696  static double GeToG4length = 2.81967*1.0e-12*1.2/1.4;
697 
698  GHepParticle* p = 0;
699  const G4ParticleDefinition* pDef = 0;
700  G4KineticTrackVector* g4secondaries = new G4KineticTrackVector;
701  G4KineticTrack* kt = 0;
702 
703  for (size_t it=0 ; it<partList.size(); it++) {
704  p= &partList.at(it);
705  pDef = PDGtoG4Particle(p->Pdg() );
706  double formationTime = p->Vt();
707  G4ThreeVector formationPosition(p->Vx()*GeToG4length,
708  p->Vy()*GeToG4length,
709  p->Vz()*GeToG4length);
710  G4LorentzVector formationMomentum(p->Px()/units::MeV, p->Py()/units::MeV,
711  p->Pz()/units::MeV, p->E()/units::MeV );
712  kt = new G4KineticTrack(pDef, formationTime, formationPosition, formationMomentum);
713  kt->SetDefinition(pDef); // defeat reset to K0L/K0S in ctor
714  kt->SetState(G4KineticTrack::inside);
715  g4secondaries->push_back(kt);
716  }
717 
718  return g4secondaries;
719 }
720 
721 //___________________________________________________________________________
722 G4KineticTrackVector* HG4BertCascIntranuke::ConvertGenieSecondariesToG4(GHepRecord* evrec) const {
723  // Get hadronic secondaries from event record and convert them to
724  // G4KineticTracks for use in cascade
725 
726  // distance units: Geant4 mm = 1.0, GENIE fermi = 1.0 => conversion 1e-12
727  // nuclear radius parameter R0: Geant4 2.81967*1.2, GENIE = 1.4 => conversion 2.41686
728  static double GeToG4length = 2.81967*1.0e-12*1.2/1.4;
729 
730  TObjArrayIter piter(evrec);
731  GHepParticle* p = 0;
732  const G4ParticleDefinition* pDef = 0;
733  G4KineticTrackVector* g4secondaries = new G4KineticTrackVector;
734  G4KineticTrack* kt = 0;
735 
736  while( (p = (GHepParticle*) piter.Next()) ) {
737  if ( p->Status() == kIStHadronInTheNucleus &&
738  this->CanRescatter(p) && ( p->RescatterCode() != 1) ) {
739  pDef = PDGtoG4Particle(p->Pdg() );
740  double formationTime = p->Vt();
741  G4ThreeVector formationPosition(p->Vx()*GeToG4length,
742  p->Vy()*GeToG4length,
743  p->Vz()*GeToG4length);
744  G4LorentzVector formationMomentum(p->Px()/units::MeV, p->Py()/units::MeV,
745  p->Pz()/units::MeV, p->E()/units::MeV );
746  kt = new G4KineticTrack(pDef, formationTime,
747  formationPosition, formationMomentum);
748  kt->SetDefinition(pDef); // defeat reset to K0L/K0S in ctor
749  kt->SetState(G4KineticTrack::inside);
750  g4secondaries->push_back(kt);
751  }
752  }
753 
754  return g4secondaries;
755 }
756 
757 //___________________________________________________________________________
758 bool HG4BertCascIntranuke::Conserve4Momentum(GHepRecord* evrec) const
759 {
760  bool failed = false;
761 
762  GHepParticle* probe = evrec->Probe(); // incoming neutrino
763  GHepParticle* targetNucleus = evrec->TargetNucleus();
764  GHepParticle* finalLepton = evrec->FinalStatePrimaryLepton();
765 
766  int Nentries = evrec->GetEntries();
767 
768  // Collect 4-momenta of final state particles
769  GHepParticle* p = 0;
770  TLorentzVector sum4mom(0.0, 0.0, 0.0, 0.0);
771  for (int j = 0; j < Nentries; j++) {
772  p = (GHepParticle*) (*evrec)[j];
773  if (p->FirstMother() == evrec->RemnantNucleusPosition() ) {
774  sum4mom += *(p->P4() );
775  /*
776  LOG("HG4BertCascIntranuke", pINFO)
777  << " Final state 4-momentum = ("
778  << p->P4()->Px() << ", " << p->P4()->Py() << ", "
779  << p->P4()->Pz() << ", " << p->P4()->E() << ") ";
780  */
781  }
782  }
783  sum4mom += *(finalLepton->P4() );
784 
785  TLorentzVector initial4mom = *(targetNucleus->P4() ) + *(probe->P4() );
786 
787  TLorentzVector diff = initial4mom - sum4mom;
788  // rwh need some threshold for acceptable differences
789  const double maxdiff = 1.0e-9; // set crazy small for now
790  double diffE = diff.E();
791  TVector3 diffp3 = diff.Vect();
792  double diffpmag = diffp3.Mag();
793  if ( TMath::Abs(diffE) > maxdiff || diffpmag > maxdiff ) {
794  failed = true;
795  LOG("HG4BertCascIntranuke", pWARN)
796  << "final state - initial state differs by > " << maxdiff << "\n"
797  << " dE = " << diffE << ", d|p3| = " << diffpmag;
798 
799  LOG("HG4BertCascIntranuke", pWARN)
800  << " Total Final state 4-momentum = (" << sum4mom.Px()
801  << ", " << sum4mom.Py()
802  << ", " << sum4mom.Pz()
803  << ", " << sum4mom.E() << ") ";
804  LOG("HG4BertCascIntranuke", pWARN)
805  << " Total Initial state 4-momentum = (" << initial4mom.Px()
806  << ", " << initial4mom.Py()
807  << ", " << initial4mom.Pz()
808  << ", " << initial4mom.E() << ") ";
809  }
810 
811  // Charge conservation
812  double Qinit = targetNucleus->Charge();
813  double Qfinal = finalLepton->Charge();
814 
815  for (int j = 0; j < Nentries; j++) {
816  p = (GHepParticle*) (*evrec)[j];
817  if ( p->FirstMother() == evrec->RemnantNucleusPosition() ) {
818  Qfinal += p->Charge();
819  }
820  }
821 
822  if (Qinit != Qfinal) {
823  failed = true;
824  LOG("HG4BertCascIntranuke", pWARN)
825  << " Charge not conserved: \n"
826  << " Qfinal = " << Qfinal
827  << " Qinit = " << Qinit;
828  }
829 
830  return ( ! failed );
831 }
832 
833 
834 //___________________________________________________________________________
835 void HG4BertCascIntranuke::LoadConfig(void)
836 {
837  // fermi momentum setup
838  // this is specifically set in Intranuke2018::Configure(string)
839  fNuclmodel = dynamic_cast<const NuclearModelI *>( this -> SubAlg("NuclearModel") ) ;
840 
841  // other intranuke config params
842  GetParam( "NUCL-R0", fR0 ); // fm
843  GetParam( "NUCL-NR", fNR );
844 
845  GetParam( "INUKE-NucRemovalE", fNucRmvE ); // GeV
846  GetParam( "INUKE-HadStep", fHadStep ) ;
847  GetParam( "INUKE-NucAbsFac", fNucAbsFac ) ;
848  GetParam( "INUKE-NucCEXFac", fNucCEXFac ) ;
849  GetParam( "INUKE-Energy_Pre_Eq", fEPreEq ) ;
850  GetParam( "INUKE-FermiFac", fFermiFac ) ;
851  GetParam( "INUKE-FermiMomentum", fFermiMomentum ) ;
852  GetParam( "INUKE-DoFermi", fDoFermi ) ;
853  GetParam( "INUKE-XsecNNCorr", fXsecNNCorr ) ;
854  GetParamDef( "UseOset", fUseOset, false ) ;
855  GetParamDef( "AltOset", fAltOset, false ) ;
856  GetParam( "HNINUKE-DelRPion", fDelRPion ) ;
857  GetParam( "HNINUKE-DelRNucleon", fDelRNucleon ) ;
858 
859  // report
860  LOG("HG4BertCascIntranuke", pNOTICE)
861  << "Settings for INTRANUKE mode: " << INukeMode::AsString(kIMdHA) << "\n"
862  << "R0 = " << fR0 << " fermi" << "\n"
863  << "NR = " << fNR << "\n"
864  << "DelRPion = " << fDelRPion << "\n"
865  << "DelRNucleon = " << fDelRNucleon << "\n"
866  << "HadStep = " << fHadStep << " fermi" << "\n"
867  << "EPreEq = " << fHadStep << " fermi" << "\n"
868  << "NucAbsFac = " << fNucAbsFac << "\n"
869  << "NucCEXFac = " << fNucCEXFac << "\n"
870  << "FermiFac = " << fFermiFac << "\n"
871  << "FermiMomtm = " << fFermiMomentum << "\n"
872  << "DoFermi? = " << ((fDoFermi)?(true):(false)) << "\n"
873  << "XsecNNCorr? = " << ((fXsecNNCorr)?(true):(false));
874 
875 }
876 
877 //___________________________________________________________________________
878 
880 {
881  Algorithm::Configure(config);
882  LoadConfig();
883 }
884 
885 //___________________________________________________________________________
886 void HG4BertCascIntranuke::Configure(string param_set)
887 {
888  Algorithm::Configure(param_set);
889  LoadConfig();
890 }
891 
892 //___________________________________________________________________________
893 #endif // __GENIE_GEANT4_INTERFACE_ENABLED__
static string AsString(INukeMode_t mode)
Definition: INukeMode.h:41
int Z(void) const
void SetFirstMother(int m)
Definition: GHepParticle.h:132
#define EE
Basic constants.
int RescatterCode(void) const
Definition: GHepParticle.h:65
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
TRandom3 & RndFsi(void) const
rnd number generator used by intranuclear cascade monte carlos
Definition: RandomGen.h:59
const int kPdgXi0
Definition: PDGCodes.h:93
double fFermiMomentum
whether or not particle collision is pauli blocked
THE MAIN GENIE PROJECT NAMESPACE
Definition: AlgCmp.h:25
double E(void) const
Get energy.
Definition: GHepParticle.h:91
const int kPdgLambda
Definition: PDGCodes.h:85
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:71
const TLorentzVector * P4(void) const
Definition: GHepParticle.h:78
Defines the EventRecordVisitorI interface. Concrete implementations of this interface use the &#39;Visito...
double fDelRPion
factor by which Pion Compton wavelength gets multiplied to become nuclear size enhancement ...
double fEPreEq
threshold for pre-equilibrium reaction
std::string string
Definition: nybbler.cc:12
virtual int RemnantNucleusPosition(void) const
Definition: GHepRecord.cxx:389
double fNR
param multiplying the nuclear radius, determining how far to track hadrons beyond the "nuclear bounda...
void StepParticle(GHepParticle *p, double step, double nuclear_radius=-1.)
Step particle.
double fFermiFac
testing parameter to modify fermi momentum
#define pFATAL
Definition: Messenger.h:56
int fRemnA
remnant nucleus A
const int kPdgSigma0
Definition: PDGCodes.h:88
Pure abstract base class. Defines the NuclearModelI interface to be implemented by any physics model ...
Definition: NuclearModelI.h:46
bool fXsecNNCorr
use nuclear medium correction for NN cross section
bool IsInNucleus(const GHepParticle *p) const
double Mass(void) const
Mass that corresponds to the PDG code.
static constexpr double MeV
Definition: Units.h:129
bool CanRescatter(const GHepParticle *p) const
double fNucAbsFac
absorption xsec correction factor (hN Mode)
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:29
double Pz(void) const
Get Pz.
Definition: GHepParticle.h:90
static QStrList * l
Definition: config.cpp:1044
GHepStatus_t Status(void) const
Definition: GHepParticle.h:64
GEvGenMode_t fGMode
event generation mode (lepton+A, hadron+A, ...)
const NuclearModelI * fNuclmodel
nuclear model used to generate fermi momentum
double Energy(void) const
Get energy.
Definition: GHepParticle.h:92
double Px(void) const
Get Px.
Definition: GHepParticle.h:88
virtual GHepParticle * Probe(void) const
Definition: GHepRecord.cxx:277
const int kPdgK0
Definition: PDGCodes.h:174
double Vt(void) const
Get production time.
Definition: GHepParticle.h:97
int Pdg(void) const
Definition: GHepParticle.h:63
int FirstMother(void) const
Definition: GHepParticle.h:66
T abs(T value)
void SetPosition(const TLorentzVector &v4)
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
static Config * config
Definition: config.cpp:1054
const int kPdgKM
Definition: PDGCodes.h:173
std::void_t< T > n
virtual GHepParticle * FinalStatePrimaryLepton(void) const
Definition: GHepRecord.cxx:326
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
GEvGenMode_t EventGenerationMode(void) const
Definition: GHepRecord.cxx:209
const int kPdgKP
Definition: PDGCodes.h:172
bool NeedsRescattering(const GHepParticle *p) const
p
Definition: test.py:223
const int kPdgPiP
Definition: PDGCodes.h:158
const int kPdgPi0
Definition: PDGCodes.h:160
double Charge(void) const
Chrg that corresponds to the PDG code.
virtual GHepParticle * TargetNucleus(void) const
Definition: GHepRecord.cxx:286
bool fDoFermi
whether or not to do fermi mom.
void SetTrackingRadius(const GHepParticle *p) const
#define pINFO
Definition: Messenger.h:62
const int kPdgK0L
Definition: PDGCodes.h:176
double fNucCEXFac
charge exchange xsec correction factor (hN Mode)
double MeanFreePath(int pdgc, const TLorentzVector &x4, const TLorentzVector &p4, double A, double Z, double nRpi=0.5, double nRnuc=1.0, const bool useOset=false, const bool altOset=false, const bool xsecNNCorr=false, string INukeMode="XX2018")
Mean free path (pions, nucleons)
Misc GENIE control constants.
void GenerateVertex(GHepRecord *ev) const
double fR0
effective nuclear size param
#define pWARN
Definition: Messenger.h:60
const int kPdgSigmaM
Definition: PDGCodes.h:89
double GenerateStep(GHepRecord *ev, GHepParticle *p) const
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:57
const int kPdgXiM
Definition: PDGCodes.h:94
bool IsIon(int pdgc)
Definition: PDGUtils.cxx:39
double Vz(void) const
Get production z.
Definition: GHepParticle.h:96
const TLorentzVector * X4(void) const
Definition: GHepParticle.h:79
void SetStatus(GHepStatus_t s)
Definition: GHepParticle.h:126
bool fAltOset
NuWro&#39;s table-based implementation (not recommended)
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:65
bool IsPseudoParticle(int pdgc)
Definition: PDGUtils.cxx:24
virtual GHepParticle * HitNucleon(void) const
Definition: GHepRecord.cxx:306
void Configure(string mesg)
Definition: gEvServ.cxx:196
double fDelRNucleon
factor by which Nucleon Compton wavelength gets multiplied to become nuclear size enhancement ...
int getA(int pdg)
double fHadStep
step size for intranuclear hadron transport
virtual GHepParticle * RemnantNucleus(void) const
Definition: GHepRecord.cxx:296
const int kPdgAntiProton
Definition: PDGCodes.h:82
#define A
Definition: memgrp.cpp:38
const int kPdgPiM
Definition: PDGCodes.h:159
const int kPdgSigmaP
Definition: PDGCodes.h:87
virtual void AddParticle(const GHepParticle &p)
Definition: GHepRecord.cxx:491
void TransportHadrons(GHepRecord *ev) const
list x
Definition: train.py:276
const int kPdgProton
Definition: PDGCodes.h:81
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:75
double Vy(void) const
Get production y.
Definition: GHepParticle.h:95
int A(void) const
#define pNOTICE
Definition: Messenger.h:61
bool GetParamDef(const RgKey &name, T &p, const T &def) const
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
double fTrackingRadius
tracking radius for the nucleus in the current event
double fNucRmvE
binding energy to subtract from cascade nucleons
const int kPdgNeutron
Definition: PDGCodes.h:83
GENIE&#39;s GHEP MC event record.
Definition: GHepRecord.h:45
bool fUseOset
Oset model for low energy pion in hN.
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:39
string Vec3AsString(const TVector3 *vec)
Definition: PrintUtils.cxx:80
Root of GENIE utility namespaces.
double Vx(void) const
Get production x.
Definition: GHepParticle.h:94
int fRemnZ
remnant nucleus Z
double Py(void) const
Get Py.
Definition: GHepParticle.h:89
const Algorithm * SubAlg(const RgKey &registry_key) const
Definition: Algorithm.cxx:345