Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ParaFissionModel.hh
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 #ifndef G4ParaFissionModel_h
27 #define G4ParaFissionModel_h 1
28 
29 #include "G4CompetitiveFission.hh"
30 #include "G4ExcitationHandler.hh"
31 #include "G4HadronicInteraction.hh"
32 #include "G4NucleiProperties.hh"
33 //#include "G4ParticleTable.hh"
34 
35 // Class Description
36 // Final state production model for (based on evaluated data
37 // libraries) description of neutron induced fission below 60 MeV;
38 // In case you need the fission fragments, use this model.
39 // To be used in your physics list in case you need this physics.
40 // In this case you want to register an object of this class with
41 // the corresponding process.
42 
43 
45 {
46 public:
47 
49  {
50  SetMinEnergy( 0.0 );
51  SetMaxEnergy( 60.*MeV );
52  }
53 
54  virtual ~G4ParaFissionModel() {};
55 
57  G4Nucleus& theNucleus)
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  }
121 private:
122 
123  G4CompetitiveFission theFission;
124  G4ExcitationHandler theHandler;
125 
126  G4HadFinalState theParticleChange;
127 };
128 #endif
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
static G4double GetNuclearMass(const G4double A, const G4double Z)
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:438
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:375
int G4int
Definition: G4Types.hh:78
G4ReactionProductVector * BreakItUp(const G4Fragment &theInitialState)
void SetStatusChange(G4HadFinalStateStatus aS)
std::vector< G4ReactionProduct * > G4ReactionProductVector
void SetMinEnergy(G4double anEnergy)
void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
Definition: G4Fragment.hh:356
const G4ParticleDefinition * GetDefinition() const
double A(double temperature)
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:307
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:63
const G4LorentzVector & Get4Momentum() const
tuple v
Definition: test.py:18
void SetEnergyChange(G4double anEnergy)
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
void SetMaxEnergy(const G4double anEnergy)
G4ThreeVector GetMomentum() const
static constexpr double MeV
Definition: G4SIunits.hh:214
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
static constexpr double keV
Definition: G4SIunits.hh:216
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:283
CLHEP::HepLorentzVector G4LorentzVector