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 | // $Id: G4Fragment.hh,v 1.16 2010/09/28 16:09:00 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-03-ref-09 $ |
---|
28 | // |
---|
29 | //--------------------------------------------------------------------- |
---|
30 | // |
---|
31 | // Geant4 header G4Fragment |
---|
32 | // |
---|
33 | // by V. Lara (May 1998) |
---|
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 |
---|
43 | // 26.09.2010 V.Ivanchenko added number of protons, neutrons, proton holes |
---|
44 | // and neutron holes as members of the class and Get/Set methods; |
---|
45 | // removed not needed 'const'; removed old debug staff and unused |
---|
46 | // private methods; add comments and reorder methods for |
---|
47 | // better reading |
---|
48 | |
---|
49 | #ifndef G4Fragment_h |
---|
50 | #define G4Fragment_h 1 |
---|
51 | |
---|
52 | #include <vector> |
---|
53 | |
---|
54 | #include "globals.hh" |
---|
55 | #include "G4LorentzVector.hh" |
---|
56 | #include "G4ThreeVector.hh" |
---|
57 | #include "G4NucleiProperties.hh" |
---|
58 | //#include "G4ParticleTable.hh" |
---|
59 | //#include "G4IonTable.hh" |
---|
60 | #include "Randomize.hh" |
---|
61 | #include "G4Proton.hh" |
---|
62 | #include "G4Neutron.hh" |
---|
63 | #include "G4HadTmpUtil.hh" |
---|
64 | |
---|
65 | class G4ParticleDefinition; |
---|
66 | |
---|
67 | class G4Fragment; |
---|
68 | typedef std::vector<G4Fragment*> G4FragmentVector; |
---|
69 | |
---|
70 | class G4Fragment |
---|
71 | { |
---|
72 | public: |
---|
73 | |
---|
74 | // ============= CONSTRUCTORS ================== |
---|
75 | |
---|
76 | // Default constructor - obsolete |
---|
77 | G4Fragment(); |
---|
78 | |
---|
79 | // Destructor |
---|
80 | ~G4Fragment(); |
---|
81 | |
---|
82 | // Copy constructor |
---|
83 | G4Fragment(const G4Fragment &right); |
---|
84 | |
---|
85 | // A,Z and 4-momentum - main constructor for fragment |
---|
86 | G4Fragment(G4int A, G4int Z, const G4LorentzVector& aMomentum); |
---|
87 | |
---|
88 | // 4-momentum and pointer to G4particleDefinition (for gammas, e-) |
---|
89 | G4Fragment(const G4LorentzVector& aMomentum, |
---|
90 | G4ParticleDefinition* aParticleDefinition); |
---|
91 | |
---|
92 | // ============= OPERATORS ================== |
---|
93 | |
---|
94 | const G4Fragment & operator=(const G4Fragment &right); |
---|
95 | G4bool operator==(const G4Fragment &right) const; |
---|
96 | G4bool operator!=(const G4Fragment &right) const; |
---|
97 | |
---|
98 | friend std::ostream& operator<<(std::ostream&, const G4Fragment*); |
---|
99 | friend std::ostream& operator<<(std::ostream&, const G4Fragment&); |
---|
100 | |
---|
101 | // ============= GENERAL METHODS ================== |
---|
102 | |
---|
103 | inline G4int GetZ_asInt() const; |
---|
104 | inline G4int GetA_asInt() const; |
---|
105 | inline void SetZandA_asInt(G4int Znew, G4int Anew); |
---|
106 | |
---|
107 | inline G4double GetExcitationEnergy() const; |
---|
108 | |
---|
109 | inline G4double GetGroundStateMass() const; |
---|
110 | |
---|
111 | inline G4double GetBindingEnergy() const; |
---|
112 | |
---|
113 | inline const G4LorentzVector& GetMomentum() const; |
---|
114 | inline void SetMomentum(const G4LorentzVector& value); |
---|
115 | |
---|
116 | inline const G4ThreeVector& GetAngularMomentum() const; |
---|
117 | inline void SetAngularMomentum(const G4ThreeVector& value); |
---|
118 | |
---|
119 | // computation of mass for any Z and A |
---|
120 | inline G4double ComputeGroundStateMass(G4int Z, G4int A) const; |
---|
121 | |
---|
122 | // obsolete methods |
---|
123 | inline G4double GetZ() const; |
---|
124 | inline G4double GetA() const; |
---|
125 | inline void SetZ(G4double value); |
---|
126 | inline void SetA(G4double value); |
---|
127 | |
---|
128 | // ============= METHODS FOR PRE-COMPOUND MODEL =============== |
---|
129 | |
---|
130 | inline G4int GetNumberOfExcitons() const; |
---|
131 | |
---|
132 | inline G4int GetNumberOfParticles() const; |
---|
133 | inline G4int GetNumberOfCharged() const; |
---|
134 | inline void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP); |
---|
135 | |
---|
136 | inline G4int GetNumberOfHoles() const; |
---|
137 | inline G4int GetNumberOfChargedHoles() const; |
---|
138 | inline void SetNumberOfHoles(G4int valueTot, G4int valueP=0); |
---|
139 | |
---|
140 | // these methods will be removed in future |
---|
141 | inline void SetNumberOfParticles(G4int value); |
---|
142 | inline void SetNumberOfCharged(G4int value); |
---|
143 | |
---|
144 | // ============= METHODS FOR PHOTON EVAPORATION =============== |
---|
145 | |
---|
146 | inline G4int GetNumberOfElectrons() const; |
---|
147 | inline void SetNumberOfElectrons(G4int value); |
---|
148 | |
---|
149 | inline G4ParticleDefinition * GetParticleDefinition() const; |
---|
150 | inline void SetParticleDefinition(G4ParticleDefinition * p); |
---|
151 | |
---|
152 | inline G4double GetCreationTime() const; |
---|
153 | inline void SetCreationTime(G4double time); |
---|
154 | |
---|
155 | // ============= PRIVATE METHODS ============================== |
---|
156 | |
---|
157 | private: |
---|
158 | |
---|
159 | void ExcitationEnergyWarning(); |
---|
160 | |
---|
161 | void NumberOfExitationWarning(const G4String&); |
---|
162 | |
---|
163 | inline void CalculateExcitationEnergy(); |
---|
164 | |
---|
165 | inline void CalculateGroundStateMass(); |
---|
166 | |
---|
167 | // ============= DATA MEMBERS ================== |
---|
168 | |
---|
169 | static G4int errCount; |
---|
170 | |
---|
171 | G4int theA; |
---|
172 | |
---|
173 | G4int theZ; |
---|
174 | |
---|
175 | G4double theExcitationEnergy; |
---|
176 | |
---|
177 | G4double theGroundStateMass; |
---|
178 | |
---|
179 | G4LorentzVector theMomentum; |
---|
180 | |
---|
181 | G4ThreeVector theAngularMomentum; |
---|
182 | |
---|
183 | // Exciton model data members |
---|
184 | |
---|
185 | G4int numberOfParticles; |
---|
186 | |
---|
187 | G4int numberOfCharged; |
---|
188 | |
---|
189 | G4int numberOfHoles; |
---|
190 | |
---|
191 | G4int numberOfChargedHoles; |
---|
192 | |
---|
193 | // Gamma evaporation data members |
---|
194 | |
---|
195 | G4int numberOfShellElectrons; |
---|
196 | |
---|
197 | G4ParticleDefinition * theParticleDefinition; |
---|
198 | |
---|
199 | G4double theCreationTime; |
---|
200 | |
---|
201 | }; |
---|
202 | |
---|
203 | // ============= INLINE METHOD IMPLEMENTATIONS =================== |
---|
204 | |
---|
205 | inline void G4Fragment::CalculateExcitationEnergy() |
---|
206 | { |
---|
207 | theExcitationEnergy = theMomentum.mag() - theGroundStateMass; |
---|
208 | if(theExcitationEnergy < 0.0) { ExcitationEnergyWarning(); } |
---|
209 | } |
---|
210 | |
---|
211 | inline void G4Fragment::CalculateGroundStateMass() |
---|
212 | { |
---|
213 | theGroundStateMass = G4NucleiProperties::GetNuclearMass(theA, theZ); |
---|
214 | } |
---|
215 | |
---|
216 | inline G4int G4Fragment::GetA_asInt() const |
---|
217 | { |
---|
218 | return theA; |
---|
219 | } |
---|
220 | |
---|
221 | inline G4int G4Fragment::GetZ_asInt() const |
---|
222 | { |
---|
223 | return theZ; |
---|
224 | } |
---|
225 | |
---|
226 | inline void G4Fragment::SetZandA_asInt(G4int Znew, G4int Anew) |
---|
227 | { |
---|
228 | theZ = Znew; |
---|
229 | theA = Anew; |
---|
230 | CalculateGroundStateMass(); |
---|
231 | } |
---|
232 | |
---|
233 | inline G4double G4Fragment::GetExcitationEnergy() const |
---|
234 | { |
---|
235 | return theExcitationEnergy; |
---|
236 | } |
---|
237 | |
---|
238 | inline G4double G4Fragment::GetGroundStateMass() const |
---|
239 | { |
---|
240 | return theGroundStateMass; |
---|
241 | } |
---|
242 | |
---|
243 | inline G4double G4Fragment::GetBindingEnergy() const |
---|
244 | { |
---|
245 | return (theA-theZ)*CLHEP::neutron_mass_c2 + theZ*CLHEP::proton_mass_c2 |
---|
246 | - theGroundStateMass; |
---|
247 | } |
---|
248 | |
---|
249 | inline const G4LorentzVector& G4Fragment::GetMomentum() const |
---|
250 | { |
---|
251 | return theMomentum; |
---|
252 | } |
---|
253 | |
---|
254 | inline void G4Fragment::SetMomentum(const G4LorentzVector& value) |
---|
255 | { |
---|
256 | theMomentum = value; |
---|
257 | CalculateExcitationEnergy(); |
---|
258 | } |
---|
259 | |
---|
260 | inline const G4ThreeVector& G4Fragment::GetAngularMomentum() const |
---|
261 | { |
---|
262 | return theAngularMomentum; |
---|
263 | } |
---|
264 | |
---|
265 | inline void G4Fragment::SetAngularMomentum(const G4ThreeVector& value) |
---|
266 | { |
---|
267 | theAngularMomentum = value; |
---|
268 | } |
---|
269 | |
---|
270 | inline G4double |
---|
271 | G4Fragment::ComputeGroundStateMass(G4int Z, G4int A) const |
---|
272 | { |
---|
273 | return G4NucleiProperties::GetNuclearMass(A, Z); |
---|
274 | } |
---|
275 | |
---|
276 | inline G4double G4Fragment::GetZ() const |
---|
277 | { |
---|
278 | return G4double(theZ); |
---|
279 | } |
---|
280 | |
---|
281 | inline G4double G4Fragment::GetA() const |
---|
282 | { |
---|
283 | return G4double(theA); |
---|
284 | } |
---|
285 | |
---|
286 | inline void G4Fragment::SetZ(const G4double value) |
---|
287 | { |
---|
288 | theZ = G4lrint(value); |
---|
289 | CalculateGroundStateMass(); |
---|
290 | } |
---|
291 | |
---|
292 | inline void G4Fragment::SetA(const G4double value) |
---|
293 | { |
---|
294 | theA = G4lrint(value); |
---|
295 | CalculateGroundStateMass(); |
---|
296 | } |
---|
297 | |
---|
298 | inline G4int G4Fragment::GetNumberOfExcitons() const |
---|
299 | { |
---|
300 | return numberOfParticles + numberOfHoles; |
---|
301 | } |
---|
302 | |
---|
303 | inline G4int G4Fragment::GetNumberOfParticles() const |
---|
304 | { |
---|
305 | return numberOfParticles; |
---|
306 | } |
---|
307 | |
---|
308 | inline G4int G4Fragment::GetNumberOfCharged() const |
---|
309 | { |
---|
310 | return numberOfCharged; |
---|
311 | } |
---|
312 | |
---|
313 | inline |
---|
314 | void G4Fragment::SetNumberOfExcitedParticle(G4int valueTot, G4int valueP) |
---|
315 | { |
---|
316 | numberOfParticles = valueTot; |
---|
317 | numberOfCharged = valueP; |
---|
318 | if(valueTot < valueP) { |
---|
319 | NumberOfExitationWarning("SetNumberOfExcitedParticle"); |
---|
320 | } |
---|
321 | } |
---|
322 | |
---|
323 | inline G4int G4Fragment::GetNumberOfHoles() const |
---|
324 | { |
---|
325 | return numberOfHoles; |
---|
326 | } |
---|
327 | |
---|
328 | inline G4int G4Fragment::GetNumberOfChargedHoles() const |
---|
329 | { |
---|
330 | return numberOfChargedHoles; |
---|
331 | } |
---|
332 | |
---|
333 | inline void G4Fragment::SetNumberOfHoles(G4int valueTot, G4int valueP) |
---|
334 | { |
---|
335 | numberOfHoles = valueTot; |
---|
336 | numberOfChargedHoles = valueP; |
---|
337 | if(valueTot < valueP) { |
---|
338 | NumberOfExitationWarning("SetNumberOfHoles"); |
---|
339 | } |
---|
340 | } |
---|
341 | |
---|
342 | inline void G4Fragment::SetNumberOfParticles(G4int value) |
---|
343 | { |
---|
344 | numberOfParticles = value; |
---|
345 | } |
---|
346 | |
---|
347 | inline void G4Fragment::SetNumberOfCharged(G4int value) |
---|
348 | { |
---|
349 | numberOfCharged = value; |
---|
350 | if(value > numberOfParticles) { |
---|
351 | NumberOfExitationWarning("SetNumberOfCharged"); |
---|
352 | } |
---|
353 | } |
---|
354 | |
---|
355 | inline G4int G4Fragment::GetNumberOfElectrons() const |
---|
356 | { |
---|
357 | return numberOfShellElectrons; |
---|
358 | } |
---|
359 | |
---|
360 | inline void G4Fragment::SetNumberOfElectrons(G4int value) |
---|
361 | { |
---|
362 | numberOfShellElectrons = value; |
---|
363 | } |
---|
364 | |
---|
365 | inline |
---|
366 | G4ParticleDefinition * G4Fragment::GetParticleDefinition(void) const |
---|
367 | { |
---|
368 | return theParticleDefinition; |
---|
369 | } |
---|
370 | |
---|
371 | inline void G4Fragment::SetParticleDefinition(G4ParticleDefinition * p) |
---|
372 | { |
---|
373 | theParticleDefinition = p; |
---|
374 | } |
---|
375 | |
---|
376 | inline G4double G4Fragment::GetCreationTime() const |
---|
377 | { |
---|
378 | return theCreationTime; |
---|
379 | } |
---|
380 | |
---|
381 | inline void G4Fragment::SetCreationTime(G4double time) |
---|
382 | { |
---|
383 | theCreationTime = time; |
---|
384 | } |
---|
385 | |
---|
386 | #endif |
---|
387 | |
---|
388 | |
---|