Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4ParaFissionModel Class Reference

#include <G4ParaFissionModel.hh>

Inheritance diagram for G4ParaFissionModel:
G4HadronicInteraction

Public Member Functions

 G4ParaFissionModel ()
 
virtual ~G4ParaFissionModel ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
 
- 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, G4double
GetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 

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 44 of file G4ParaFissionModel.hh.

Constructor & Destructor Documentation

G4ParaFissionModel::G4ParaFissionModel ( )
inline

Definition at line 48 of file G4ParaFissionModel.hh.

References python.hepunit::MeV, G4HadronicInteraction::SetMaxEnergy(), and G4HadronicInteraction::SetMinEnergy().

49  {
50  SetMinEnergy( 0.0 );
51  SetMaxEnergy( 60.*MeV );
52  }
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
virtual G4ParaFissionModel::~G4ParaFissionModel ( )
inlinevirtual

Definition at line 54 of file G4ParaFissionModel.hh.

54 {};

Member Function Documentation

virtual G4HadFinalState* G4ParaFissionModel::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus theNucleus 
)
inlinevirtual

Implements G4HadronicInteraction.

Definition at line 56 of file G4ParaFissionModel.hh.

References G4HadFinalState::AddSecondary(), G4ExcitationHandler::BreakItUp(), G4CompetitiveFission::BreakUp(), G4HadFinalState::Clear(), G4HadProjectile::Get4Momentum(), G4Nucleus::GetA_asInt(), G4ParticleDefinition::GetBaryonNumber(), G4HadProjectile::GetDefinition(), G4ReactionProduct::GetDefinition(), G4Fragment::GetExcitationEnergy(), G4Fragment::GetMomentum(), G4ReactionProduct::GetMomentum(), G4NucleiProperties::GetNuclearMass(), G4Fragment::GetParticleDefinition(), G4ParticleDefinition::GetPDGCharge(), G4Nucleus::GetZ_asInt(), python.hepunit::keV, G4HadFinalState::SetEnergyChange(), G4Fragment::SetNumberOfExcitedParticle(), G4Fragment::SetNumberOfHoles(), G4HadFinalState::SetStatusChange(), stopAndKill, and test::v.

58  {
59  theParticleChange.Clear();
60  theParticleChange.SetStatusChange( stopAndKill );
61  theParticleChange.SetEnergyChange( 0.0 );
62 
63  // prepare the fragment
64 
65  G4int A = theNucleus.GetA_asInt();
66  G4int Z = theNucleus.GetZ_asInt();
68 
69  G4int numberOfEx = aTrack.GetDefinition()->GetBaryonNumber();
70  G4int numberOfCh = G4int(aTrack.GetDefinition()->GetPDGCharge() + 0.5);
71  G4int numberOfHoles = 0;
72 
73  A += numberOfEx;
74  Z += numberOfCh;
75 
76  G4LorentzVector v = aTrack.Get4Momentum() + G4LorentzVector(0.0,0.0,0.0,nucMass);
77  G4Fragment anInitialState(A,Z,v);
78  anInitialState.SetNumberOfExcitedParticle(numberOfEx,numberOfCh);
79  anInitialState.SetNumberOfHoles(0,0);
80 
81  // do the fission
82  G4FragmentVector * theFissionResult = theFission.BreakUp(anInitialState);
83 
84  // deexcite the fission fragments and fill result
85 
86  G4int ll = theFissionResult->size();
87  for(G4int i=0; i<ll; i++)
88  {
89  G4ReactionProductVector* theExcitationResult = 0;
90  G4Fragment* aFragment = (*theFissionResult)[i];
91  if(aFragment->GetExcitationEnergy() > keV)
92  {
93  theExcitationResult = theHandler.BreakItUp(*aFragment);
94 
95  // add secondaries
96  for(G4int j = 0; j < G4int(theExcitationResult->size()); j++)
97  {
98  G4ReactionProduct* rp0 = (*theExcitationResult)[j];
99  G4DynamicParticle* p0 =
100  new G4DynamicParticle(rp0->GetDefinition(),rp0->GetMomentum());
101  theParticleChange.AddSecondary(p0);
102  delete rp0;
103  }
104  delete theExcitationResult;
105  }
106  else
107  {
108  // add secondary
109  G4DynamicParticle* p0 =
110  new G4DynamicParticle(aFragment->GetParticleDefinition(),
111  aFragment->GetMomentum());
112  theParticleChange.AddSecondary(p0);
113  }
114  delete aFragment;
115  }
116 
117  delete theFissionResult;
118 
119  return &theParticleChange;
120  }
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
static G4double GetNuclearMass(const G4double A, const G4double Z)
int G4int
Definition: G4Types.hh:78
virtual G4FragmentVector * BreakUp(const G4Fragment &theNucleus)
G4ReactionProductVector * BreakItUp(const G4Fragment &theInitialState) const
G4ParticleDefinition * GetDefinition() const
void SetStatusChange(G4HadFinalStateStatus aS)
std::vector< G4ReactionProduct * > G4ReactionProductVector
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:271
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
const G4LorentzVector & Get4Momentum() const
void SetEnergyChange(G4double anEnergy)
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:388
G4ThreeVector GetMomentum() const
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
void AddSecondary(G4DynamicParticle *aP)
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:255
CLHEP::HepLorentzVector G4LorentzVector

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