source: trunk/source/processes/hadronic/models/util/include/G4Fragment.hh@ 1330

Last change on this file since 1330 was 1315, checked in by garnier, 15 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 9.3 KB
RevLine 
[819]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//
[1315]26// $Id: G4Fragment.hh,v 1.11 2010/05/19 10:23:00 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-cand-01 $
[819]28//
[1315]29//---------------------------------------------------------------------
[819]30//
[1315]31// Geant4 header G4Fragment
32//
[819]33// by V. Lara (May 1998)
[1315]34//
35// Modifications:
36// 03.05.2010 V.Ivanchenko General cleanup of inline functions: objects
37// are accessed by reference; remove double return
38// tolerance of excitation energy at modent it is computed;
39// safe computation of excitation for exotic fragments
40// 18.05.2010 V.Ivanchenko added member theGroundStateMass and inline
41// method which allowing to compute this value once and use
42// many times
[819]43
44#ifndef G4Fragment_h
45#define G4Fragment_h 1
46
47#include "G4ios.hh"
48#include <iomanip>
49#include <vector>
50
51#include "globals.hh"
52#include "G4LorentzVector.hh"
[1315]53//#include "G4ParticleMomentum.hh"
[819]54#include "G4ThreeVector.hh"
55#include "G4NucleiProperties.hh"
56#include "G4ParticleTable.hh"
57#include "G4IonTable.hh"
58#include "Randomize.hh"
59#include "G4Proton.hh"
60#include "G4Neutron.hh"
61#include "G4HadronicException.hh"
62#include "G4HadTmpUtil.hh"
63
64
65class G4ParticleDefinition;
66
67class G4Fragment; // Forward deckaration
68typedef std::vector<G4Fragment*> G4FragmentVector;
69
70class G4Fragment
71{
72public:
73
74 // ============= CONSTRUCTORS ==================
75
[1315]76 // Default constructor - obsolete
[819]77 G4Fragment();
78
79 // Destructor
80 ~G4Fragment();
81
82 // Copy constructor
83 G4Fragment(const G4Fragment &right);
84
[1315]85 // A,Z and 4-momentum - main constructor for fragment
86 G4Fragment(const G4int A, const G4int Z, const G4LorentzVector& aMomentum);
[819]87
88 // 4-momentum and pointer to G4particleDefinition (for gammas)
[1315]89 G4Fragment(const G4LorentzVector& aMomentum, G4ParticleDefinition * aParticleDefinition);
[819]90
91 // ============= OPERATORS ==================
92
93 const G4Fragment & operator=(const G4Fragment &right);
94 G4bool operator==(const G4Fragment &right) const;
95 G4bool operator!=(const G4Fragment &right) const;
96
97 friend std::ostream& operator<<(std::ostream&, const G4Fragment*);
98 friend std::ostream& operator<<(std::ostream&, const G4Fragment&);
99
100 // ============= METHODS ==================
101
[1315]102 inline G4double GetA() const;
103 inline void SetA(const G4double value);
[819]104
[1315]105 inline G4double GetZ() const;
106 inline void SetZ(const G4double value);
107
108 inline G4int GetZ_asInt() const;
109 inline G4int GetA_asInt() const;
110 inline void SetZandA_asInt(G4int Znew, G4int Anew);
[819]111
[1315]112 inline G4double GetExcitationEnergy() const;
[819]113 void SetExcitationEnergy(const G4double value);
114
[1315]115 inline const G4LorentzVector& GetMomentum() const;
116 inline void SetMomentum(const G4LorentzVector& value);
[819]117
[1315]118 inline const G4ThreeVector& GetAngularMomentum() const;
119 inline void SetAngularMomentum(const G4ThreeVector& value);
[819]120
[1315]121 inline G4int GetNumberOfExcitons() const;
[819]122
[1315]123 inline G4int GetNumberOfHoles() const;
124 inline void SetNumberOfHoles(const G4int value);
[819]125
[1315]126 inline G4int GetNumberOfCharged() const;
[819]127 void SetNumberOfCharged(const G4int value);
128
[1315]129 inline G4int GetNumberOfParticles() const;
130 inline void SetNumberOfParticles(const G4int value);
[819]131
[1315]132 inline G4ParticleDefinition * GetParticleDefinition() const;
133 inline void SetParticleDefinition(G4ParticleDefinition * aParticleDefinition);
[819]134
[1315]135 inline G4double GetCreationTime() const;
136 inline void SetCreationTime(const G4double time);
[819]137
[1315]138 inline G4double GetGroundStateMass() const;
139
140 inline G4double GetBindingEnergy() const;
[819]141
[1315]142 // computation of mass for any Z and A
143 inline G4double ComputeGroundStateMass(const G4int Z, const G4int A) const;
[819]144
145#ifdef PRECOMPOUND_TEST
[1315]146 G4String GetCreatorModel() const;
147 void SetCreatorModel(const G4String & aModel);
[819]148#endif
149
150private:
151
[1315]152 void ExcitationEnegryWarning();
[819]153
[1315]154 inline void CalculateGroundStateMass();
155
156 inline void CalculateExcitationEnergy();
157
[819]158 G4ThreeVector IsotropicRandom3Vector(const G4double Magnitude = 1.0) const;
159
160 // ============= DATA MEMBERS ==================
161
[1315]162 static G4int errCount;
163
164 G4int theA;
[819]165
[1315]166 G4int theZ;
[819]167
168 G4double theExcitationEnergy;
169
[1315]170 G4double theGroundStateMass;
171
[819]172 G4LorentzVector theMomentum;
173
174 G4ThreeVector theAngularMomentum;
175
176 G4int numberOfParticles;
177
178 G4int numberOfHoles;
179
180 G4int numberOfCharged;
181
182 // Gamma evaporation requeriments
183
184 G4ParticleDefinition * theParticleDefinition;
185
186 G4double theCreationTime;
187
188#ifdef PRECOMPOUND_TEST
189 G4String theCreatorModel;
190#endif
191};
192
193// Class G4Fragment
[1315]194inline void G4Fragment::CalculateGroundStateMass()
195{
196 theGroundStateMass = G4NucleiProperties::GetNuclearMass(theA, theZ);
197}
[819]198
199inline G4double G4Fragment::GetA() const
200{
[1315]201 return G4double(theA);
[819]202}
203
204inline void G4Fragment::SetA(const G4double value)
205{
[1315]206 theA = G4lrint(value);
207 CalculateGroundStateMass();
[819]208}
209
210inline G4double G4Fragment::GetZ() const
211{
[1315]212 return G4double(theZ);
[819]213}
214
215inline void G4Fragment::SetZ(const G4double value)
216{
[1315]217 theZ = G4lrint(value);
218 CalculateGroundStateMass();
[819]219}
220
[1315]221inline G4int G4Fragment::GetA_asInt() const
222{
223 return theA;
224}
225
226inline G4int G4Fragment::GetZ_asInt() const
227{
228 return theZ;
229}
230
231inline void G4Fragment::SetZandA_asInt(G4int Znew, G4int Anew)
232{
233 theZ = Znew;
234 theA = Anew;
235 CalculateGroundStateMass();
236}
237
[819]238inline G4double G4Fragment::GetExcitationEnergy() const
239{
240 // temporary fix for what seems to be
241 // a problem with rounding errors for on-shell lorentz-vectors in CLHEP.
242 // HPW Apr 1999 @@@@@@@
243
[1315]244 //VI if(std::abs(theExcitationEnergy)<10*eV) return 0;
[819]245 return theExcitationEnergy;
246}
247
[1315]248inline const G4LorentzVector& G4Fragment::GetMomentum() const
[819]249{
250 return theMomentum;
251}
252
[1315]253inline const G4ThreeVector& G4Fragment::GetAngularMomentum() const
[819]254{
255 return theAngularMomentum;
256}
257
[1315]258inline void G4Fragment::SetAngularMomentum(const G4ThreeVector& value)
[819]259{
260 theAngularMomentum = value;
261}
262
263inline G4int G4Fragment::GetNumberOfExcitons() const
264{
265 return numberOfParticles + numberOfHoles;
266}
267
268inline void G4Fragment::SetNumberOfParticles(const G4int value)
269{
270 numberOfParticles = value;
271}
272
273inline G4int G4Fragment::GetNumberOfHoles() const
274{
275 return numberOfHoles;
276}
277
278inline void G4Fragment::SetNumberOfHoles(const G4int value)
279{
280 numberOfHoles = value;
281}
282
283inline G4int G4Fragment::GetNumberOfCharged() const
284{
285 return numberOfCharged;
286}
287
288inline void G4Fragment::SetNumberOfCharged(const G4int value)
289{
[1315]290 if (value <= numberOfParticles) { numberOfCharged = value; }
[819]291 else
292 {
[1315]293 G4String text = "G4Fragment::SetNumberOfCharged: Number of charged particles can't be greater than number of particles";
294 throw G4HadronicException(__FILE__, __LINE__, text);
[819]295 }
296}
297
298inline G4int G4Fragment::GetNumberOfParticles() const
299{
300 return numberOfParticles;
301}
302
303inline G4ParticleDefinition * G4Fragment::GetParticleDefinition(void) const
304{
305 return theParticleDefinition;
306}
307
308inline void G4Fragment::SetParticleDefinition(G4ParticleDefinition * aParticleDefinition)
309{
310 theParticleDefinition = aParticleDefinition;
311}
312
[1315]313inline G4double G4Fragment::GetCreationTime() const
[819]314{
315 return theCreationTime;
316}
317
318inline void G4Fragment::SetCreationTime(const G4double time)
319{
320 theCreationTime = time;
321}
322
[1315]323inline G4double G4Fragment::GetGroundStateMass() const
[819]324{
[1315]325 return theGroundStateMass;
[819]326}
[1315]327
328inline G4double
329G4Fragment::ComputeGroundStateMass(const G4int Z, const G4int A) const
330{
331 return G4NucleiProperties::GetNuclearMass(A, Z);
332}
333
334inline void G4Fragment::CalculateExcitationEnergy()
335{
336 theExcitationEnergy = theMomentum.mag() - theGroundStateMass;
337 if(theExcitationEnergy < 0.0) { ExcitationEnegryWarning(); }
338}
[819]339
[1315]340inline G4double G4Fragment::GetBindingEnergy() const
[819]341{
[1315]342 return (theA-theZ)*CLHEP::neutron_mass_c2 + theZ*CLHEP::proton_mass_c2
343 - theGroundStateMass;
[819]344}
345
[1315]346inline void G4Fragment::SetMomentum(const G4LorentzVector& value)
347{
348 theMomentum = value;
349 CalculateExcitationEnergy();
350}
[819]351
352#endif
353
354
Note: See TracBrowser for help on using the repository browser.