Geant4  10.02.p03
G4AblaInterface Class Reference

#include <G4AblaInterface.hh>

Inheritance diagram for G4AblaInterface:
Collaboration diagram for G4AblaInterface:

Public Member Functions

 G4AblaInterface ()
 
virtual ~G4AblaInterface ()
 
virtual G4ReactionProductVectorDeExcite (G4Fragment &aFragment)
 
virtual G4HadFinalStateApplyYourself (G4HadProjectile const &, G4Nucleus &)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void DeExciteModelDescription (std::ostream &outFile) const
 
- Public Member Functions inherited from G4VPreCompoundModel
 G4VPreCompoundModel (G4ExcitationHandler *ptr=0, const G4String &modelName="PrecompoundModel")
 
virtual ~G4VPreCompoundModel ()
 
void SetExcitationHandler (G4ExcitationHandler *ptr)
 
G4ExcitationHandlerGetExcitationHandler () const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
const G4HadronicInteractionGetMyPointer () const
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 

Private Member Functions

G4ReactionProducttoG4Particle (G4int A, G4int Z, G4double kinE, G4double px, G4double py, G4double pz) const
 Convert an Abla particle to a G4ReactionProduct. More...
 
G4ParticleDefinitiontoG4ParticleDefinition (G4int A, G4int Z) const
 Convert A and Z to a G4ParticleDefinition. More...
 

Private Attributes

G4VarNtpablaResult
 
G4Volantvolant
 
G4AblatheABLAModel
 
G4long eventNumber
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 49 of file G4AblaInterface.hh.

Constructor & Destructor Documentation

◆ G4AblaInterface()

G4AblaInterface::G4AblaInterface ( )

Definition at line 49 of file G4AblaInterface.cc.

49  :
50  G4VPreCompoundModel(NULL, "ABLA"),
51  ablaResult(new G4VarNtp),
52  volant(new G4Volant),
54  eventNumber(0)
55 {
57 }
Definition: G4Abla.hh:54
G4VarNtp * ablaResult
void initEvapora()
Definition: G4Abla.cc:682
G4VPreCompoundModel(G4ExcitationHandler *ptr=0, const G4String &modelName="PrecompoundModel")
Here is the call graph for this function:

◆ ~G4AblaInterface()

G4AblaInterface::~G4AblaInterface ( )
virtual

Definition at line 59 of file G4AblaInterface.cc.

59  {
60  delete volant;
61  delete ablaResult;
62  delete theABLAModel;
63 }
G4VarNtp * ablaResult

Member Function Documentation

◆ ApplyYourself()

virtual G4HadFinalState* G4AblaInterface::ApplyYourself ( G4HadProjectile const &  ,
G4Nucleus  
)
inlinevirtual

Implements G4VPreCompoundModel.

Definition at line 56 of file G4AblaInterface.hh.

56  {
57  return NULL;
58  }
Here is the call graph for this function:

◆ DeExcite()

G4ReactionProductVector * G4AblaInterface::DeExcite ( G4Fragment aFragment)
virtual

Implements G4VPreCompoundModel.

Definition at line 65 of file G4AblaInterface.cc.

65  {
66  volant->clear();
67  ablaResult->clear();
68 
69  const G4int ARem = aFragment.GetA_asInt();
70  const G4int ZRem = aFragment.GetZ_asInt();
71  const G4double nuclearMass = aFragment.GetGroundStateMass() / MeV;
72  const G4double eStarRem = aFragment.GetExcitationEnergy() / MeV;
73  const G4double jRem = aFragment.GetAngularMomentum().mag() / hbar_Planck;
74  const G4LorentzVector &pRem = aFragment.GetMomentum();
75  const G4double eTotRem = pRem.e();
76  const G4double eKinRem = (eTotRem - pRem.invariantMass()) / MeV;
77  const G4double pxRem = pRem.x() / MeV;
78  const G4double pyRem = pRem.y() / MeV;
79  const G4double pzRem = pRem.z() / MeV;
80 
81  eventNumber++;
82 
83  theABLAModel->breakItUp(ARem, ZRem,
84  nuclearMass,
85  eStarRem,
86  jRem,
87  eKinRem,
88  pxRem,
89  pyRem,
90  pzRem,
91  eventNumber);
92 
94 
95  for(int j = 0; j < ablaResult->ntrack; ++j) { // Copy ABLA result to the EventInfo
97  ablaResult->zvv[j],
98  ablaResult->enerj[j],
99  ablaResult->plab[j]*std::sin(ablaResult->tetlab[j]*pi/180.0)*std::cos(ablaResult->philab[j]*pi/180.0),
100  ablaResult->plab[j]*std::sin(ablaResult->tetlab[j]*pi/180.0)*std::sin(ablaResult->philab[j]*pi/180.0),
101  ablaResult->plab[j]*std::cos(ablaResult->tetlab[j]*pi/180.0));
102  if(product)
103  result->push_back(product);
104  }
105  return result;
106 }
void clear()
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:273
static const double MeV
Definition: G4SIunits.hh:211
G4int GetA_asInt() const
Definition: G4Fragment.hh:256
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:278
G4double plab[VARNTPSIZE]
double invariantMass(const HepLorentzVector &w) const
G4int avv[VARNTPSIZE]
int G4int
Definition: G4Types.hh:78
G4double enerj[VARNTPSIZE]
std::vector< G4ReactionProduct * > G4ReactionProductVector
void breakItUp(G4int nucleusA, G4int nucleusZ, G4double nucleusMass, G4double excitationEnergy, G4double angularMomentum, G4double recoilEnergy, G4double momX, G4double momY, G4double momZ, G4int eventnumber)
Definition: G4Abla.cc:102
G4VarNtp * ablaResult
G4int GetZ_asInt() const
Definition: G4Fragment.hh:261
double mag() const
G4double tetlab[VARNTPSIZE]
static const double pi
Definition: G4SIunits.hh:74
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:289
float hbar_Planck
Definition: hepunit.py:264
G4ReactionProduct * toG4Particle(G4int A, G4int Z, G4double kinE, G4double px, G4double py, G4double pz) const
Convert an Abla particle to a G4ReactionProduct.
G4ThreeVector GetAngularMomentum() const
Definition: G4Fragment.cc:275
void clear()
G4double philab[VARNTPSIZE]
G4int zvv[VARNTPSIZE]
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ DeExciteModelDescription()

