Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CascadeCheckBalance.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 G4CASCADE_CHECK_BALANCE_HH
27#define G4CASCADE_CHECK_BALANCE_HH
28//
29// Verify and report four-momentum conservation for collision output; uses
30// same interface as collision generators.
31//
32// 20100624 M. Kelsey -- Add baryon conservation check and kinetic energy
33// 20100628 M. Kelsey -- Add interface to take list of particles directly
34// 20100711 M. Kelsey -- Add name of parent Collider for reporting messages,
35// allow changing parent name, add interface for nuclear fragments
36// 20100713 M. Kelsey -- Add (publicly adjustable) tolerance for zeroes
37// 20100715 M. Kelsey -- FPE! Need to check initial values before computing
38// relative error.
39// 20100715 M. Kelsey -- Add G4CascadParticle interface for G4NucleiModel;
40// do momentum check on direction, not just magnitude. Move
41// temporary G4CollisionOutput buffer here, for thread-safety
42// 20100909 M. Kelsey -- Add interface to get four-vector difference, and
43// to supply both kinds of particle lists (G4IntraNucleiCascader)
44// 20100923 M. Kelsey -- Baryon and charge deltas should have been integer
45// 20110328 M. Kelsey -- Add default ctor and explicit limit setting
46// 20110722 M. Kelsey -- For IntraNucleiCascader, take G4CollOut as argument
47// 20121002 M. Kelsey -- Add strangeness check (useful for Omega- beam)
48// 20130620 Address Coverity complaint about missing copy actions
49// 20130621 Add interface to take G4Fragment input instead of G4InuclNuclei.
50// 20140930 Change name from "const char*" to "const G4String"
51
52#include "G4VCascadeCollider.hh"
53#include "globals.hh"
54#include "G4CollisionOutput.hh"
55#include "G4LorentzVector.hh"
56#include <cmath>
57#include <vector>
58
61class G4InuclNuclei;
62class G4InuclParticle;
63
65public:
66 static const G4double tolerance; // Don't do floating zero!
67
68 explicit G4CascadeCheckBalance(const G4String& owner="G4CascadeCheckBalance");
69
70 G4CascadeCheckBalance(G4double relative, G4double absolute,
71 const G4String& owner="G4CascadeCheckBalance");
73
74 void setOwner(const G4String& owner) { setName(owner); }
75
76 void setLimits(G4double relative, G4double absolute) {
77 setRelativeLimit(relative);
78 setAbsoluteLimit(absolute);
79 }
80
81 void setRelativeLimit(G4double limit) { relativeLimit = limit; }
82 void setAbsoluteLimit(G4double limit) { absoluteLimit = limit; }
83
84 void collide(G4InuclParticle* bullet, G4InuclParticle* target,
85 G4CollisionOutput& output);
86
87 // This is for use with G4VCascadeDeexcitation modules
88 void collide(const G4Fragment& fragment, G4CollisionOutput& output);
89
90 // This is for use with G4EPCollider internal checks
91 void collide(G4InuclParticle* bullet, G4InuclParticle* target,
92 const std::vector<G4InuclElementaryParticle>& particles);
93
94 // This is for use with G4NucleiModel internal checks
95 void collide(G4InuclParticle* bullet, G4InuclParticle* target,
96 const std::vector<G4CascadParticle>& particles);
97
98 // This is for use with G4IntraNucleiCascader
99 void collide(G4InuclParticle* bullet, G4InuclParticle* target,
100 G4CollisionOutput& output,
101 const std::vector<G4CascadParticle>& cparticles);
102
103 // This is for use with G4BigBanger internal checks
104 void collide(const G4Fragment& target,
105 const std::vector<G4InuclElementaryParticle>& particles);
106
107 // This is for use with G4Fissioner internal checks
108 void collide(const G4Fragment& target,
109 const std::vector<G4InuclNuclei>& fragments);
110
111 // Checks on conservation laws (kinematics, baryon number, charge, hyperons)
112 G4bool energyOkay() const;
113 G4bool ekinOkay() const;
114 G4bool momentumOkay() const;
115 G4bool baryonOkay() const;
116 G4bool chargeOkay() const;
117 G4bool strangeOkay() const;
118
119 // Global check, used by G4CascadeInterface validation loop
120 // NOTE: Strangeness is not required to be conserved in final state
121 G4bool okay() const { return (energyOkay() && momentumOkay() &&
122 baryonOkay() && chargeOkay()); }
123
124 // Calculations of conserved quantities from initial and final state
125 // FIXME: Relative comparisons don't work for zero!
126 G4double deltaE() const { return (final.e() - initial.e()); }
128 return ( (std::abs(deltaE())<tolerance) ? 0. :
129 (initial.e()<tolerance) ? 1. : deltaE()/initial.e() );
130 }
131
132 G4double deltaKE() const { return (ekin(final) - ekin(initial)); }
134 return ( (std::abs(deltaKE())<tolerance) ? 0. :
135 (ekin(initial)<tolerance) ? 1. : deltaKE()/ekin(initial) );
136 }
137
138 G4double deltaP() const { return deltaLV().rho(); }
140 return ( (std::abs(deltaP())<tolerance) ? 0. :
141 (initial.rho()<tolerance) ? 1. : deltaP()/initial.rho() );
142 }
143
144 G4LorentzVector deltaLV() const { return final - initial; }
145
146 // Baryon number, charge, S are discrete; no bounds and no "relative" scale
147 G4int deltaB() const { return (finalBaryon - initialBaryon); }
148 G4int deltaQ() const { return (finalCharge - initialCharge); }
149 G4int deltaS() const { return (finalStrange- initialStrange); }
150
151protected:
152 // Utility function for kinetic energy
153 G4double ekin(const G4LorentzVector& p) const { return (p.e() - p.m()); }
154
155private:
156 G4double relativeLimit; // Fractional bound on conservation
157 G4double absoluteLimit; // Absolute (GeV) bound on conservation
158
159 G4LorentzVector initial; // Four-vectors for computing violations
160 G4LorentzVector final;
161
162 G4int initialBaryon; // Total baryon number
163 G4int finalBaryon;
164
165 G4int initialCharge; // Total charge
166 G4int finalCharge;
167
168 G4int initialStrange; // Total strangeness (s-quark content)
169 G4int finalStrange;
170
171 G4CollisionOutput tempOutput; // Buffer for direct-list interfaces
172
173private:
174 // Copying of modules is forbidden
177};
178
179#endif /* G4CASCADE_CHECK_BALANCE_HH */
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void setRelativeLimit(G4double limit)
void setLimits(G4double relative, G4double absolute)
void setOwner(const G4String &owner)
G4double ekin(const G4LorentzVector &p) const
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
G4LorentzVector deltaLV() const
static const G4double tolerance
void setAbsoluteLimit(G4double limit)
virtual void setName(const G4String &name)