Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4CascadeDeexciteBase.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//
27// Semi-concrete base class for de-excitation modules, analogous to
28// G4CascadeColliderBase.
29//
30// 20130806 M. Kelsey -- Per A. Dotti, move zero vector to file scope to
31// address thread-collision problem.
32
36#include "G4CollisionOutput.hh"
37#include "G4Fragment.hh"
38#include "G4InuclNuclei.hh"
40#include "G4SystemOfUnits.hh"
41#include <vector>
42
43using namespace G4InuclSpecialFunctions;
44
45
46// Constructor and destructor
47
49 : G4VCascadeDeexcitation(name), balance(0), A(0), Z(0), EEXS(0) {
52}
53
55 delete balance;
56}
57
60 if (balance) balance->setVerboseLevel(verbose);
61}
62
63
64// Copy pertinent information from G4Fragment for modules
65
67 A = target.GetA_asInt();
68 Z = target.GetZ_asInt();
69 PEX = target.GetMomentum()/GeV; // Convert from G4 to Bertini units
70 EEXS = target.GetExcitationEnergy();
71}
72
73
74// Create (fill) new G4Fragment with proper momentum/energy handling
75
76namespace {
77 static const G4LorentzVector zero(0.,0.,0.,0.); // File scope avoids churn
78}
79
80const G4Fragment&
82 return makeFragment(zero, fragA, fragZ, EX);
83}
84
85const G4Fragment&
87 G4int fragZ, G4double EX) {
88 if (verboseLevel>2) {
89 G4cout << " >>> " << theName << "::makeFragment " << mom << " " << fragA
90 << " " << fragZ << " " << EX << G4endl;
91 }
92
93 // Adjust four-momentum so that mass is nucleus + excitation
94 G4double mass =
95 G4InuclNuclei::getNucleiMass(fragA,fragZ) + EX/GeV;
96 mom.setVectM(mom.vect(), mass);
97
98 // Overwrite previous fragment contents, zeroing out excitons
99 aFragment.SetZandA_asInt(fragZ, fragA);
100 aFragment.SetMomentum(mom*GeV); // Bertini uses GeV!
103
104 return aFragment;
105}
106
107// Decide wether nuclear fragment is candidate for G4BigBanger
108
110 return explosion(fragment.GetA_asInt(), fragment.GetZ_asInt(),
111 fragment.GetExcitationEnergy()); // in MeV
112}
113
115 G4double excitation) const {
116 if (verboseLevel) G4cout << " >>> " << theName << "::explosion ?" << G4endl;
117
118 const G4int a_cut = 20;
119 const G4double be_cut = 3.0;
120
121 // Neutron balls, or small fragments with high excitations can explode
122 return ((fragA <= a_cut || fragZ==0) &&
123 (excitation >= be_cut * bindingEnergy(fragA,fragZ))
124 );
125}
126
127
128// Validate output for energy, momentum conservation, etc.
129
131 G4CollisionOutput& output) {
132 if (!balance) return true; // Skip checks unless requested
133
134 if (verboseLevel > 1)
135 G4cout << " >>> " << theName << "::validateOutput" << G4endl;
136
138 balance->collide(target, output);
139 return balance->okay(); // Returns false if violations
140}
141
143 const std::vector<G4InuclElementaryParticle>& particles) {
144 if (!balance) return true; // Skip checks unless requested
145
146 if (verboseLevel > 1)
147 G4cout << " >>> " << theName << "::validateOutput" << G4endl;
148
150 balance->collide(target, particles);
151 return balance->okay(); // Returns false if violations
152}
153
155 const std::vector<G4InuclNuclei>& fragments) {
156 if (!balance) return true; // Skip checks unless requested
157
158 if (verboseLevel > 1)
159 G4cout << " >>> " << theName << "::validateOutput" << G4endl;
160
162 balance->collide(target, fragments);
163 return balance->okay(); // Returns false if violations
164}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void setVectM(const Hep3Vector &spatial, double mass)
Hep3Vector vect() const
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
G4CascadeDeexciteBase(const char *name)
void getTargetData(const G4Fragment &target)
const G4Fragment & makeFragment(G4LorentzVector mom, G4int A, G4int Z, G4double EX=0.)
virtual void setVerboseLevel(G4int verbose=0)
virtual G4bool explosion(const G4Fragment &target) const
G4CascadeCheckBalance * balance
virtual G4bool validateOutput(const G4Fragment &target, G4CollisionOutput &output)
static G4bool checkConservation()
void SetZandA_asInt(G4int Znew, G4int Anew, G4int Lnew=0)
Definition: G4Fragment.hh:295
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:312
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:322
G4int GetZ_asInt() const
Definition: G4Fragment.hh:289
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:396
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:327
void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
Definition: G4Fragment.hh:377
G4int GetA_asInt() const
Definition: G4Fragment.hh:284
G4double getNucleiMass() const
virtual void setVerboseLevel(G4int verbose=0)
G4double bindingEnergy(G4int A, G4int Z)