void G4AblaInterface::DeExciteModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VPreCompoundModel.

Definition at line 147 of file G4AblaInterface.cc.

147  {
148  outFile << "ABLA V3 is a statistical model for nuclear de-excitation. It simulates\n"
149  << "evaporation of neutrons, protons and alpha particles, as well as fission\n"
150  << "where applicable. The code included in Geant4 is a C++ translation of the\n"
151  << "original Fortran code. More details about the physics are available in the\n"
152  << "the Geant4 Physics Reference Manual and in the reference articles.\n\n"
153  << "References: A.R. Junghans et al., Nucl. Phys. A629 (1998) 635;\n"
154  << " J. Benlliure et al., Nucl. Phys. A628 (1998) 458.\n\n"; }
Here is the caller graph for this function:

◆ ModelDescription()

void G4AblaInterface::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 143 of file G4AblaInterface.cc.

143  {
144  outFile << "ABLA V3 does not provide an implementation of the ApplyYourself method!\n\n";
145 }
Here is the caller graph for this function:

◆ toG4Particle()

G4ReactionProduct * G4AblaInterface::toG4Particle ( G4int  A,
G4int  Z,
G4double  kinE,
G4double  px,
G4double  py,
G4double  pz 
) const
private

Convert an Abla particle to a G4ReactionProduct.

Definition at line 126 of file G4AblaInterface.cc.

129  {
131  if(def == 0) { // Check if we have a valid particle definition
132  return 0;
133  }
134  const double energy = kinE * MeV;
135  const G4ThreeVector momentum(px, py, pz);
136  const G4ThreeVector momentumDirection = momentum.unit();
137  G4DynamicParticle p(def, momentumDirection, energy);
139  (*r) = p;
140  return r;
141 }
static const double MeV
Definition: G4SIunits.hh:211
G4ParticleDefinition * toG4ParticleDefinition(G4int A, G4int Z) const
Convert A and Z to a G4ParticleDefinition.
double energy
Definition: plottest35.C:25
double A(double temperature)
Float_t Z
Hep3Vector unit() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ toG4ParticleDefinition()

G4ParticleDefinition * G4AblaInterface::toG4ParticleDefinition ( G4int  A,
G4int  Z 
) const
private

Convert A and Z to a G4ParticleDefinition.

Definition at line 108 of file G4AblaInterface.cc.

108  {
109  if (A == 1 && Z == 1) return G4Proton::Proton();
110  else if(A == 1 && Z == 0) return G4Neutron::Neutron();
111  else if(A == 0 && Z == 1) return G4PionPlus::PionPlus();
112  else if(A == 0 && Z == -1) return G4PionMinus::PionMinus();
113  else if(A == 0 && Z == 0) return G4PionZero::PionZero();
114  else if(A == 2 && Z == 1) return G4Deuteron::Deuteron();
115  else if(A == 3 && Z == 1) return G4Triton::Triton();
116  else if(A == 3 && Z == 2) return G4He3::He3();
117  else if(A == 4 && Z == 2) return G4Alpha::Alpha();
118  else if(A > 0 && Z > 0 && A >= Z) { // Returns ground state ion definition
119  return G4IonTable::GetIonTable()->GetIon(Z, A, 0);
120  } else { // Error, unrecognized particle
121  G4cout << "Can't convert particle with A=" << A << ", Z=" << Z << " to G4ParticleDefinition, trouble ahead" << G4endl;
122  return 0;
123  }
124 }
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
Float_t Z
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4PionZero * PionZero()
Definition: G4PionZero.cc:108
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
#define G4endl
Definition: G4ios.hh:61
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4He3 * He3()
Definition: G4He3.cc:94
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ ablaResult

G4VarNtp* G4AblaInterface::ablaResult
private

Definition at line 64 of file G4AblaInterface.hh.

◆ eventNumber

G4long G4AblaInterface::eventNumber
private

Definition at line 67 of file G4AblaInterface.hh.

◆ theABLAModel

G4Abla* G4AblaInterface::theABLAModel
private

Definition at line 66 of file G4AblaInterface.hh.

◆ volant

G4Volant* G4AblaInterface::volant
private

Definition at line 65 of file G4AblaInterface.hh.


The documentation for this class was generated from the following files: