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
G4CascadeRecoilMaker.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// Collects generated cascade data (using Collider::collide() interface)
28// and computes the nuclear recoil kinematics needed to balance the event.
29//
30// 20100909 M. Kelsey -- Inspired by G4CascadeCheckBalance
31// 20100909 M. Kelsey -- Move G4IntraNucleiCascader::goodCase() here, add
32// tolerance for "almost zero" excitation energy
33// 20100910 M. Kelsey -- Drop getRecoilFragment() in favor of user calling
34// makeRecoilFragment() with returned non-const pointer. Drop
35// handling of excitons.
36// 20100921 M. Kelsey -- Return G4InuclNuclei using "makeRecoilNuclei()".
37// Repurpose "makeRecoilFragment()" to return G4Fragment.
38// 20100924 M. Kelsey -- Remove unusable G4Fragment::SetExcitationEnergy().
39// Add deltaM to compute mass difference, use excitationEnergy
40// to force G4Fragment four-vector to match.
41// 20110214 M. Kelsey -- Follow G4InuclParticle::Model enumerator migration
42// 20110303 M. Kelsey -- Add diagnostic messages to goodNucleus().
43// 20110308 M. Kelsey -- Follow new G4Fragment interface for hole types
44// 20110722 M. Kelsey -- For IntraNucleiCascader, take G4CollOut as argument
45
46#include <vector>
47
49#include "globals.hh"
50#include "G4SystemOfUnits.hh"
51#include "G4CascadParticle.hh"
53#include "G4CollisionOutput.hh"
54#include "G4Fragment.hh"
56#include "G4InuclNuclei.hh"
57#include "G4InuclParticle.hh"
59#include "G4LorentzVector.hh"
60
61using namespace G4InuclSpecialFunctions;
62
63
64// Constructor and destructor
65
67 : G4VCascadeCollider("G4CascadeRecoilMaker"),
68 excTolerance(tolerance), inputEkin(0.),
69 recoilA(0), recoilZ(0), excitationEnergy(0.) {
70 balance = new G4CascadeCheckBalance(tolerance, tolerance, theName);
71}
72
74 delete balance;
75}
76
77
78// Standard Collider interface (non-const output "buffer")
79
81 G4InuclParticle* target,
82 G4CollisionOutput& output) {
83 if (verboseLevel > 1)
84 G4cout << " >>> G4CascadeRecoilMaker::collide" << G4endl;
85
86 // Available energy needed for "goodNucleus()" test at end
87 inputEkin = bullet ? bullet->getKineticEnergy() : 0.;
88
90 balance->collide(bullet, target, output);
91 fillRecoil();
92}
93
94// This is for use with G4IntraNucleiCascader
95
97 G4InuclParticle* target,
98 G4CollisionOutput& output,
99 const std::vector<G4CascadParticle>& cparticles) {
100 if (verboseLevel > 1)
101 G4cout << " >>> G4CascadeRecoilMaker::collide(<EP>,<CP>)" << G4endl;
102
103 // Available energy needed for "goodNucleus()" test at end
104 inputEkin = bullet ? bullet->getKineticEnergy() : 0.;
105
107 balance->collide(bullet, target, output, cparticles);
108 fillRecoil();
109}
110
111
112// Used current event configuration to construct recoil nucleus
113// NOTE: CheckBalance uses "final-initial", we want "initial-final"
114
116 recoilZ = -(balance->deltaQ()); // Charge "non-conservation"
117 recoilA = -(balance->deltaB()); // Baryon "non-conservation"
118 recoilMomentum = -(balance->deltaLV());
119
120 theExcitons.clear(); // Discard previous exciton configuraiton
121
122 // Bertini uses MeV for excitation energy
123 if (!goodFragment()) excitationEnergy = 0.;
124 else excitationEnergy = deltaM() * GeV;
125
126 // Allow for very small negative mass difference, and round to zero
127 if (std::abs(excitationEnergy) < excTolerance) excitationEnergy = 0.;
128
129 if (verboseLevel > 2) {
130 G4cout << " recoil px " << recoilMomentum.px()
131 << " py " << recoilMomentum.py() << " pz " << recoilMomentum.pz()
132 << " E " << recoilMomentum.e() << " baryon " << recoilA
133 << " charge " << recoilZ
134 << "\n recoil mass " << recoilMomentum.m()
135 << " 'excitation' energy " << excitationEnergy << G4endl;
136 }
137}
138
139
140// Construct physical nucleus from recoil parameters, if reasonable
141
144 if (verboseLevel > 1)
145 G4cout << " >>> G4CascadeRecoilMaker::makeRecoilNuclei" << G4endl;
146
147 if (!goodRecoil()) {
148 if (verboseLevel > 2 && !wholeEvent())
149 G4cout << theName << ": event recoil is not a physical nucleus" << G4endl;
150
151 return 0; // Null pointer means no fragment
152 }
153
154 theRecoilNuclei.fill(recoilMomentum, recoilA, recoilZ,
155 excitationEnergy, model);
156 theRecoilNuclei.setExitonConfiguration(theExcitons);
157
158 return &theRecoilNuclei;
159}
160
161
162// Construct pre-compound nuclear fragment from recoil parameters
163
165 if (verboseLevel > 1)
166 G4cout << " >>> G4CascadeRecoilMaker::makeRecoilFragment" << G4endl;
167
168 if (!goodRecoil()) {
169 if (verboseLevel > 2 && !wholeEvent())
170 G4cout << theName << ": event recoil is not a physical nucleus" << G4endl;
171
172 return 0; // Null pointer means no fragment
173 }
174
175 theRecoilFragment.SetZandA_asInt(recoilZ, recoilA); // Note convention!
176
177 // User may have overridden excitation energy; force four-momentum to match
178 G4double fragMass =
179 G4InuclNuclei::getNucleiMass(recoilA,recoilZ) + excitationEnergy/GeV;
180
181 G4LorentzVector fragMom; fragMom.setVectM(recoilMomentum.vect(), fragMass);
182 theRecoilFragment.SetMomentum(fragMom*GeV); // Bertini uses GeV!
183
184 // Note: exciton configuration has to be set piece by piece
185 // (arguments are Ntotal,Nproton in both cases)
186 G4int nholes = theExcitons.protonHoles+theExcitons.neutronHoles;
187 theRecoilFragment.SetNumberOfHoles(nholes, theExcitons.protonHoles);
188
189 G4int nexcit = (theExcitons.protonQuasiParticles
190 + theExcitons.neutronQuasiParticles);
191 theRecoilFragment.SetNumberOfExcitedParticle(nexcit,
192 theExcitons.protonQuasiParticles);
193
194 return &theRecoilFragment;
195}
196
197
198// Compute raw mass difference from recoil parameters
199
201 G4double nucMass = G4InuclNuclei::getNucleiMass(recoilA,recoilZ);
202 return (recoilMomentum.m() - nucMass);
203}
204
205
206// Data quality checks
207
209 return (recoilA>0 && recoilZ>=0 && recoilA >= recoilZ);
210}
211
213 return (goodFragment() && excitationEnergy > -excTolerance);
214}
215
217 if (verboseLevel > 2) {
218 G4cout << " >>> G4CascadeRecoilMaker::wholeEvent:"
219 << " A " << recoilA << " Z " << recoilZ
220 << " P " << recoilMomentum.rho() << " E " << recoilMomentum.e()
221 << "\n wholeEvent returns "
222 << (recoilA==0 && recoilZ==0 &&
223 recoilMomentum.rho() < excTolerance/GeV &&
224 std::abs(recoilMomentum.e()) < excTolerance/GeV) << G4endl;
225 }
226
227 return (recoilA==0 && recoilZ==0 &&
228 recoilMomentum.rho() < excTolerance/GeV &&
229 std::abs(recoilMomentum.e()) < excTolerance/GeV);
230}
231
232// Determine whether desired nuclear fragment is constructable outcome
233
235 if (verboseLevel > 2) {
236 G4cout << " >>> G4CascadeRecoilMaker::goodNucleus" << G4endl;
237 }
238
239 const G4double minExcitation = 0.1*keV;
240 const G4double reasonableExcitation = 7.0; // Multiple of binding energy
241 const G4double fractionalExcitation = 0.2; // Fraction of input to excite
242
243 if (!goodRecoil()) {
244 if (verboseLevel>2) {
245 if (!goodFragment()) G4cerr << " goodNucleus: invalid A/Z" << G4endl;
246 else if (excitationEnergy < -excTolerance)
247 G4cerr << " goodNucleus: negative excitation" << G4endl;
248 }
249 return false; // Not a sensible nucleus
250 }
251
252 if (excitationEnergy <= minExcitation) return true; // Effectively zero
253
254 // Maximum possible excitation energy determined by initial energy
255 G4double dm = bindingEnergy(recoilA,recoilZ);
256 G4double exc_max0z = fractionalExcitation * inputEkin*GeV;
257 G4double exc_dm = reasonableExcitation * dm;
258 G4double exc_max = (exc_max0z > exc_dm) ? exc_max0z : exc_dm;
259
260 if (verboseLevel > 3) {
261 G4cout << " eexs " << excitationEnergy << " max " << exc_max
262 << " dm " << dm << G4endl;
263 }
264
265 if (verboseLevel > 2 && excitationEnergy >= exc_max)
266 G4cerr << " goodNucleus: too much excitation" << G4endl;
267
268 return (excitationEnergy < exc_max); // Below maximum possible
269}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#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)
G4LorentzVector deltaLV() const
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
G4InuclNuclei * makeRecoilNuclei(G4InuclParticle::Model model=G4InuclParticle::DefaultModel)
G4Fragment * makeRecoilFragment()
G4CascadeRecoilMaker(G4double tolerance=0.001 *CLHEP::MeV)
void SetZandA_asInt(G4int Znew, G4int Anew, G4int Lnew=0)
Definition: G4Fragment.hh:295
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
void setExitonConfiguration(const G4ExitonConfiguration &config)
G4double getNucleiMass() const
void fill(G4int a, G4int z, G4double exc=0., Model model=DefaultModel)
G4double getKineticEnergy() const
virtual void setVerboseLevel(G4int verbose=0)
G4double bindingEnergy(G4int A, G4int Z)