Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4Fragment.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 // $Id: G4Fragment.hh 110267 2018-05-17 14:34:00Z gcosmo $
27 //
28 //---------------------------------------------------------------------
29 //
30 // Geant4 header G4Fragment
31 //
32 // by V. Lara (May 1998)
33 //
34 // Modifications:
35 // 03.05.2010 V.Ivanchenko General cleanup of inline functions: objects
36 // are accessed by reference; remove double return
37 // tolerance of excitation energy at modent it is computed;
38 // safe computation of excitation for exotic fragments
39 // 18.05.2010 V.Ivanchenko added member theGroundStateMass and inline
40 // method which allowing to compute this value once and use
41 // many times
42 // 26.09.2010 V.Ivanchenko added number of protons, neutrons, proton holes
43 // and neutron holes as members of the class and Get/Set methods;
44 // removed not needed 'const'; removed old debug staff and unused
45 // private methods; add comments and reorder methods for
46 // better reading
47 
48 #ifndef G4Fragment_h
49 #define G4Fragment_h 1
50 
51 #include "globals.hh"
52 #include "G4Allocator.hh"
53 #include "G4LorentzVector.hh"
54 #include "G4ThreeVector.hh"
55 #include "G4NuclearPolarization.hh"
56 #include "G4NucleiProperties.hh"
57 #include "G4Proton.hh"
58 #include "G4Neutron.hh"
59 #include <vector>
60 
62 
63 class G4Fragment;
64 typedef std::vector<G4Fragment*> G4FragmentVector;
65 
66 class G4Fragment
67 {
68 public:
69 
70  // ============= CONSTRUCTORS ==================
71 
72  // Default constructor - obsolete
73  G4Fragment();
74 
75  // Destructor
76  ~G4Fragment();
77 
78  // Copy constructor
79  G4Fragment(const G4Fragment &right);
80 
81  // A,Z and 4-momentum - main constructor for fragment
82  G4Fragment(G4int A, G4int Z, const G4LorentzVector& aMomentum);
83 
84  // 4-momentum and pointer to G4particleDefinition (for gammas, e-)
85  G4Fragment(const G4LorentzVector& aMomentum,
86  const G4ParticleDefinition* aParticleDefinition);
87 
88  // ============= OPERATORS ==================
89 
90  G4Fragment & operator=(const G4Fragment &right);
91  G4bool operator==(const G4Fragment &right) const;
92  G4bool operator!=(const G4Fragment &right) const;
93 
94  friend std::ostream& operator<<(std::ostream&, const G4Fragment&);
95 
96  // new/delete operators are overloded to use G4Allocator
97  inline void *operator new(size_t);
98  inline void operator delete(void *aFragment);
99 
100  // ============= GENERAL METHODS ==================
101 
102  inline G4int GetZ_asInt() const;
103  inline G4int GetA_asInt() const;
104  inline void SetZandA_asInt(G4int Znew, G4int Anew);
105 
106  inline G4double GetExcitationEnergy() const;
107  inline void SetExcEnergyAndMomentum(G4double eexc, const G4LorentzVector&);
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  // computation of mass for any Z and A
117  inline G4double ComputeGroundStateMass(G4int Z, G4int A) const;
118 
119  // extra methods
120  inline G4double GetSpin() const;
121  inline void SetSpin(G4double value);
122 
123  inline G4int GetCreatorModelType() const;
124  inline void SetCreatorModelType(G4int value);
125 
126  // obsolete methods
127 
128  inline G4double GetZ() const;
129  inline G4double GetA() const;
130  inline void SetZ(G4double value);
131  inline void SetA(G4double value);
132 
133  // ============= METHODS FOR PRE-COMPOUND MODEL ===============
134 
135  inline G4int GetNumberOfExcitons() const;
136 
137  inline G4int GetNumberOfParticles() const;
138  inline G4int GetNumberOfCharged() const;
139  inline void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP);
140 
141  inline G4int GetNumberOfHoles() const;
142  inline G4int GetNumberOfChargedHoles() const;
143  inline void SetNumberOfHoles(G4int valueTot, G4int valueP=0);
144 
145  // these methods will be removed in future
146  inline void SetNumberOfParticles(G4int value);
147  inline void SetNumberOfCharged(G4int value);
148 
149  // ============= METHODS FOR PHOTON EVAPORATION ===============
150 
151  inline G4int GetNumberOfElectrons() const;
152  inline void SetNumberOfElectrons(G4int value);
153 
154  inline G4int GetFloatingLevelNumber() const;
155  inline void SetFloatingLevelNumber(G4int value);
156 
157  inline const G4ParticleDefinition * GetParticleDefinition() const;
158  inline void SetParticleDefinition(const G4ParticleDefinition * p);
159 
160  inline G4double GetCreationTime() const;
161  inline void SetCreationTime(G4double time);
162 
163  // G4Fragment class is not responsible for creation and delition of
164  // G4NuclearPolarization object
168 
169  void SetAngularMomentum(const G4ThreeVector&);
171 
172  // ============= PRIVATE METHODS ==============================
173 
174 private:
175 
177 
178  void NumberOfExitationWarning(const G4String&);
179 
180  inline void CalculateExcitationEnergy();
181 
182  inline void CalculateGroundStateMass();
183 
184  // ============= DATA MEMBERS ==================
185 
187 
189 
191 
193 
195 
196  // Nuclear polarisation by default is nullptr
198 
199  // creator model type
201 
202  // Exciton model data members
207 
208  // Gamma evaporation data members
211 
213 
216 
218 };
219 
220 // ============= INLINE METHOD IMPLEMENTATIONS ===================
221 
222 #if defined G4HADRONIC_ALLOC_EXPORT
224 #else
226 #endif
227 
228 inline void * G4Fragment::operator new(size_t)
229 {
231  return (void*) pFragmentAllocator()->MallocSingle();
232 }
233 
234 inline void G4Fragment::operator delete(void * aFragment)
235 {
236  pFragmentAllocator()->FreeSingle((G4Fragment *) aFragment);
237 }
238 
240 {
244  theExcitationEnergy = 0.0;
245  }
246 }
247 
248 inline G4double
250 {
252 }
253 
255 {
257 }
258 
260 {
261  return theA;
262 }
263 
265 {
266  return theZ;
267 }
268 
269 inline void G4Fragment::SetZandA_asInt(G4int Znew, G4int Anew)
270 {
271  theZ = Znew;
272  theA = Anew;
274 }
275 
277 {
278  return theExcitationEnergy;
279 }
280 
282 {
283  return theGroundStateMass;
284 }
285 
287  const G4LorentzVector& v)
288 {
289  theExcitationEnergy = eexc;
290  theMomentum.set(0.0, 0.0, 0.0, theGroundStateMass + eexc);
292 }
293 
295 {
298 }
299 
301 {
302  return theMomentum;
303 }
304 
306 {
307  theMomentum = value;
309 }
310 
312 {
313  return G4double(theZ);
314 }
315 
317 {
318  return G4double(theA);
319 }
320 
321 inline void G4Fragment::SetZ(const G4double value)
322 {
323  theZ = G4lrint(value);
325 }
326 
327 inline void G4Fragment::SetA(const G4double value)
328 {
329  theA = G4lrint(value);
331 }
332 
334 {
336 }
337 
339 {
340  return numberOfParticles;
341 }
342 
344 {
345  return numberOfCharged;
346 }
347 
348 inline
350 {
351  numberOfParticles = valueTot;
352  numberOfCharged = valueP;
353  if(valueTot < valueP) {
354  NumberOfExitationWarning("SetNumberOfExcitedParticle");
355  }
356 }
357 
359 {
360  return numberOfHoles;
361 }
362 
364 {
365  return numberOfChargedHoles;
366 }
367 
368 inline void G4Fragment::SetNumberOfHoles(G4int valueTot, G4int valueP)
369 {
370  numberOfHoles = valueTot;
371  numberOfChargedHoles = valueP;
372  if(valueTot < valueP) {
373  NumberOfExitationWarning("SetNumberOfHoles");
374  }
375 }
376 
378 {
380 }
381 
383 {
385  if(value > numberOfParticles) {
386  NumberOfExitationWarning("SetNumberOfCharged");
387  }
388 }
389 
391 {
392  return numberOfShellElectrons;
393 }
394 
396 {
398 }
399 
401 {
402  return creatorModel;
403 }
404 
406 {
408 }
409 
411 {
412  return spin;
413 }
414 
416 {
417  spin = value;
418 }
419 
421 {
422  return xLevel;
423 }
424 
426 {
427  xLevel = value;
428 }
429 
430 inline
432 {
433  return theParticleDefinition;
434 }
435 
437 {
439 }
440 
442 {
443  return theCreationTime;
444 }
445 
447 {
448  theCreationTime = time;
449 }
450 
452 {
453  return thePolarization;
454 }
455 
457 {
458  return thePolarization;
459 }
460 
462 {
463  thePolarization = p;
464 }
465 
466 #endif
467 
468 
G4double GetA() const
Definition: G4Fragment.hh:316
G4double GetCreationTime() const
Definition: G4Fragment.hh:441
static const G4double minFragExcitation
Definition: G4Fragment.hh:217
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:300
void SetNuclearPolarization(G4NuclearPolarization *)
Definition: G4Fragment.hh:461
void SetZ(G4double value)
Definition: G4Fragment.hh:321
G4double theExcitationEnergy
Definition: G4Fragment.hh:190
G4LorentzVector theMomentum
Definition: G4Fragment.hh:194
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:368
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:446
G4NuclearPolarization * NuclearPolarization()
Definition: G4Fragment.hh:451
G4int numberOfParticles
Definition: G4Fragment.hh:203
G4NuclearPolarization * GetNuclearPolarization() const
Definition: G4Fragment.hh:456
void SetZandA_asInt(G4int Znew, G4int Anew)
Definition: G4Fragment.hh:269
const char * p
Definition: xmltok.h:285
G4int GetCreatorModelType() const
Definition: G4Fragment.hh:400
G4double GetSpin() const
Definition: G4Fragment.hh:410
G4double GetBindingEnergy() const
Definition: G4Fragment.hh:294
G4int GetNumberOfHoles() const
Definition: G4Fragment.hh:358
G4int GetA_asInt() const
Definition: G4Fragment.hh:259
G4int theA
Definition: G4Fragment.hh:186
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:281
G4int GetNumberOfParticles() const
Definition: G4Fragment.hh:338
G4int numberOfChargedHoles
Definition: G4Fragment.hh:206
void NumberOfExitationWarning(const G4String &)
Definition: G4Fragment.cc:241
G4int theZ
Definition: G4Fragment.hh:188
void CalculateExcitationEnergy()
Definition: G4Fragment.hh:239
static constexpr double neutron_mass_c2
G4double theCreationTime
Definition: G4Fragment.hh:215
static constexpr double proton_mass_c2
Float_t Z
void SetExcEnergyAndMomentum(G4double eexc, const G4LorentzVector &)
Definition: G4Fragment.hh:286
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4NuclearPolarization * thePolarization
Definition: G4Fragment.hh:197
const XML_Char int const XML_Char * value
Definition: expat.h:331
G4double GetZ() const
Definition: G4Fragment.hh:311
void SetSpin(G4double value)
Definition: G4Fragment.hh:415
G4double spin
Definition: G4Fragment.hh:214
G4int GetNumberOfCharged() const
Definition: G4Fragment.hh:343
G4int numberOfHoles
Definition: G4Fragment.hh:205
G4double ComputeGroundStateMass(G4int Z, G4int A) const
Definition: G4Fragment.hh:249
void CalculateGroundStateMass()
Definition: G4Fragment.hh:254
G4Fragment & operator=(const G4Fragment &right)
Definition: G4Fragment.cc:149
friend std::ostream & operator<<(std::ostream &, const G4Fragment &)
Definition: G4Fragment.cc:182
double A(double temperature)
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:305
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int creatorModel
Definition: G4Fragment.hh:200
G4DLLIMPORT G4Allocator< G4Fragment > *& pFragmentAllocator()
Definition: G4Fragment.cc:46
void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
Definition: G4Fragment.hh:349
void SetNumberOfCharged(G4int value)
Definition: G4Fragment.hh:382
void set(double x, double y, double z, double t)
void ExcitationEnergyWarning()
Definition: G4Fragment.cc:233
G4int GetZ_asInt() const
Definition: G4Fragment.hh:264
void SetFloatingLevelNumber(G4int value)
Definition: G4Fragment.hh:425
double mag() const
void SetCreatorModelType(G4int value)
Definition: G4Fragment.hh:405
int G4lrint(double ad)
Definition: templates.hh:151
int G4int
Definition: G4Types.hh:78
void SetA(G4double value)
Definition: G4Fragment.hh:327
G4int GetNumberOfElectrons() const
Definition: G4Fragment.hh:390
G4int numberOfShellElectrons
Definition: G4Fragment.hh:209
const G4ParticleDefinition * theParticleDefinition
Definition: G4Fragment.hh:212
void SetNumberOfParticles(G4int value)
Definition: G4Fragment.hh:377
G4int xLevel
Definition: G4Fragment.hh:210
G4bool operator==(const G4Fragment &right) const
Definition: G4Fragment.cc:172
#define G4DLLEXPORT
Definition: G4Types.hh:62
Hep3Vector boostVector() const
G4int GetFloatingLevelNumber() const
Definition: G4Fragment.hh:420
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:431
G4int GetNumberOfChargedHoles() const
Definition: G4Fragment.hh:363
G4int GetNumberOfExcitons() const
Definition: G4Fragment.hh:333
G4double theGroundStateMass
Definition: G4Fragment.hh:192
void SetAngularMomentum(const G4ThreeVector &)
Definition: G4Fragment.cc:250
#define G4DLLIMPORT
Definition: G4Types.hh:63
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:276
G4ThreeVector GetAngularMomentum() const
Definition: G4Fragment.cc:255
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:63
void SetNumberOfElectrons(G4int value)
Definition: G4Fragment.hh:395
G4int numberOfCharged
Definition: G4Fragment.hh:204
HepLorentzVector & boost(double, double, double)
void SetParticleDefinition(const G4ParticleDefinition *p)
Definition: G4Fragment.hh:436
G4bool operator!=(const G4Fragment &right) const
Definition: G4Fragment.cc:177