Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4GeneratorPrecompoundInterface.cc
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// $Id$
27//
28// -----------------------------------------------------------------------------
29// GEANT 4 class file
30//
31// History: first implementation
32// HPW, 10DEC 98, the decay part originally written by Gunter Folger
33// in his FTF-test-program.
34//
35// M.Kelsey, 28 Jul 2011 -- Replace loop to decay input secondaries
36// with new utility class, simplify cleanup loops
37// -----------------------------------------------------------------------------
38
39#include <algorithm>
40#include <vector>
41
44#include "G4SystemOfUnits.hh"
47#include "G4Proton.hh"
48#include "G4Neutron.hh"
49#include "G4V3DNucleus.hh"
50#include "G4Nucleon.hh"
51#include "G4FragmentVector.hh"
52#include "G4ReactionProduct.hh"
54#include "G4PreCompoundModel.hh"
58
60: CaptureThreshold(10*MeV)
61{
62 proton = G4Proton::Proton();
63 neutron = G4Neutron::Neutron();
64 if(preModel) { SetDeExcitation(preModel); }
65 else {
68 G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(hadi);
69 if(!pre) { pre = new G4PreCompoundModel(); }
70 SetDeExcitation(pre);
71 }
72}
73
75{
76}
77
78//---------------------------------------------------------------------
79// choose to calculate excitation energy from energy balance
80#define exactExcitationEnergy
81
82//#define G4GPI_debug_excitation
83
85Propagate(G4KineticTrackVector* theSecondaries, G4V3DNucleus* theNucleus)
86{
88
89 // decay the strong resonances
90 G4DecayKineticTracks decay(theSecondaries);
91
92 // prepare the fragment
93 G4int anA=theNucleus->GetMassNumber();
94 G4int aZ=theNucleus->GetCharge();
95 G4int numberOfEx = 0;
96 G4int numberOfCh = 0;
97 G4int numberOfHoles = 0;
98 G4double exEnergy = 0.0;
99 G4double R = theNucleus->GetNuclearRadius();
100 G4ThreeVector exciton3Momentum(0.,0.,0.);
101 G4ThreeVector captured3Momentum(0.,0.,0.);
102 G4ThreeVector wounded3Momentum(0.,0.,0.);
103
104 // loop over secondaries
105#ifdef exactExcitationEnergy
106 G4LorentzVector secondary4Momemtum(0,0,0,0);
107#endif
108 G4KineticTrackVector::iterator iter;
109 for(iter=theSecondaries->begin(); iter !=theSecondaries->end(); ++iter)
110 {
111 G4ParticleDefinition* part = (*iter)->GetDefinition();
112 G4double e = (*iter)->Get4Momentum().e();
113 G4double mass = (*iter)->Get4Momentum().mag();
114 G4ThreeVector mom = (*iter)->Get4Momentum().vect();
115 if((part != proton && part != neutron) ||
116 (e > mass + CaptureThreshold) ||
117 ((*iter)->GetPosition().mag() > R)) {
118 G4ReactionProduct * theNew = new G4ReactionProduct(part);
119 theNew->SetMomentum(mom);
120 theNew->SetTotalEnergy(e);
121 theTotalResult->push_back(theNew);
122#ifdef exactExcitationEnergy
123 secondary4Momemtum += (*iter)->Get4Momentum();
124#endif
125 } else {
126 // within the nucleus, neutron or proton
127 // now calculate A, Z of the fragment, momentum, number of exciton states
128 ++anA;
129 ++numberOfEx;
130 G4int Z = G4int(part->GetPDGCharge()/eplus + 0.1);
131 aZ += Z;
132 numberOfCh += Z;
133 captured3Momentum += mom;
134 exEnergy += (e - mass);
135 }
136 delete (*iter);
137 }
138 delete theSecondaries;
139
140 // loop over wounded nucleus
141 G4Nucleon * theCurrentNucleon =
142 theNucleus->StartLoop() ? theNucleus->GetNextNucleon() : 0;
143 while(theCurrentNucleon) {
144 if(theCurrentNucleon->AreYouHit()) {
145 ++numberOfHoles;
146 ++numberOfEx;
147 --anA;
148 aZ -= G4int(theCurrentNucleon->GetDefinition()->GetPDGCharge()/eplus + 0.1);
149 wounded3Momentum += theCurrentNucleon->Get4Momentum().vect();
150 //G4cout << "hit nucleon " << theCurrentNucleon->Get4Momentum() << G4endl;
151 exEnergy += theCurrentNucleon->GetBindingEnergy();
152 }
153 theCurrentNucleon = theNucleus->GetNextNucleon();
154 }
155 exciton3Momentum = captured3Momentum - wounded3Momentum;
156
157 if(anA>0 && aZ>0) {
159#ifdef exactExcitationEnergy
160 // recalculate exEnergy from Energy balance....
161 const G4HadProjectile * primary = GetPrimaryProjectile();
162 G4double Einitial= primary->Get4Momentum().e()
163 + G4NucleiProperties::GetNuclearMass(theNucleus->GetMassNumber(),theNucleus->GetCharge());
164 G4double Efinal = fMass + secondary4Momemtum.e();
165 if ( (Einitial - Efinal) > 0 ) {
166 // G4cout << "G4GPI::Propagate() : positive exact excitation Energy "
167 // << (Einitial - Efinal)/MeV << " MeV, exciton estimate " << exEnergy/MeV << " MeV" << G4endl;
168 exEnergy=Einitial - Efinal;
169 }
170 else {
171 // G4cout << "G4GeneratorPrecompoundInterface::Propagate() : negative exact excitation Energy "
172 // << (Einitial - Efinal)/MeV << " MeV, setting excitation to 0 MeV" << G4endl;
173 exEnergy=0.;
174 }
175#endif
176 fMass += exEnergy;
177 G4ThreeVector balance=primary->Get4Momentum().vect() - secondary4Momemtum.vect() - exciton3Momentum;
178 #ifdef G4GPI_debug_excitation
179 G4cout << "momentum balance init/final " << balance << " value " << balance.mag() << G4endl
180 << "primary / secondaries "<< primary->Get4Momentum() << " / "
181 << secondary4Momemtum << " captured/wounded: " << captured3Momentum << " / " << wounded3Momentum
182 << " exciton " << exciton3Momentum << G4endl
183 << secondary4Momemtum.vect() + exciton3Momentum << G4endl;
184 #endif
185#ifdef exactExcitationEnergy
186 G4LorentzVector exciton4Momentum(exciton3Momentum, fMass);
187#else
188 G4LorentzVector exciton4Momentum(exciton3Momentum,
189 std::sqrt(exciton3Momentum.mag2() + fMass*fMass));
190#endif
191 if ( exEnergy > 0.0 ) { // Need to de-excite the remnant nucleus only if excitation energy > 0.
192 G4Fragment anInitialState(anA, aZ, exciton4Momentum);
193 anInitialState.SetNumberOfParticles(numberOfEx-numberOfHoles);
194 anInitialState.SetNumberOfCharged(numberOfCh);
195 anInitialState.SetNumberOfHoles(numberOfHoles);
196
197 G4ReactionProductVector * aPrecoResult = theDeExcitation->DeExcite(anInitialState);
198 // fill pre-compound part into the result, and return
199 theTotalResult->insert(theTotalResult->end(),aPrecoResult->begin(),aPrecoResult->end() );
200 delete aPrecoResult;
201
202 } else { // No/negative excitation energy, we only need to create the remnant nucleus
203 // energy is not conserved, ignore exciton momentum, i.e. remnant nucleus will be at rest
204 G4ParticleDefinition* theKindOfFragment = 0;
205 if (anA == 1 && aZ == 0) {
206 theKindOfFragment = G4Neutron::NeutronDefinition();
207 } else if (anA == 1 && aZ == 1) {
208 theKindOfFragment = G4Proton::ProtonDefinition();
209 } else if (anA == 2 && aZ == 1) {
210 theKindOfFragment = G4Deuteron::DeuteronDefinition();
211 } else if (anA == 3 && aZ == 1) {
212 theKindOfFragment = G4Triton::TritonDefinition();
213 } else if (anA == 3 && aZ == 2) {
214 theKindOfFragment = G4He3::He3Definition();
215 } else if (anA == 4 && aZ == 2) {
216 theKindOfFragment = G4Alpha::AlphaDefinition();;
217 } else {
218 theKindOfFragment =
220 }
221 if (theKindOfFragment != 0) {
222 G4ReactionProduct * theNew = new G4ReactionProduct(theKindOfFragment);
223 theNew->SetMomentum(G4ThreeVector(0.,0.,0.));
224 theNew->SetTotalEnergy(fMass);
225 //theNew->SetFormationTime(??0.??);
226 theTotalResult->push_back(theNew);
227 }
228 }
229 }
230
231 return theTotalResult;
232}
233
236{
237 G4cout << "G4GeneratorPrecompoundInterface: ApplyYourself interface called stand-allone."
238 << G4endl;
239 G4cout << "This class is only a mediator between generator and precompound"<<G4endl;
240 G4cout << "Please remove from your physics list."<<G4endl;
241 throw G4HadronicException(__FILE__, __LINE__, "SEVERE: G4GeneratorPrecompoundInterface model interface called stand-allone.");
242 return new G4HadFinalState;
243}
245{
246 outFile << "G4GeneratorPrecompoundInterface interfaces a high\n"
247 << "energy model through the wounded nucleus to precompound de-excition.\n"
248 << "Low energy protons and neutron present among secondaries produced by \n"
249 << "the high energy generator and within the nucleus are captured. The wounded\n"
250 << "nucleus and the captured particles form an excited nuclear fragment. This\n"
251 << "fragment is passed to the Geant4 pre-compound model for de-excitation.\n"
252 << "Nuclear de-excitation:\n";
253 // preco
254
255}
std::vector< G4ReactionProduct * > G4ReactionProductVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
double mag2() const
double mag() const
Hep3Vector vect() const
static G4Alpha * AlphaDefinition()
Definition: G4Alpha.cc:84
static G4Deuteron * DeuteronDefinition()
Definition: G4Deuteron.cc:89
void SetNumberOfCharged(G4int value)
Definition: G4Fragment.hh:349
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:335
void SetNumberOfParticles(G4int value)
Definition: G4Fragment.hh:344
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
virtual G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)
virtual void PropagateModelDescription(std::ostream &) const
G4GeneratorPrecompoundInterface(G4VPreCompoundModel *p=0)
const G4LorentzVector & Get4Momentum() const
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
static G4He3 * He3Definition()
Definition: G4He3.cc:89
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int J=0)
Definition: G4IonTable.cc:267
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4bool AreYouHit() const
Definition: G4Nucleon.hh:97
virtual const G4LorentzVector & Get4Momentum() const
Definition: G4Nucleon.hh:72
G4double GetBindingEnergy() const
Definition: G4Nucleon.hh:75
virtual G4ParticleDefinition * GetDefinition() const
Definition: G4Nucleon.hh:85
G4double GetPDGCharge() const
static G4ParticleTable * GetParticleTable()
G4IonTable * GetIonTable()
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
static G4Triton * TritonDefinition()
Definition: G4Triton.cc:90
virtual G4Nucleon * GetNextNucleon()=0
virtual G4int GetCharge()=0
virtual G4bool StartLoop()=0
virtual G4double GetNuclearRadius()=0
virtual G4int GetMassNumber()=0
const G4HadProjectile * GetPrimaryProjectile() const
void SetDeExcitation(G4VPreCompoundModel *ptr)
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment)=0