Geant4  10.01.p03
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 85824 2014-11-05 15:26: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 <vector>
52 
53 #include "globals.hh"
54 #include "G4Allocator.hh"
55 #include "G4LorentzVector.hh"
56 #include "G4ThreeVector.hh"
57 #include "G4NucleiProperties.hh"
58 #include "Randomize.hh"
59 #include "G4Proton.hh"
60 #include "G4Neutron.hh"
61 #include "G4HadTmpUtil.hh"
62 
64 
65 class G4Fragment;
66 typedef std::vector<G4Fragment*> G4FragmentVector;
67 
68 class G4Fragment
69 {
70 public:
71 
72  // ============= CONSTRUCTORS ==================
73 
74  // Default constructor - obsolete
75  G4Fragment();
76 
77  // Destructor
78  ~G4Fragment();
79 
80  // Copy constructor
81  G4Fragment(const G4Fragment &right);
82 
83  // A,Z and 4-momentum - main constructor for fragment
84  G4Fragment(G4int A, G4int Z, const G4LorentzVector& aMomentum);
85 
86  // 4-momentum and pointer to G4particleDefinition (for gammas, e-)
87  G4Fragment(const G4LorentzVector& aMomentum,
88  const G4ParticleDefinition* aParticleDefinition);
89 
90  // ============= OPERATORS ==================
91 
92  G4Fragment & operator=(const G4Fragment &right);
93  G4bool operator==(const G4Fragment &right) const;
94  G4bool operator!=(const G4Fragment &right) const;
95 
96  friend std::ostream& operator<<(std::ostream&, const G4Fragment*);
97  friend std::ostream& operator<<(std::ostream&, const G4Fragment&);
98 
99  // new/delete operators are overloded to use G4Allocator
100  inline void *operator new(size_t);
101  inline void operator delete(void *aFragment);
102 
103  // ============= GENERAL METHODS ==================
104 
105  inline G4int GetZ_asInt() const;
106  inline G4int GetA_asInt() const;
107  inline void SetZandA_asInt(G4int Znew, G4int Anew);
108 
109  inline G4double GetExcitationEnergy() const;
110 
111  inline G4double GetGroundStateMass() const;
112 
113  inline G4double GetBindingEnergy() const;
114 
115  inline const G4LorentzVector& GetMomentum() const;
116  inline void SetMomentum(const G4LorentzVector& value);
117 
118  inline const G4ThreeVector& GetAngularMomentum() const;
119  inline void SetAngularMomentum(const G4ThreeVector& value);
120 
121  // computation of mass for any Z and A
122  inline G4double ComputeGroundStateMass(G4int Z, G4int A) const;
123 
124  inline G4int GetCreatorModelType() const;
125  inline void SetCreatorModelType(G4int value);
126 
127  // obsolete methods
128 
129  inline G4double GetZ() const;
130  inline G4double GetA() const;
131  inline void SetZ(G4double value);
132  inline void SetA(G4double value);
133 
134  // ============= METHODS FOR PRE-COMPOUND MODEL ===============
135 
136  inline G4int GetNumberOfExcitons() const;
137 
138  inline G4int GetNumberOfParticles() const;
139  inline G4int GetNumberOfCharged() const;
140  inline void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP);
141 
142  inline G4int GetNumberOfHoles() const;
143  inline G4int GetNumberOfChargedHoles() const;
144  inline void SetNumberOfHoles(G4int valueTot, G4int valueP=0);
145 
146  // these methods will be removed in future
147  inline void SetNumberOfParticles(G4int value);
148  inline void SetNumberOfCharged(G4int value);
149 
150  // ============= METHODS FOR PHOTON EVAPORATION ===============
151 
152  inline G4int GetNumberOfElectrons() const;
153  inline void SetNumberOfElectrons(G4int value);
154 
155  inline const G4ParticleDefinition * GetParticleDefinition() const;
156  inline void SetParticleDefinition(const G4ParticleDefinition * p);
157 
158  inline G4double GetCreationTime() const;
159  inline void SetCreationTime(G4double time);
160 
161  inline G4bool IsStable() const;
162  inline void SetStable(G4bool val);
163 
164  // ============= PRIVATE METHODS ==============================
165 
166 private:
167 
169 
170  void NumberOfExitationWarning(const G4String&);
171 
172  inline void CalculateExcitationEnergy();
173 
174  inline void CalculateGroundStateMass();
175 
176  // ============= DATA MEMBERS ==================
177 
179 
181 
183 
185 
187 
189 
191 
192  // Exciton model data members
194 
196 
198 
200 
201  // Gamma evaporation data members
203 
205 
207 
209 
210 };
211 
212 // ============= INLINE METHOD IMPLEMENTATIONS ===================
213 
214 #if defined G4HADRONIC_ALLOC_EXPORT
216 #else
218 #endif
219 
220 inline void * G4Fragment::operator new(size_t)
221 {
223  return (void*) pFragmentAllocator->MallocSingle();
224 }
225 
226 inline void G4Fragment::operator delete(void * aFragment)
227 {
228  pFragmentAllocator->FreeSingle((G4Fragment *) aFragment);
229 }
230 
232 {
235  if(theExcitationEnergy < CLHEP::keV) { isStable = true; }
236 }
237 
239 {
241 }
242 
244 {
245  return theA;
246 }
247 
249 {
250  return theZ;
251 }
252 
253 inline void G4Fragment::SetZandA_asInt(G4int Znew, G4int Anew)
254 {
255  theZ = Znew;
256  theA = Anew;
258 }
259 
261 {
262  return theExcitationEnergy;
263 }
264 
266 {
267  return theGroundStateMass;
268 }
269 
271 {
272  return (theA-theZ)*CLHEP::neutron_mass_c2 + theZ*CLHEP::proton_mass_c2
274 }
275 
277 {
278  return theMomentum;
279 }
280 
281 inline void G4Fragment::SetMomentum(const G4LorentzVector& value)
282 {
283  theMomentum = value;
285 }
286 
288 {
289  return theAngularMomentum;
290 }
291 
293 {
294  theAngularMomentum = value;
295 }
296 
297 inline G4double
299 {
301 }
302 
304 {
305  return G4double(theZ);
306 }
307 
309 {
310  return G4double(theA);
311 }
312 
313 inline void G4Fragment::SetZ(const G4double value)
314 {
315  theZ = G4lrint(value);
317 }
318 
319 inline void G4Fragment::SetA(const G4double value)
320 {
321  theA = G4lrint(value);
323 }
324 
326 {
328 }
329 
331 {
332  return numberOfParticles;
333 }
334 
336 {
337  return numberOfCharged;
338 }
339 
340 inline
342 {
343  numberOfParticles = valueTot;
344  numberOfCharged = valueP;
345  if(valueTot < valueP) {
346  NumberOfExitationWarning("SetNumberOfExcitedParticle");
347  }
348 }
349 
351 {
352  return numberOfHoles;
353 }
354 
356 {
357  return numberOfChargedHoles;
358 }
359 
360 inline void G4Fragment::SetNumberOfHoles(G4int valueTot, G4int valueP)
361 {
362  numberOfHoles = valueTot;
363  numberOfChargedHoles = valueP;
364  if(valueTot < valueP) {
365  NumberOfExitationWarning("SetNumberOfHoles");
366  }
367 }
368 
370 {
371  numberOfParticles = value;
372 }
373 
375 {
376  numberOfCharged = value;
377  if(value > numberOfParticles) {
378  NumberOfExitationWarning("SetNumberOfCharged");
379  }
380 }
381 
383 {
384  return numberOfShellElectrons;
385 }
386 
388 {
389  numberOfShellElectrons = value;
390 }
391 
393 {
394  return creatorModel;
395 }
396 
398 {
399  creatorModel = value;
400 }
401 
402 inline
404 {
405  return theParticleDefinition;
406 }
407 
409 {
411 }
412 
414 {
415  return theCreationTime;
416 }
417 
419 {
420  theCreationTime = time;
421 }
422 
424 {
425  return isStable;
426 }
427 
429 {
430  isStable = val;
431 }
432 
433 #endif
434 
435 
G4int theZ
Definition: G4Fragment.hh:180
G4double theCreationTime
Definition: G4Fragment.hh:206
G4double theGroundStateMass
Definition: G4Fragment.hh:184
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int numberOfCharged
Definition: G4Fragment.hh:195
#define G4DLLEXPORT
Definition: G4Types.hh:62
friend std::ostream & operator<<(std::ostream &, const G4Fragment *)
Definition: G4Fragment.cc:180
CLHEP::Hep3Vector G4ThreeVector
void ExcitationEnergyWarning()
Definition: G4Fragment.cc:233
G4double GetA() const
Definition: G4Fragment.hh:308
#define G4DLLIMPORT
Definition: G4Types.hh:63
const G4ThreeVector & GetAngularMomentum() const
Definition: G4Fragment.hh:287
G4int GetNumberOfElectrons() const
Definition: G4Fragment.hh:382
void SetAngularMomentum(const G4ThreeVector &value)
Definition: G4Fragment.hh:292
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:403
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:360
void SetNumberOfElectrons(G4int value)
Definition: G4Fragment.hh:387
G4int GetNumberOfParticles() const
Definition: G4Fragment.hh:330
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
void SetA(G4double value)
Definition: G4Fragment.hh:319
void SetCreatorModelType(G4int value)
Definition: G4Fragment.hh:397
G4int numberOfShellElectrons
Definition: G4Fragment.hh:202
G4int GetNumberOfHoles() const
Definition: G4Fragment.hh:350
void SetParticleDefinition(const G4ParticleDefinition *p)
Definition: G4Fragment.hh:408
void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
Definition: G4Fragment.hh:341
G4int GetA_asInt() const
Definition: G4Fragment.hh:243
G4bool operator==(const G4Fragment &right) const
Definition: G4Fragment.cc:170
G4double GetCreationTime() const
Definition: G4Fragment.hh:413
G4DLLIMPORT G4ThreadLocal G4Allocator< G4Fragment > * pFragmentAllocator
Definition: G4Fragment.cc:50
bool G4bool
Definition: G4Types.hh:79
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:276
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:281
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:265
void CalculateExcitationEnergy()
Definition: G4Fragment.hh:231
void SetNumberOfParticles(G4int value)
Definition: G4Fragment.hh:369
G4int theA
Definition: G4Fragment.hh:178
G4ThreeVector theAngularMomentum
Definition: G4Fragment.hh:188
static const G4double A[nN]
G4LorentzVector theMomentum
Definition: G4Fragment.hh:186
G4double GetBindingEnergy() const
Definition: G4Fragment.hh:270
G4int numberOfHoles
Definition: G4Fragment.hh:197
G4int creatorModel
Definition: G4Fragment.hh:190
G4int numberOfParticles
Definition: G4Fragment.hh:193
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:418
G4bool IsStable() const
Definition: G4Fragment.hh:423
G4int GetNumberOfExcitons() const
Definition: G4Fragment.hh:325
int G4lrint(double ad)
Definition: templates.hh:163
void SetZandA_asInt(G4int Znew, G4int Anew)
Definition: G4Fragment.hh:253
G4bool operator!=(const G4Fragment &right) const
Definition: G4Fragment.cc:175
G4int GetZ_asInt() const
Definition: G4Fragment.hh:248
void SetStable(G4bool val)
Definition: G4Fragment.hh:428
void SetZ(G4double value)
Definition: G4Fragment.hh:313
static const double keV
Definition: G4SIunits.hh:195
void NumberOfExitationWarning(const G4String &)
Definition: G4Fragment.cc:252
G4int GetNumberOfChargedHoles() const
Definition: G4Fragment.hh:355
void SetNumberOfCharged(G4int value)
Definition: G4Fragment.hh:374
double G4double
Definition: G4Types.hh:76
void CalculateGroundStateMass()
Definition: G4Fragment.hh:238
G4Fragment & operator=(const G4Fragment &right)
Definition: G4Fragment.cc:148
G4int GetCreatorModelType() const
Definition: G4Fragment.hh:392
G4int GetNumberOfCharged() const
Definition: G4Fragment.hh:335
G4double ComputeGroundStateMass(G4int Z, G4int A) const
Definition: G4Fragment.hh:298
G4bool isStable
Definition: G4Fragment.hh:208
G4int numberOfChargedHoles
Definition: G4Fragment.hh:199
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:260
G4double GetZ() const
Definition: G4Fragment.hh:303
const G4ParticleDefinition * theParticleDefinition
Definition: G4Fragment.hh:204
G4double theExcitationEnergy
Definition: G4Fragment.hh:182
CLHEP::HepLorentzVector G4LorentzVector