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