Geant4  10.02.p02
G4Fragment.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 // $Id: G4Fragment.hh 97617 2016-06-06 12:47:17Z 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  friend std::ostream& operator<<(std::ostream&, const G4Fragment&);
96 
97  // new/delete operators are overloded to use G4Allocator
98  inline void *operator new(size_t);
99  inline void operator delete(void *aFragment);
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  // 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 const G4ParticleDefinition * GetParticleDefinition() const;
155  inline void SetParticleDefinition(const G4ParticleDefinition * p);
156 
157  inline G4double GetCreationTime() const;
158  inline void SetCreationTime(G4double time);
159 
162 
165 
166  // ============= PRIVATE METHODS ==============================
167 
168 private:
169 
171 
172  void NumberOfExitationWarning(const G4String&);
173 
174  inline void CalculateExcitationEnergy();
175 
176  inline void CalculateGroundStateMass();
177 
178  // ============= DATA MEMBERS ==================
179 
181 
183 
185 
187 
189 
190  // Nuclear polarisation by default is nullptr
192 
193  // creator model type
195 
196  // Exciton model data members
201 
202  // Gamma evaporation data members
204 
206 
209 };
210 
211 // ============= INLINE METHOD IMPLEMENTATIONS ===================
212 
214 {
215  if(p != thePolarization) {
216  delete thePolarization;
217  thePolarization = p;
218  }
219 }
220 
221 #if defined G4HADRONIC_ALLOC_EXPORT
223 #else
225 #endif
226 
227 inline void * G4Fragment::operator new(size_t)
228 {
230  return (void*) pFragmentAllocator->MallocSingle();
231 }
232 
233 inline void G4Fragment::operator delete(void * aFragment)
234 {
235  ((G4Fragment *)aFragment)->SetNuclearPolarization(nullptr);
236  pFragmentAllocator->FreeSingle((G4Fragment *) aFragment);
237 }
238 
240 {
243 }
244 
245 inline G4double
247 {
249 }
250 
252 {
254 }
255 
257 {
258  return theA;
259 }
260 
262 {
263  return theZ;
264 }
265 
266 inline void G4Fragment::SetZandA_asInt(G4int Znew, G4int Anew)
267 {
268  theZ = Znew;
269  theA = Anew;
271 }
272 
274 {
275  return theExcitationEnergy;
276 }
277 
279 {
280  return theGroundStateMass;
281 }
282 
284 {
285  return (theA-theZ)*CLHEP::neutron_mass_c2 + theZ*CLHEP::proton_mass_c2
287 }
288 
290 {
291  return theMomentum;
292 }
293 
294 inline void G4Fragment::SetMomentum(const G4LorentzVector& value)
295 {
296  theMomentum = value;
298 }
299 
301 {
302  return G4double(theZ);
303 }
304 
306 {
307  return G4double(theA);
308 }
309 
310 inline void G4Fragment::SetZ(const G4double value)
311 {
312  theZ = G4lrint(value);
314 }
315 
316 inline void G4Fragment::SetA(const G4double value)
317 {
318  theA = G4lrint(value);
320 }
321 
323 {
325 }
326 
328 {
329  return numberOfParticles;
330 }
331 
333 {
334  return numberOfCharged;
335 }
336 
337 inline
339 {
340  numberOfParticles = valueTot;
341  numberOfCharged = valueP;
342  if(valueTot < valueP) {
343  NumberOfExitationWarning("SetNumberOfExcitedParticle");
344  }
345 }
346 
348 {
349  return numberOfHoles;
350 }
351 
353 {
354  return numberOfChargedHoles;
355 }
356 
357 inline void G4Fragment::SetNumberOfHoles(G4int valueTot, G4int valueP)
358 {
359  numberOfHoles = valueTot;
360  numberOfChargedHoles = valueP;
361  if(valueTot < valueP) {
362  NumberOfExitationWarning("SetNumberOfHoles");
363  }
364 }
365 
367 {
368  numberOfParticles = value;
369 }
370 
372 {
373  numberOfCharged = value;
374  if(value > numberOfParticles) {
375  NumberOfExitationWarning("SetNumberOfCharged");
376  }
377 }
378 
380 {
381  return numberOfShellElectrons;
382 }
383 
385 {
386  numberOfShellElectrons = value;
387 }
388 
390 {
391  return creatorModel;
392 }
393 
395 {
396  creatorModel = value;
397 }
398 
400 {
401  return spin;
402 }
403 
404 inline void G4Fragment::SetSpin(G4double value)
405 {
406  spin = value;
407 }
408 
409 inline
411 {
412  return theParticleDefinition;
413 }
414 
416 {
418 }
419 
421 {
422  return theCreationTime;
423 }
424 
426 {
427  theCreationTime = time;
428 }
429 
431 {
432  return thePolarization;
433 }
434 
435 #endif
436 
437 
G4int theZ
Definition: G4Fragment.hh:182
G4double theCreationTime
Definition: G4Fragment.hh:208
G4double theGroundStateMass
Definition: G4Fragment.hh:186
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int numberOfCharged
Definition: G4Fragment.hh:198
G4NuclearPolarization * thePolarization
Definition: G4Fragment.hh:191
#define G4DLLEXPORT
Definition: G4Types.hh:62
friend std::ostream & operator<<(std::ostream &, const G4Fragment *)
Definition: G4Fragment.cc:181
CLHEP::Hep3Vector G4ThreeVector
void ExcitationEnergyWarning()
Definition: G4Fragment.cc:241
G4double GetA() const
Definition: G4Fragment.hh:305
#define G4DLLIMPORT
Definition: G4Types.hh:63
G4int GetNumberOfElectrons() const
Definition: G4Fragment.hh:379
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:410
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:357
void SetNumberOfElectrons(G4int value)
Definition: G4Fragment.hh:384
G4int GetNumberOfParticles() const
Definition: G4Fragment.hh:327
G4ThreeVector GetAngularMomentum() const
Definition: G4Fragment.cc:275
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
void SetA(G4double value)
Definition: G4Fragment.hh:316
void SetCreatorModelType(G4int value)
Definition: G4Fragment.hh:394
G4int numberOfShellElectrons
Definition: G4Fragment.hh:203
G4int GetNumberOfHoles() const
Definition: G4Fragment.hh:347
void SetParticleDefinition(const G4ParticleDefinition *p)
Definition: G4Fragment.hh:415
void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
Definition: G4Fragment.hh:338
double A(double temperature)
G4int GetA_asInt() const
Definition: G4Fragment.hh:256
G4bool operator==(const G4Fragment &right) const
Definition: G4Fragment.cc:171
void SetNuclearPolarization(G4NuclearPolarization *)
Definition: G4Fragment.hh:213
G4double GetCreationTime() const
Definition: G4Fragment.hh:420
G4DLLIMPORT G4ThreadLocal G4Allocator< G4Fragment > * pFragmentAllocator
Definition: G4Fragment.cc:48
bool G4bool
Definition: G4Types.hh:79
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:289
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:294
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:63
G4double spin
Definition: G4Fragment.hh:207
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:278
void CalculateExcitationEnergy()
Definition: G4Fragment.hh:239
void SetNumberOfParticles(G4int value)
Definition: G4Fragment.hh:366
G4int theA
Definition: G4Fragment.hh:180
G4LorentzVector theMomentum
Definition: G4Fragment.hh:188
G4double GetBindingEnergy() const
Definition: G4Fragment.hh:283
G4double GetSpin() const
Definition: G4Fragment.hh:399
G4int numberOfHoles
Definition: G4Fragment.hh:199
G4int creatorModel
Definition: G4Fragment.hh:194
G4int numberOfParticles
Definition: G4Fragment.hh:197
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:425
G4int GetNumberOfExcitons() const
Definition: G4Fragment.hh:322
int G4lrint(double ad)
Definition: templates.hh:163
void SetZandA_asInt(G4int Znew, G4int Anew)
Definition: G4Fragment.hh:266
G4bool operator!=(const G4Fragment &right) const
Definition: G4Fragment.cc:176
G4int GetZ_asInt() const
Definition: G4Fragment.hh:261
void SetAngularMomentum(G4ThreeVector &)
Definition: G4Fragment.cc:270
void SetZ(G4double value)
Definition: G4Fragment.hh:310
void SetSpin(G4double value)
Definition: G4Fragment.hh:404
void NumberOfExitationWarning(const G4String &)
Definition: G4Fragment.cc:261
G4int GetNumberOfChargedHoles() const
Definition: G4Fragment.hh:352
void SetNumberOfCharged(G4int value)
Definition: G4Fragment.hh:371
double G4double
Definition: G4Types.hh:76
void CalculateGroundStateMass()
Definition: G4Fragment.hh:251
G4Fragment & operator=(const G4Fragment &right)
Definition: G4Fragment.cc:145
G4int GetCreatorModelType() const
Definition: G4Fragment.hh:389
G4int GetNumberOfCharged() const
Definition: G4Fragment.hh:332
G4double ComputeGroundStateMass(G4int Z, G4int A) const
Definition: G4Fragment.hh:246
G4int numberOfChargedHoles
Definition: G4Fragment.hh:200
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:273
G4NuclearPolarization * GetNuclearPolarization() const
Definition: G4Fragment.hh:430
G4double GetZ() const
Definition: G4Fragment.hh:300
const G4ParticleDefinition * theParticleDefinition
Definition: G4Fragment.hh:205
G4double theExcitationEnergy
Definition: G4Fragment.hh:184
CLHEP::HepLorentzVector G4LorentzVector