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
G4CascadeCheckBalance.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// Verify and report four-momentum conservation for collision output; uses
28// same interface as collision generators.
29//
30// 20100624 M. Kelsey -- Add baryon number, charge, and kinetic energy
31// 20100624 M. Kelsey -- Bug fix: All checks should be |delta| !
32// 20100627 M. Kelsey -- Always report violations on cerr, regardless of
33// verbosity.
34// 20100628 M. Kelsey -- Add interface to take list of particles directly,
35// bug fix reporting of charge conservation error.
36// 20100630 M. Kelsey -- for nuclei, include excitation energies in total.
37// 20100701 M. Kelsey -- Undo previous change, handled by G4InuclNuclei.
38// 20100711 M. Kelsey -- Use name of parent collider for reporting messages
39// 20100713 M. kelsey -- Hide conservation errors behind verbosity
40// 20100715 M. Kelsey -- Use new G4CollisionOutput totals instead of loops,
41// move temporary buffer to be data member
42// 20100719 M. Kelsey -- Change zero tolerance to 10 keV instead of 1 keV.
43// 20100909 M. Kelsey -- Add interface for both kinds of particle lists
44// 20101019 M. Kelsey -- CoVerity report: unitialized constructor
45// 20110328 M. Kelsey -- Add default ctor and explicit limit setting
46// 20110722 M. Kelsey -- For IntraNucleiCascader, take G4CollOut as argument
47// 20120525 M. Kelsey -- Follow process-level checking: allow _either_ rel.
48// or abs. to pass (instead of requiring both)
49// 20121002 M. Kelsey -- Add strangeness check (useful for Omega- beam)
50// 20130621 Add interface to take G4Fragment input instead of G4InuclNuclei.
51// 20140930 Change name from "const char*" to "const G4String"
52// 20150626 Fix logical condition for report relative or absolute violations.
53
55#include "globals.hh"
56#include "G4CascadParticle.hh"
58#include "G4InuclNuclei.hh"
59#include "G4InuclParticle.hh"
60#include "G4CollisionOutput.hh"
61#include "G4LorentzVector.hh"
62#include "G4Electron.hh"
63#include "G4SystemOfUnits.hh"
64#include <vector>
65
66
67// Constructor sets acceptance limits
68
69 const G4double G4CascadeCheckBalance::tolerance = 1e-6; // How small is zero?
70
72 : G4VCascadeCollider(owner), relativeLimit(G4CascadeCheckBalance::tolerance),
73 absoluteLimit(G4CascadeCheckBalance::tolerance), initialBaryon(0),
74 finalBaryon(0), initialCharge(0), finalCharge(0), initialStrange(0),
75 finalStrange(0) {}
76
78 G4double absolute,
79 const G4String& owner)
80 : G4VCascadeCollider(owner), relativeLimit(relative),
81 absoluteLimit(absolute), initialBaryon(0), finalBaryon(0),
82 initialCharge(0), finalCharge(0), initialStrange(0),
83 finalStrange(0) {}
84
85
86// Pseudo-collision just computes input and output four-vectors
87
89 G4InuclParticle* target,
90 G4CollisionOutput& output) {
91 if (verboseLevel)
92 G4cout << " >>> G4CascadeCheckBalance(" << theName << ")::collide"
93 << G4endl;
94
95 initial *= 0.; // Fast reset; some colliders only have one pointer
96 if (bullet) initial += bullet->getMomentum();
97 if (target) initial += target->getMomentum();
98
99 // Baryon number, charge and strangeness must be computed "by hand"
100 initialCharge = 0;
101 if (bullet) initialCharge += G4int(bullet->getCharge());
102 if (target) initialCharge += G4int(target->getCharge());
103
105 dynamic_cast<G4InuclElementaryParticle*>(bullet);
107 dynamic_cast<G4InuclElementaryParticle*>(target);
108
109 G4InuclNuclei* nbullet = dynamic_cast<G4InuclNuclei*>(bullet);
110 G4InuclNuclei* ntarget = dynamic_cast<G4InuclNuclei*>(target);
111
112 initialBaryon =
113 ((pbullet ? pbullet->baryon() : nbullet ? nbullet->getA() : 0) +
114 (ptarget ? ptarget->baryon() : ntarget ? ntarget->getA() : 0) );
115
116 // NOTE: Currently we ignore possibility of hypernucleus target
117 initialStrange = 0;
118 if (pbullet) initialStrange += pbullet->getStrangeness();
119 if (ptarget) initialStrange += ptarget->getStrangeness();
120
121 G4int nelec = 0;
122 G4double elMass = 0.;
123 std::vector<G4InuclElementaryParticle>& outParts =
124 output.getOutgoingParticles();
125 for (G4int i = 0; i < output.numberOfOutgoingParticles(); i++) {
126 if (outParts[i].getDefinition() == G4Electron::Electron() ) {
127 elMass += outParts[i].getDefinition()->GetPDGMass();
128 ++nelec;
129 }
130 }
131 if(nelec > 0) {
132 initial += G4LorentzVector(0.,0.,0.,elMass/GeV);
133 initialCharge -= nelec;
134 }
135
136 // Final state totals are computed for us
137 final = output.getTotalOutputMomentum();
138 finalBaryon = output.getTotalBaryonNumber();
139 finalCharge = output.getTotalCharge();
140 finalStrange = output.getTotalStrangeness();
141
142 // Report results
143 if (verboseLevel) {
144 G4cout << " initial px " << initial.px() << " py " << initial.py()
145 << " pz " << initial.pz() << " E " << initial.e()
146 << " baryon " << initialBaryon << " charge " << initialCharge
147 << " strange " << initialStrange << G4endl
148 << " final px " << final.px() << " py " << final.py()
149 << " pz " << final.pz() << " E " << final.e()
150 << " baryon " << finalBaryon << " charge " << finalCharge
151 << " strange " << finalStrange << G4endl;
152 }
153}
154
155// For de-excitation, take G4Fragment as initial state
157 G4CollisionOutput& output) {
158 if (verboseLevel)
159 G4cout << " >>> G4CascadeCheckBalance(" << theName << ")::collide(<FRAG>)"
160 << G4endl;
161
162 // Copy initial state directly from fragment (no bullet/target sums)
163 initial = fragment.GetMomentum()/GeV; // G4Fragment returns MeV
164 initialCharge = fragment.GetZ_asInt();
165 initialBaryon = fragment.GetA_asInt();
166 initialStrange = 0; // No hypernuclei at present
167
168 // Final state totals are computed for us
169 final = output.getTotalOutputMomentum();
170
171 // Remove electron masses when internal conversion occurs
172 G4double elMass = 0.;
173 std::vector<G4InuclElementaryParticle>& outParts =
174 output.getOutgoingParticles();
175 G4int nelec = 0;
176 for (G4int i = 0; i < output.numberOfOutgoingParticles(); i++) {
177 if (outParts[i].getDefinition() == G4Electron::Electron() ) {
178 elMass += outParts[i].getDefinition()->GetPDGMass();
179 ++nelec;
180 }
181 }
182 if(nelec > 0) {
183 initial += G4LorentzVector(0.,0.,0.,elMass/GeV);
184 initialCharge -= nelec;
185 }
186
187 finalBaryon = output.getTotalBaryonNumber();
188 finalCharge = output.getTotalCharge();
189 finalStrange = output.getTotalStrangeness();
190
191 // Report results
192 if (verboseLevel) {
193 G4cout << " initial px " << initial.px() << " py " << initial.py()
194 << " pz " << initial.pz() << " E " << initial.e()
195 << " baryon " << initialBaryon << " charge " << initialCharge
196 << " strange " << initialStrange << G4endl
197 << " final px " << final.px() << " py " << final.py()
198 << " pz " << final.pz() << " E " << final.e()
199 << " baryon " << finalBaryon << " charge " << finalCharge
200 << " strange " << finalStrange << G4endl;
201 }
202}
203
204
205// Take list of output particles directly (e.g., from G4EPCollider internals)
206
208 G4InuclParticle* target,
209 const std::vector<G4InuclElementaryParticle>& particles) {
210 if (verboseLevel)
211 G4cout << " >>> G4CascadeCheckBalance(" << theName << ")::collide(<vector>)"
212 << G4endl;
213
214 tempOutput.reset(); // Buffer for processing
215 tempOutput.addOutgoingParticles(particles);
216 collide(bullet, target, tempOutput);
217}
218
220 const std::vector<G4InuclElementaryParticle>& particles) {
221 if (verboseLevel)
222 G4cout << " >>> G4CascadeCheckBalance(" << theName
223 << ")::collide(<FRAG>,<vector>)" << G4endl;
224
225 tempOutput.reset(); // Buffer for processing
226 tempOutput.addOutgoingParticles(particles);
227 collide(target, tempOutput);
228}
229
230
231// Take list of nuclear fragments directly (e.g., from G4Fissioner internals)
233 const std::vector<G4InuclNuclei>& fragments) {
234 if (verboseLevel)
235 G4cout << " >>> G4CascadeCheckBalance(" << theName << ")::collide(<vector>)"
236 << G4endl;
237
238 tempOutput.reset(); // Buffer for processing
239 tempOutput.addOutgoingNuclei(fragments);
240 collide(target, tempOutput);
241}
242
243
244// Take list of "cparticles" (e.g., from G4NucleiModel internals)
246 G4InuclParticle* target,
247 const std::vector<G4CascadParticle>& particles) {
248 if (verboseLevel)
249 G4cout << " >>> G4CascadeCheckBalance(" << theName
250 << ")::collide(<cparticles>)" << G4endl;
251
252 tempOutput.reset(); // Buffer for processing
253 tempOutput.addOutgoingParticles(particles);
254 collide(bullet, target, tempOutput);
255}
256
257
258// Take lists of both G4InuclEP & G4CP (e.g., from G4IntraNucleiCascader)
260 G4InuclParticle* target,
261 G4CollisionOutput& output,
262 const std::vector<G4CascadParticle>& cparticles) {
263 if (verboseLevel)
264 G4cout << " >>> G4CascadeCheckBalance(" << theName
265 << ")::collide(<EP>,<CP>)" << G4endl;
266
267 tempOutput.reset(); // Buffer for processing
268 tempOutput.add(output);
269 tempOutput.addOutgoingParticles(cparticles);
270 collide(bullet, target, tempOutput);
271}
272
273// Compare relative and absolute violations to limits, and report
274
276 G4bool relokay = (std::abs(relativeE()) < relativeLimit);
277 G4bool absokay = (std::abs(deltaE()) < absoluteLimit);
278
279 if (verboseLevel && (!relokay || !absokay)) {
280 G4cerr << theName << ": Energy conservation: relative " << relativeE()
281 << (relokay ? " conserved" : " VIOLATED")
282 << " absolute " << deltaE()
283 << (absokay ? " conserved" : " VIOLATED") << G4endl;
284 } else if (verboseLevel > 1) {
285 G4cout << theName << ": Energy conservation: relative " << relativeE()
286 << " conserved absolute " << deltaE() << " conserved" << G4endl;
287 }
288
289 return (relokay && absokay);
290}
291
292
294 G4bool relokay = (std::abs(relativeKE()) < relativeLimit);
295 G4bool absokay = (std::abs(deltaKE()) < absoluteLimit);
296
297 if (verboseLevel && (!relokay || !absokay)) {
298 G4cerr << theName << ": Kinetic energy balance: relative "
299 << relativeKE() << (relokay ? " conserved" : " VIOLATED")
300 << " absolute " << deltaKE()
301 << (absokay ? " conserved" : " VIOLATED") << G4endl;
302 } else if (verboseLevel > 1) {
303 G4cout << theName << ": Kinetic energy balance: relative "
304 << relativeKE() << " conserved absolute " << deltaKE()
305 << " conserved" << G4endl;
306 }
307
308 return (relokay && absokay);
309}
310
311
313 G4bool relokay = (std::abs(relativeP()) < 10.*relativeLimit);
314 G4bool absokay = (std::abs(deltaP()) < 10.*absoluteLimit);
315
316 if (verboseLevel && (!relokay || !absokay)) {
317 G4cerr << theName << ": Momentum conservation: relative " << relativeP()
318 << (relokay ? " conserved" : " VIOLATED")
319 << " absolute " << deltaP()
320 << (absokay ? " conserved" : " VIOLATED") << G4endl;
321 } else if (verboseLevel > 1) {
322 G4cout << theName << ": Momentum conservation: relative " << relativeP()
323 << " conserved absolute " << deltaP() << " conserved" << G4endl;
324 }
325
326 return (relokay && absokay);
327}
328
329
331 G4bool bokay = (deltaB() == 0); // Must be perfect!
332
333 if (verboseLevel && !bokay)
334 G4cerr << theName << ": Baryon number VIOLATED " << deltaB() << G4endl;
335 return bokay;
336}
337
339 G4bool qokay = (deltaQ() == 0); // Must be perfect!
340
341 if (verboseLevel && !qokay)
342 G4cerr << theName << ": Charge conservation VIOLATED " << deltaQ()
343 << G4endl;
344 return qokay;
345}
346
348 G4bool sokay = (deltaS() == 0); // Must be perfect!
349
350 if (verboseLevel && !sokay)
351 G4cerr << theName << ": Strangeness conservation VIOLATED " << deltaS()
352 << G4endl;
353
354 return sokay;
355}
CLHEP::HepLorentzVector G4LorentzVector
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 collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
static const G4double tolerance
G4CascadeCheckBalance(const G4String &owner="G4CascadeCheckBalance")
G4int numberOfOutgoingParticles() const
G4int getTotalStrangeness() const
G4LorentzVector getTotalOutputMomentum() const
G4int getTotalBaryonNumber() const
G4int getTotalCharge() const
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
void add(const G4CollisionOutput &right)
void addOutgoingNuclei(const std::vector< G4InuclNuclei > &nuclea)
void addOutgoingParticles(const std::vector< G4InuclElementaryParticle > &particles)
static G4Electron * Electron()
Definition: G4Electron.cc:93
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:322
G4int GetZ_asInt() const
Definition: G4Fragment.hh:289
G4int GetA_asInt() const
Definition: G4Fragment.hh:284
G4int getA() const
G4LorentzVector getMomentum() const
G4double getCharge() const