Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4CascadeCheckBalance.hh
이 파일의 문서화 페이지로 가기
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: G4CascadeCheckBalance.hh 71942 2013-06-28 19:08:11Z mkelsey $
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 // 20130620 Address Coverity complaint about missing copy actions
50 // 20130621 Add interface to take G4Fragment input instead of G4InuclNuclei.
51 // 20140930 Change name from "const char*" to "const G4String"
52 
53 #include "G4VCascadeCollider.hh"
54 #include "globals.hh"
55 #include "G4CollisionOutput.hh"
56 #include "G4LorentzVector.hh"
57 #include <cmath>
58 #include <vector>
59 
60 class G4CascadParticle;
62 class G4InuclNuclei;
63 class G4InuclParticle;
64 
66 public:
67  static const G4double tolerance; // Don't do floating zero!
68 
69  explicit G4CascadeCheckBalance(const G4String& owner="G4CascadeCheckBalance");
70 
71  G4CascadeCheckBalance(G4double relative, G4double absolute,
72  const G4String& owner="G4CascadeCheckBalance");
73  virtual ~G4CascadeCheckBalance() {};
74 
75  void setOwner(const G4String& owner) { setName(owner); }
76 
77  void setLimits(G4double relative, G4double absolute) {
78  setRelativeLimit(relative);
79  setAbsoluteLimit(absolute);
80  }
81 
82  void setRelativeLimit(G4double limit) { relativeLimit = limit; }
83  void setAbsoluteLimit(G4double limit) { absoluteLimit = limit; }
84 
86  G4CollisionOutput& output);
87 
88  // This is for use with G4VCascadeDeexcitation modules
89  void collide(const G4Fragment& fragment, G4CollisionOutput& output);
90 
91  // This is for use with G4EPCollider internal checks
93  const std::vector<G4InuclElementaryParticle>& particles);
94 
95  // This is for use with G4NucleiModel internal checks
97  const std::vector<G4CascadParticle>& particles);
98 
99  // This is for use with G4IntraNucleiCascader
101  G4CollisionOutput& output,
102  const std::vector<G4CascadParticle>& cparticles);
103 
104  // This is for use with G4BigBanger internal checks
105  void collide(const G4Fragment& target,
106  const std::vector<G4InuclElementaryParticle>& particles);
107 
108  // This is for use with G4Fissioner internal checks
109  void collide(const G4Fragment& target,
110  const std::vector<G4InuclNuclei>& fragments);
111 
112  // Checks on conservation laws (kinematics, baryon number, charge, hyperons)
113  G4bool energyOkay() const;
114  G4bool ekinOkay() const;
115  G4bool momentumOkay() const;
116  G4bool baryonOkay() const;
117  G4bool chargeOkay() const;
118  G4bool strangeOkay() const;
119 
120  // Global check, used by G4CascadeInterface validation loop
121  // NOTE: Strangeness is not required to be conserved in final state
122  G4bool okay() const { return (energyOkay() && momentumOkay() &&
123  baryonOkay() && chargeOkay()); }
124 
125  // Calculations of conserved quantities from initial and final state
126  // FIXME: Relative comparisons don't work for zero!
127  G4double deltaE() const { return (final.e() - initial.e()); }
128  G4double relativeE() const {
129  return ( (std::abs(deltaE())<tolerance) ? 0. :
130  (initial.e()<tolerance) ? 1. : deltaE()/initial.e() );
131  }
132 
133  G4double deltaKE() const { return (ekin(final) - ekin(initial)); }
135  return ( (std::abs(deltaKE())<tolerance) ? 0. :
136  (ekin(initial)<tolerance) ? 1. : deltaKE()/ekin(initial) );
137  }
138 
139  G4double deltaP() const { return deltaLV().rho(); }
140  G4double relativeP() const {
141  return ( (std::abs(deltaP())<tolerance) ? 0. :
142  (initial.rho()<tolerance) ? 1. : deltaP()/initial.rho() );
143  }
144 
145  G4LorentzVector deltaLV() const { return final - initial; }
146 
147  // Baryon number, charge, S are discrete; no bounds and no "relative" scale
148  G4int deltaB() const { return (finalBaryon - initialBaryon); }
149  G4int deltaQ() const { return (finalCharge - initialCharge); }
150  G4int deltaS() const { return (finalStrange- initialStrange); }
151 
152 protected:
153  // Utility function for kinetic energy
154  G4double ekin(const G4LorentzVector& p) const { return (p.e() - p.m()); }
155 
156 private:
157  G4double relativeLimit; // Fractional bound on conservation
158  G4double absoluteLimit; // Absolute (GeV) bound on conservation
159 
160  G4LorentzVector initial; // Four-vectors for computing violations
162 
163  G4int initialBaryon; // Total baryon number
165 
166  G4int initialCharge; // Total charge
168 
169  G4int initialStrange; // Total strangeness (s-quark content)
171 
172  G4CollisionOutput tempOutput; // Buffer for direct-list interfaces
173 
174 private:
175  // Copying of modules is forbidden
178 };
179 
180 #endif /* G4CASCADE_CHECK_BALANCE_HH */
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
virtual void setName(const G4String &name)
const XML_Char * target
Definition: expat.h:268
void setOwner(const G4String &owner)
void setAbsoluteLimit(G4double limit)
const char * p
Definition: xmltok.h:285
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4double ekin(const G4LorentzVector &p) const
void setLimits(G4double relative, G4double absolute)
static const G4double tolerance
G4CascadeCheckBalance(const G4String &owner="G4CascadeCheckBalance")
double rho() const
int G4int
Definition: G4Types.hh:78
void setRelativeLimit(G4double limit)
G4CascadeCheckBalance & operator=(const G4CascadeCheckBalance &)
G4LorentzVector deltaLV() const