1#include "Framework/Conventions/GBuild.h"
2#ifdef __GENIE_GEANT4_INTERFACE_ENABLED__
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#include "Geant4/G4Diproton.hh"
55#include "Geant4/G4Dineutron.hh"
56#include "Geant4/G4UnboundPN.hh"
59using std::ostringstream;
68HG4BertCascIntranuke::HG4BertCascIntranuke()
75HG4BertCascIntranuke::HG4BertCascIntranuke(
string config)
81HG4BertCascIntranuke::~HG4BertCascIntranuke()
85int HG4BertCascIntranuke::G4BertCascade(
GHepRecord * evrec)
const{
90 const G4ParticleDefinition* incidentDef = PDGtoG4Particle(probe->
Pdg() );
92 int Zinit = tgtNucl->
Z();
93 int Ainit = tgtNucl->
A();
95 G4InuclNuclei *theNucleus =
new G4InuclNuclei(Ainit,Zinit);
98 const G4ThreeVector tgtmomentum(0.,0.,0.);
99 const G4double tgtEnergy =tgtNucl->
Energy()*1000;
100 G4LorentzVector tgt4momentum(tgtmomentum,tgtEnergy);
103 TLorentzVector *p4Genie = probe->
P4();
106 G4ThreeVector incidentDir(p4Genie->Vect().Unit().Px(),p4Genie->Vect().Unit().Py(),p4Genie->Vect().Unit().Pz());
108 double dynamicMass = std::sqrt(p4Genie->M2() );
109 double incidentKE = p4Genie->E() - dynamicMass;
110 const G4DynamicParticle p(incidentDef, incidentDir, incidentKE/
units::MeV, dynamicMass/
units::MeV);
111 G4InuclElementaryParticle* incident =
new G4InuclElementaryParticle(p,G4InuclParticle::bullet);
114 this->SetTrackingRadius(tgtNucl);
115 this->GenerateVertex(evrec);
120 bool has_interacted =
false;
121 while ( this-> IsInNucleus(sp) ) {
126 double d = this->GenerateStep(evrec,sp);
127 has_interacted = (d<fHadStep);
128 if(has_interacted)
break;
131 if ( has_interacted ) {
134 G4CollisionOutput cascadeOutput;
135 G4InuclCollider bertCollider;
136 bertCollider.useCascadeDeexcitation();
138 bertCollider.collide(incident,theNucleus,cascadeOutput);
146 TLorentzVector remX(0.,0.,0.,0.);
147 int Nfrag = cascadeOutput.numberOfOutgoingNuclei();
148 const std::vector<G4InuclNuclei>& outgoingFragments =
149 cascadeOutput.getOutgoingNuclei();
156 int Nhad = cascadeOutput.numberOfOutgoingParticles();
157 const std::vector<G4InuclElementaryParticle>& outgoingHadrons =
158 cascadeOutput.getOutgoingParticles();
159 for (
int l = 0; l < Nhad; l++) {
160 npdg = outgoingHadrons[l].getDefinition()->GetPDGEncoding();
161 G4LorentzVector hadP = outgoingHadrons[l].getMomentum();
162 TLorentzVector tempP(hadP.px(), hadP.py(), hadP.pz(), hadP.e() );
171 for (
int j = 0; j < Nfrag; j++) {
172 if (outgoingFragments[j].
getA() > maxA) {
173 maxA = outgoingFragments[j].getA();
178 fRemnZ = outgoingFragments[rem_index].getZ();
179 fRemnA = outgoingFragments[rem_index].getA();
182 G4LorentzVector nucP = outgoingFragments[rem_index].getMomentum();
183 TLorentzVector remP(nucP.px(), nucP.py(), nucP.pz(), nucP.e() );
184 npdg = outgoingFragments[rem_index].getDefinition()->GetPDGEncoding();
191 <<
"NO Particle with pdg = " << npdg <<
" in PDGLibrary!";
194 1,0,-1,-1, remP, remX);
203 for (G4int k = 0; k < Nfrag; k++) {
204 if (k != rem_index) {
205 npdg = outgoingFragments[k].getDefinition()->GetPDGEncoding();
206 nucP = outgoingFragments[k].getMomentum();
207 TLorentzVector tempP(nucP.px(), nucP.py(), nucP.pz(), nucP.e() );
216 TLorentzVector p4h (0.,0.,probe->
Pz(),probe->
E());
217 TLorentzVector x4null(0.,0.,0.,0.);
218 TLorentzVector p4tgt (0.,0.,0.,tgtNucl->
Mass());
235 fDelRPion, fDelRNucleon,
239 double d = -1.*LL * TMath::Log(rnd->
RndFsi().Rndm());
244void HG4BertCascIntranuke::ProcessEventRecord(
GHepRecord* evrec)
const
249 LOG(
"HG4BertCascIntranuke",
pINFO) <<
" Lepton-nucleus event generation mode ";
253 SetTrackingRadius(nucltgt);
255 LOG(
"HG4BertCascIntranuke",
pINFO) <<
"No nuclear target found - exit ";
261 TransportHadrons(evrec);
263 G4BertCascade(evrec);
265 LOG(
"HG4BertCascIntranuke",
pINFO) <<
" Inappropriate event generation mode - exit ";
269 LOG(
"HG4BertCascIntranuke",
pINFO) <<
"Done with this event";
273void HG4BertCascIntranuke::InitG4Particles()
const
275 G4LeptonConstructor::ConstructParticle();
276 G4MesonConstructor::ConstructParticle();
277 G4BaryonConstructor::ConstructParticle();
280 G4ParticleTable* pTable = G4ParticleTable::GetParticleTable();
281 pTable->SetReadiness(
true);
282 G4GenericIon* gIon = G4GenericIon::Definition();
283 gIon->SetProcessManager(
new G4ProcessManager(gIon) );
286 const G4ParticleDefinition* electron = PDGtoG4Particle(11);
287 const G4ParticleDefinition* proton = PDGtoG4Particle(2212);
288 const G4ParticleDefinition* neutron = PDGtoG4Particle(2112);
289 const G4ParticleDefinition* piplus = PDGtoG4Particle(211);
291 <<
"testing initialization of G4 particles \n"
292 <<
" e 0x" << electron <<
"\n"
293 <<
" p 0x" << proton <<
"\n"
294 <<
" n 0x" << neutron <<
"\n"
295 <<
" pi+ 0x" << piplus <<
"\n"
296 <<
"...InitG4Particles complete";
297 if ( electron == 0 || proton == 0 || neutron == 0 || piplus == 0 ) {
299 <<
"something is seriously wrong with g4 particle lookup";
305void HG4BertCascIntranuke::GenerateVertex(
GHepRecord * evrec)
const
314 TVector3 vtx(999999.,999999.,999999.);
319 double x=999999., y=999999.,
epsilon = 0.001;
320 double R2 = TMath::Power(fTrackingRadius,2.);
321 double rp2 = TMath::Power(x,2.) + TMath::Power(y,2.);
324 y = -fTrackingRadius + 2*fTrackingRadius * rnd->
RndFsi().Rndm();
326 rp2 = TMath::Power(x,2.) + TMath::Power(y,2.);
328 vtx.SetXYZ(x,y, -1.*TMath::Sqrt(TMath::Max(0.,R2-rp2)) +
epsilon);
331 TVector3 direction = evrec->
Probe()->
P4()->Vect().Unit();
334 vtx.RotateUz(direction);
337 <<
"Generated vtx @ R = " << vtx.Mag() <<
" fm / "
340 TObjArrayIter piter(evrec);
351void HG4BertCascIntranuke::SetTrackingRadius(
const GHepParticle * p)
const
355 fTrackingRadius = fR0*TMath::Power(A, 1./3.);
360 fTrackingRadius *= fNR;
363 <<
"Setting tracking radius to R = " << fTrackingRadius;
367bool HG4BertCascIntranuke::IsInNucleus(
const GHepParticle * p)
const
370 return (p->
X4()->Vect().Mag() < fTrackingRadius + fHadStep);
373void HG4BertCascIntranuke::TransportHadrons(
GHepRecord * evrec)
const
389 const G4ParticleDefinition* incidentDef = 0;
391 TObjArrayIter piter(evrec);
392 TObjArrayIter pitter(evrec);
393 TObjArrayIter pittter(evrec);
395 bool has_incidentBaryon(
false), has_secondaries(
false);
396 bool has_remnant(
false), has_incidentparticle(
false);
403 bool has_interacted(
false);
404 if ( ! this->NeedsRescattering(p) )
continue;
409 if ( ! this->CanRescatter(sp) ) {
417 while ( this-> IsInNucleus(sp) ) {
420 double d = this->GenerateStep(evrec,sp);
421 has_interacted = (d<fHadStep);
422 if ( has_interacted ) {
423 has_secondaries=
true;
427 if ( ! has_interacted ) {
434 if ( IsBaryon(sp) ) {
436 incidentDef = PDGtoG4Particle(sp->
Pdg() );
437 has_incidentBaryon=
true;
445 LOG(
"G4BertCascInterface::TransportHadrons",
pWARN)
446 <<
"IsBaryon failed to tag " << *sp;
454 std::vector<int> Postion_evrec;
455 if ( has_secondaries ) {
456 if ( ! incidentBaryon )
LOG(
"G4BertCascInterface::TransportHadrons",
pINFO)
457 <<
"Unrecognized baryon in nucleus";
459 delete incidentBaryon;
460 G4Fancy3DNucleus* g4Nucleus =
new G4Fancy3DNucleus();
462 TLorentzVector pIncident;
464 g4Nucleus->Init(remNucl->
A(),remNucl->
Z());
465 double EE = struckNucleon->
E() - tgtNucl->
Mass() + g4Nucleus->GetMass()*
units::MeV;
466 TLorentzVector struckMomentum(struckNucleon->
Px(), struckNucleon->
Py(), struckNucleon->
Pz(), EE);
467 Double_t PxI(0),PyI(0),PzI(0),EEI(0), Q(0), P(0), N(0);
478 if ( pos_in_evrec==0 ) pos_in_evrec = icccur;
479 Postion_evrec.push_back(icccur);
482 if ( pos_in_evrec==0 ) pos_in_evrec = icccur;
483 if ( ! has_incidentparticle ) {
485 incidentDef = PDGtoG4Particle(p->
Pdg() );
486 has_incidentparticle=
true;
491 if ( ! has_incidentparticle) {
493 incidentDef=PDGtoG4Particle(pinN->
Pdg());
496 else if (N == 2) incidentDef = PDGtoG4Particle(
kPdgClusterNN);
497 else if (P == 1 && N == 1) incidentDef = PDGtoG4Particle(
kPdgClusterNP);
500 pIncident.SetPxPyPzE(PxI,PyI,PzI,EEI);
503 G4ThreeVector incidentDir(pIncident.Vect().Unit().Px(),
504 pIncident.Vect().Unit().Py(),
505 pIncident.Vect().Unit().Pz());
507 double dynamicMass = std::sqrt(pIncident.M2() );
508 double incidentKE = pIncident.E() - dynamicMass;
514 G4InuclElementaryParticle* incident =
new G4InuclElementaryParticle(dp,G4InuclParticle::bullet);
518 G4KineticTrackVector* g4secondaries = ConvertGenieSecondariesToG4(evrec);
520 Nsec = g4secondaries->size();
523 G4CollisionOutput cascadeOutput;
524 G4InuclCollider bertCollider;
525 bertCollider.useCascadeDeexcitation();
527 bertCollider.rescatter(incident, g4secondaries, g4Nucleus, cascadeOutput);
530 for (
int n = 0; n < Nsec; n++)
delete (*g4secondaries)[n];
531 delete g4secondaries;
536 TLorentzVector remX(tgtNucl->
Vx(), tgtNucl->
Vy(), tgtNucl->
Vz(), tgtNucl->
Vt() );
539 int Nfrag = cascadeOutput.numberOfOutgoingNuclei();
540 const std::vector<G4InuclNuclei>& outgoingFragments = cascadeOutput.getOutgoingNuclei();
547 int Nhad = cascadeOutput.numberOfOutgoingParticles();
548 const std::vector<G4InuclElementaryParticle>& outgoingHadrons = cascadeOutput.getOutgoingParticles();
550 int mother1=Postion_evrec.at(0);
552 if(Nsec==1)mother2=-1;
553 else if(Nsec>1)mother2=Postion_evrec.at(Nsec-1);
554 for (
int l = 0; l < Nhad; l++) {
555 npdg = outgoingHadrons[l].getDefinition()->GetPDGEncoding();
557 G4LorentzVector hadP = outgoingHadrons[l].getMomentum();
558 TLorentzVector tempP(hadP.px(), hadP.py(), hadP.pz(), hadP.e() );
567 for (
int j = 0; j < Nfrag; j++) {
568 if (outgoingFragments[j].
getA() > maxA) {
569 maxA = outgoingFragments[j].getA();
574 fRemnZ = outgoingFragments[rem_index].getZ();
575 fRemnA = outgoingFragments[rem_index].getA();
578 G4LorentzVector nucP = outgoingFragments[rem_index].getMomentum();
579 TLorentzVector remP(nucP.px(), nucP.py(), nucP.pz(), nucP.e() );
582 for (G4int k = 0; k < Nfrag; k++) {
583 if (k != rem_index) {
584 npdg = outgoingFragments[k].getDefinition()->GetPDGEncoding();
585 nucP = outgoingFragments[k].getMomentum();
586 TLorentzVector tempP(nucP.px(), nucP.py(), nucP.pz(), nucP.e() );
594 npdg = outgoingFragments[rem_index].getDefinition()->GetPDGEncoding();
595 remP.SetPx(remP.Px()+remNucl->
P4()->Px());
596 remP.SetPy(remP.Py()+remNucl->
P4()->Py());
597 remP.SetPz(remP.Pz()+remNucl->
P4()->Pz());
604 <<
"NO Particle with pdg = " << npdg <<
" in PDGLibrary!";
607 rem_nucl,-1,-1,-1, remP, remX);
628 int dau1(0), dau2(0);
633 for(
int ii=1;ii<Nsec;ii++){
651bool HG4BertCascIntranuke::NeedsRescattering(
const GHepParticle * p)
const {
659bool HG4BertCascIntranuke::CanRescatter(
const GHepParticle * p)
const
683bool HG4BertCascIntranuke::IsBaryon(
const GHepParticle * p)
const
691 <<
"no entry for PDG " << p->
Pdg() <<
" in PDGLibrary";
693 if ( std::string(ppdg->ParticleClass()) == std::string(
"Baryon") ) {
700const G4ParticleDefinition* HG4BertCascIntranuke::PDGtoG4Particle(
int pdg)
const
702 const G4ParticleDefinition* pDef = 0;
708 if ( abs(
pdg) < 1000000000 ) {
709 pDef = G4ParticleTable::GetParticleTable()->FindParticle(
pdg);
710 }
else if (
pdg < 2000000000 ) {
711 pDef = G4IonTable::GetIonTable()->GetIon(
pdg);
716 <<
"Unrecognized Bertini particle type: " <<
pdg;
723G4KineticTrackVector* HG4BertCascIntranuke::ConvertGenieSecondariesToG4(std::vector<GHepParticle> partList)
const
725 static double GeToG4length = 2.81967*1.0e-12*1.2/1.4;
728 const G4ParticleDefinition* pDef = 0;
729 G4KineticTrackVector* g4secondaries =
new G4KineticTrackVector;
730 G4KineticTrack* kt = 0;
732 for (
size_t it=0 ; it<partList.size(); it++) {
734 pDef = PDGtoG4Particle(p->
Pdg() );
735 double formationTime = p->
Vt();
736 G4ThreeVector formationPosition(p->
Vx()*GeToG4length,
737 p->
Vy()*GeToG4length,
738 p->
Vz()*GeToG4length);
741 kt =
new G4KineticTrack(pDef, formationTime, formationPosition, formationMomentum);
742 kt->SetDefinition(pDef);
743 kt->SetState(G4KineticTrack::inside);
744 g4secondaries->push_back(kt);
747 return g4secondaries;
751G4KineticTrackVector* HG4BertCascIntranuke::ConvertGenieSecondariesToG4(
GHepRecord* evrec)
const {
757 static double GeToG4length = 2.81967*1.0e-12*1.2/1.4;
759 TObjArrayIter piter(evrec);
761 const G4ParticleDefinition* pDef = 0;
762 G4KineticTrackVector* g4secondaries =
new G4KineticTrackVector;
763 G4KineticTrack* kt = 0;
768 pDef = PDGtoG4Particle(p->
Pdg() );
769 double formationTime = p->
Vt();
770 G4ThreeVector formationPosition(p->
Vx()*GeToG4length,
771 p->
Vy()*GeToG4length,
772 p->
Vz()*GeToG4length);
775 kt =
new G4KineticTrack(pDef, formationTime,
776 formationPosition, formationMomentum);
777 kt->SetDefinition(pDef);
778 kt->SetState(G4KineticTrack::inside);
779 g4secondaries->push_back(kt);
783 return g4secondaries;
787bool HG4BertCascIntranuke::Conserve4Momentum(
GHepRecord* evrec)
const
795 int Nentries = evrec->GetEntries();
799 TLorentzVector sum4mom(0.0, 0.0, 0.0, 0.0);
800 for (
int j = 0; j < Nentries; j++) {
803 sum4mom += *(p->
P4() );
805 <<
" Final state 4-momentum = ("
806 << p->
P4()->Px() <<
", " << p->
P4()->Py() <<
", "
807 << p->
P4()->Pz() <<
", " << p->
P4()->E() <<
") ";
810 sum4mom += *(finalLepton->
P4() );
812 TLorentzVector initial4mom = *(targetNucleus->
P4() ) + *(probe->
P4() );
814 TLorentzVector diff = initial4mom - sum4mom;
816 const double maxdiff = 1.0e-9;
817 double diffE = diff.E();
818 TVector3 diffp3 = diff.Vect();
819 double diffpmag = diffp3.Mag();
820 if ( TMath::Abs(diffE) > maxdiff || diffpmag > maxdiff ) {
823 <<
"final state - initial state differs by > " << maxdiff <<
"\n"
824 <<
" dE = " << diffE <<
", d|p3| = " << diffpmag;
827 <<
" Total Final state 4-momentum = (" << sum4mom.Px()
828 <<
", " << sum4mom.Py()
829 <<
", " << sum4mom.Pz()
830 <<
", " << sum4mom.E() <<
") ";
832 <<
" Total Initial state 4-momentum = (" << initial4mom.Px()
833 <<
", " << initial4mom.Py()
834 <<
", " << initial4mom.Pz()
835 <<
", " << initial4mom.E() <<
") ";
839 double Qinit = targetNucleus->
Charge();
840 double Qfinal = finalLepton->
Charge();
842 for (
int j = 0; j < Nentries; j++) {
849 if (Qinit != Qfinal) {
852 <<
" Charge not conserved: \n"
853 <<
" Qfinal = " << Qfinal
854 <<
" Qinit = " << Qinit;
862void HG4BertCascIntranuke::LoadConfig(
void)
866 fNuclmodel =
dynamic_cast<const NuclearModelI *
>(
this -> SubAlg(
"NuclearModel") ) ;
869 GetParam(
"NUCL-R0", fR0 );
870 GetParam(
"NUCL-NR", fNR );
872 GetParam(
"INUKE-NucRemovalE", fNucRmvE );
873 GetParam(
"INUKE-HadStep", fHadStep ) ;
874 GetParam(
"INUKE-NucAbsFac", fNucAbsFac ) ;
875 GetParam(
"INUKE-NucCEXFac", fNucCEXFac ) ;
876 GetParam(
"INUKE-Energy_Pre_Eq", fEPreEq ) ;
877 GetParam(
"INUKE-FermiFac", fFermiFac ) ;
878 GetParam(
"INUKE-FermiMomentum", fFermiMomentum ) ;
879 GetParam(
"INUKE-DoFermi", fDoFermi ) ;
880 GetParam(
"INUKE-XsecNNCorr", fXsecNNCorr ) ;
881 GetParamDef(
"UseOset", fUseOset,
false ) ;
882 GetParamDef(
"AltOset", fAltOset,
false ) ;
883 GetParam(
"HNINUKE-DelRPion", fDelRPion ) ;
884 GetParam(
"HNINUKE-DelRNucleon", fDelRNucleon ) ;
889 <<
"R0 = " << fR0 <<
" fermi" <<
"\n"
890 <<
"NR = " << fNR <<
"\n"
891 <<
"DelRPion = " << fDelRPion <<
"\n"
892 <<
"DelRNucleon = " << fDelRNucleon <<
"\n"
893 <<
"HadStep = " << fHadStep <<
" fermi" <<
"\n"
894 <<
"EPreEq = " << fHadStep <<
" fermi" <<
"\n"
895 <<
"NucAbsFac = " << fNucAbsFac <<
"\n"
896 <<
"NucCEXFac = " << fNucCEXFac <<
"\n"
897 <<
"FermiFac = " << fFermiFac <<
"\n"
898 <<
"FermiMomtm = " << fFermiMomentum <<
"\n"
899 <<
"DoFermi? = " << ((fDoFermi)?(
true):(
false)) <<
"\n"
900 <<
"XsecNNCorr? = " << ((fXsecNNCorr)?(
true):(
false));
913void HG4BertCascIntranuke::Configure(
string param_set)
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils.
virtual void Configure(const Registry &config)
Defines the EventRecordVisitorI interface. Concrete implementations of this interface use the 'Visito...
STDHEP-like event record entry that can fit a particle or a nucleus.
int FirstMother(void) const
void SetPosition(const TLorentzVector &v4)
void SetLastDaughter(int d)
void SetRescatterCode(int code)
void SetFirstMother(int m)
double Vy(void) const
Get production y.
const TLorentzVector * P4(void) const
double Charge(void) const
Chrg that corresponds to the PDG code.
double Mass(void) const
Mass that corresponds to the PDG code.
const TLorentzVector * X4(void) const
int LastDaughter(void) const
void SetStatus(GHepStatus_t s)
double Px(void) const
Get Px.
double E(void) const
Get energy.
double Pz(void) const
Get Pz.
double Py(void) const
Get Py.
double Vz(void) const
Get production z.
double Energy(void) const
Get energy.
int RescatterCode(void) const
void SetFirstDaughter(int d)
GHepStatus_t Status(void) const
double Vx(void) const
Get production x.
double Vt(void) const
Get production time.
int FirstDaughter(void) const
GENIE's GHEP MC event record.
virtual GHepParticle * Probe(void) const
virtual GHepParticle * TargetNucleus(void) const
GEvGenMode_t EventGenerationMode(void) const
virtual GHepParticle * RemnantNucleus(void) const
virtual void AddParticle(const GHepParticle &p)
virtual int RemnantNucleusPosition(void) const
virtual GHepParticle * Particle(int position) const
virtual GHepParticle * FinalStatePrimaryLepton(void) const
virtual GHepParticle * HitNucleon(void) const
static string AsString(INukeMode_t mode)
Pure abstract base class. Defines the NuclearModelI interface to be implemented by any physics model ...
static PDGLibrary * Instance(void)
A singleton holding random number generator classes. All random number generation in GENIE should tak...
static RandomGen * Instance()
Access instance.
TRandom3 & RndFsi(void) const
rnd number generator used by intranuclear cascade monte carlos
A registry. Provides the container for algorithm configuration parameters.
Misc GENIE control constants.
Utilities for improving the code readability when using PDG codes.
bool IsPseudoParticle(int pdgc)
static constexpr double MeV
Simple functions for loading and reading nucleus dependent keys from config files.
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)
void StepParticle(GHepParticle *p, double step, double nuclear_radius=-1.)
Step particle.
string Vec3AsString(const TVector3 *vec)
Root of GENIE utility namespaces.
THE MAIN GENIE PROJECT NAMESPACE
@ kIStFinalStateNuclearRemnant
const int kPdgHadronicBlob