Geant4-11
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//
27//---------------------------------------------------------------------
28//
29// Geant4 header G4Fragment
30//
31// by V. Lara (May 1998)
32//
33// Modifications:
34// 03.05.2010 V.Ivanchenko General cleanup of inline functions: objects
35// are accessed by reference; remove double return
36// tolerance of excitation energy at modent it is computed;
37// safe computation of excitation for exotic fragments
38// 18.05.2010 V.Ivanchenko added member theGroundStateMass and inline
39// method which allowing to compute this value once and use
40// many times
41// 26.09.2010 V.Ivanchenko added number of protons, neutrons, proton holes
42// and neutron holes as members of the class and Get/Set methods;
43// removed not needed 'const'; removed old debug staff and unused
44// private methods; add comments and reorder methods for
45// better reading
46// 27.10.2021 A.Ribon extension for hypernuclei.
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"
56#include "G4NucleiProperties.hh"
58#include "G4Proton.hh"
59#include "G4Neutron.hh"
60#include <vector>
61
63
64class G4Fragment;
65typedef std::vector<G4Fragment*> G4FragmentVector;
66
68{
69public:
70
71 // ============= CONSTRUCTORS ==================
72
73 // Default constructor - obsolete
74 G4Fragment();
75
76 // Destructor
78
79 // Copy constructor
80 G4Fragment(const G4Fragment &right);
81
82 // A,Z and 4-momentum - main constructor for fragment
83 G4Fragment(G4int A, G4int Z, const G4LorentzVector& aMomentum, G4bool warning=true);
84
85 // A,Z,numberOfLambdas and 4-momentum
86 G4Fragment(G4int A, G4int Z, G4int numberOfLambdas, const G4LorentzVector& aMomentum, G4bool warning=true);
87
88 // 4-momentum and pointer to G4particleDefinition (for gammas, e-)
89 G4Fragment(const G4LorentzVector& aMomentum,
90 const G4ParticleDefinition* aParticleDefinition);
91
92 // ============= OPERATORS ==================
93
94 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
100 // new/delete operators are overloded to use G4Allocator
101 inline void *operator new(size_t);
102 inline void operator delete(void *aFragment);
103
104 // ============= GENERAL METHODS ==================
105
106 inline G4int GetZ_asInt() const;
107 inline G4int GetA_asInt() const;
108 inline void SetZandA_asInt(G4int Znew, G4int Anew);
109
110 inline G4int GetNumberOfLambdas() const;
111 inline void SetNumberOfLambdas(G4int numberOfLambdas);
112 // Non-negative number of lambdas/anti-lambdas inside the nucleus/anti-nucleus
113
114 inline G4double GetExcitationEnergy() const;
115 inline void SetExcEnergyAndMomentum(G4double eexc, const G4LorentzVector&);
116
117 inline G4double GetGroundStateMass() const;
118
119 inline G4double GetBindingEnergy() const;
120
121 inline const G4LorentzVector& GetMomentum() const;
122 inline void SetMomentum(const G4LorentzVector& value);
123
124 // computation of mass for any Z, A and numberOfLambdas
125 inline G4double ComputeGroundStateMass(G4int Z, G4int A, G4int numberOfLambdas = 0) const;
126
127 // extra methods
128 inline G4double GetSpin() const;
129 inline void SetSpin(G4double value);
130
131 inline G4int GetCreatorModelID() const;
132 inline void SetCreatorModelID(G4int value);
133
134 // obsolete methods
135
136 inline G4double GetZ() const;
137 inline G4double GetA() const;
138 inline void SetZ(G4double value);
139 inline void SetA(G4double value);
140
141 // ============= METHODS FOR PRE-COMPOUND MODEL ===============
142
143 inline G4int GetNumberOfExcitons() const;
144
145 inline G4int GetNumberOfParticles() const;
146 inline G4int GetNumberOfCharged() const;
147 inline void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP);
148
149 inline G4int GetNumberOfHoles() const;
150 inline G4int GetNumberOfChargedHoles() const;
151 inline void SetNumberOfHoles(G4int valueTot, G4int valueP=0);
152
153 // these methods will be removed in future
154 inline void SetNumberOfParticles(G4int value);
155 inline void SetNumberOfCharged(G4int value);
156
157 // ============= METHODS FOR PHOTON EVAPORATION ===============
158
159 inline G4int GetNumberOfElectrons() const;
160 inline void SetNumberOfElectrons(G4int value);
161
162 inline G4int GetFloatingLevelNumber() const;
163 inline void SetFloatingLevelNumber(G4int value);
164
165 inline const G4ParticleDefinition * GetParticleDefinition() const;
166 inline void SetParticleDefinition(const G4ParticleDefinition * p);
167
168 inline G4double GetCreationTime() const;
169 inline void SetCreationTime(G4double time);
170
171 // G4Fragment class is not responsible for creation and delition of
172 // G4NuclearPolarization object
176
179
180 // ============= PRIVATE METHODS ==============================
181
182private:
183
185
187
188 inline void CalculateExcitationEnergy(G4bool warning=true);
189
190 inline void CalculateGroundStateMass();
191
192 // ============= DATA MEMBERS ==================
193
195
197
198 G4int theL; // Non-negative number of lambdas/anti-lambdas inside the nucleus/anti-nucleus
199
201
203
205
206 // Nuclear polarisation by default is nullptr
208
209 // creator model type
211
212 // Exciton model data members
217
218 // Gamma evaporation data members
221
223
226
228};
229
230// ============= INLINE METHOD IMPLEMENTATIONS ===================
231
232#if defined G4HADRONIC_ALLOC_EXPORT
234#else
236#endif
237
238inline void * G4Fragment::operator new(size_t)
239{
241 return (void*) pFragmentAllocator()->MallocSingle();
242}
243
244inline void G4Fragment::operator delete(void * aFragment)
245{
246 pFragmentAllocator()->FreeSingle((G4Fragment *) aFragment);
247}
248
250{
255 }
256}
257
258inline G4double
260{
261 if ( numberOfLambdas <= 0 ) return G4NucleiProperties::GetNuclearMass(A, Z);
262 else return G4HyperNucleiProperties::GetNuclearMass(A, Z, numberOfLambdas);
263}
264
266{
269}
270
272{
273 return theA;
274}
275
277{
278 return theZ;
279}
280
282{
283 theZ = Znew;
284 theA = Anew;
286}
287
289{
290 return theL;
291}
292
293inline void G4Fragment::SetNumberOfLambdas(G4int numberOfLambdas)
294{
295 theL = std::max( numberOfLambdas, 0 ); // Cannot be negative
297}
298
300{
301 return theExcitationEnergy;
302}
303
305{
306 return theGroundStateMass;
307}
308
310 const G4LorentzVector& v)
311{
312 theExcitationEnergy = eexc;
313 theMomentum.set(0.0, 0.0, 0.0, theGroundStateMass + eexc);
315}
316
318{
321}
322
324{
325 return theMomentum;
326}
327
329{
330 theMomentum = value;
332}
333
335{
336 return G4double(theZ);
337}
338
340{
341 return G4double(theA);
342}
343
344inline void G4Fragment::SetZ(const G4double value)
345{
346 theZ = G4lrint(value);
348}
349
350inline void G4Fragment::SetA(const G4double value)
351{
352 theA = G4lrint(value);
354}
355
357{
359}
360
362{
363 return numberOfParticles;
364}
365
367{
368 return numberOfCharged;
369}
370
371inline
373{
374 numberOfParticles = valueTot;
375 numberOfCharged = valueP;
376 if(valueTot < valueP) {
377 NumberOfExitationWarning("SetNumberOfExcitedParticle");
378 }
379}
380
382{
383 return numberOfHoles;
384}
385
387{
389}
390
391inline void G4Fragment::SetNumberOfHoles(G4int valueTot, G4int valueP)
392{
393 numberOfHoles = valueTot;
394 numberOfChargedHoles = valueP;
395 if(valueTot < valueP) {
396 NumberOfExitationWarning("SetNumberOfHoles");
397 }
398}
399
401{
402 numberOfParticles = value;
403}
404
406{
407 numberOfCharged = value;
408 if(value > numberOfParticles) {
409 NumberOfExitationWarning("SetNumberOfCharged");
410 }
411}
412
414{
416}
417
419{
421}
422
424{
425 return creatorModel;
426}
427
429{
430 creatorModel = value;
431}
432
434{
435 return spin;
436}
437
439{
440 spin = value;
441}
442
444{
445 return xLevel;
446}
447
449{
450 xLevel = value;
451}
452
453inline
455{
457}
458
460{
462}
463
465{
466 return theCreationTime;
467}
468
470{
471 theCreationTime = time;
472}
473
475{
476 return thePolarization;
477}
478
480{
481 return thePolarization;
482}
483
485{
486 thePolarization = p;
487}
488
489#endif
490
491
G4DLLIMPORT G4Allocator< G4Fragment > *& pFragmentAllocator()
Definition: G4Fragment.cc:46
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:64
double G4double
Definition: G4Types.hh:83
#define G4DLLIMPORT
Definition: G4Types.hh:69
#define G4DLLEXPORT
Definition: G4Types.hh:68
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
void set(double x, double y, double z, double t)
G4double theCreationTime
Definition: G4Fragment.hh:225
G4double theGroundStateMass
Definition: G4Fragment.hh:202
void SetZ(G4double value)
Definition: G4Fragment.hh:344
void SetFloatingLevelNumber(G4int value)
Definition: G4Fragment.hh:448
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:304
void SetNumberOfLambdas(G4int numberOfLambdas)
Definition: G4Fragment.hh:293
void ExcitationEnergyWarning()
Definition: G4Fragment.cc:264
G4int GetNumberOfParticles() const
Definition: G4Fragment.hh:361
G4int GetCreatorModelID() const
Definition: G4Fragment.hh:423
G4int GetNumberOfHoles() const
Definition: G4Fragment.hh:381
G4NuclearPolarization * GetNuclearPolarization() const
Definition: G4Fragment.hh:479
G4double ComputeGroundStateMass(G4int Z, G4int A, G4int numberOfLambdas=0) const
Definition: G4Fragment.hh:259
void NumberOfExitationWarning(const G4String &)
Definition: G4Fragment.cc:272
void SetNumberOfCharged(G4int value)
Definition: G4Fragment.hh:405
G4int numberOfChargedHoles
Definition: G4Fragment.hh:216
friend std::ostream & operator<<(std::ostream &, const G4Fragment &)
Definition: G4Fragment.cc:212
G4int GetNumberOfChargedHoles() const
Definition: G4Fragment.hh:386
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:299
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:323
void SetCreatorModelID(G4int value)
Definition: G4Fragment.hh:428
G4double GetCreationTime() const
Definition: G4Fragment.hh:464
void SetNuclearPolarization(G4NuclearPolarization *)
Definition: G4Fragment.hh:484
G4double GetSpin() const
Definition: G4Fragment.hh:433
G4double GetBindingEnergy() const
Definition: G4Fragment.hh:317
G4int theZ
Definition: G4Fragment.hh:196
G4NuclearPolarization * thePolarization
Definition: G4Fragment.hh:207
G4double GetZ() const
Definition: G4Fragment.hh:334
G4bool operator!=(const G4Fragment &right) const
Definition: G4Fragment.cc:207
G4int GetZ_asInt() const
Definition: G4Fragment.hh:276
G4double GetA() const
Definition: G4Fragment.hh:339
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:469
G4int GetFloatingLevelNumber() const
Definition: G4Fragment.hh:443
void SetNumberOfElectrons(G4int value)
Definition: G4Fragment.hh:418
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:391
G4int GetNumberOfExcitons() const
Definition: G4Fragment.hh:356
G4int GetNumberOfLambdas() const
Definition: G4Fragment.hh:288
G4int numberOfParticles
Definition: G4Fragment.hh:213
G4ThreeVector GetAngularMomentum() const
Definition: G4Fragment.cc:286
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:454
void CalculateExcitationEnergy(G4bool warning=true)
Definition: G4Fragment.hh:249
void SetA(G4double value)
Definition: G4Fragment.hh:350
G4int numberOfHoles
Definition: G4Fragment.hh:215
G4int theA
Definition: G4Fragment.hh:194
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:328
const G4ParticleDefinition * theParticleDefinition
Definition: G4Fragment.hh:222
void SetExcEnergyAndMomentum(G4double eexc, const G4LorentzVector &)
Definition: G4Fragment.hh:309
G4int GetNumberOfElectrons() const
Definition: G4Fragment.hh:413
void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
Definition: G4Fragment.hh:372
G4int theL
Definition: G4Fragment.hh:198
G4double theExcitationEnergy
Definition: G4Fragment.hh:200
G4Fragment & operator=(const G4Fragment &right)
Definition: G4Fragment.cc:178
G4int creatorModel
Definition: G4Fragment.hh:210
void SetNumberOfParticles(G4int value)
Definition: G4Fragment.hh:400
G4int numberOfCharged
Definition: G4Fragment.hh:214
G4int GetNumberOfCharged() const
Definition: G4Fragment.hh:366
G4int numberOfShellElectrons
Definition: G4Fragment.hh:219
G4bool operator==(const G4Fragment &right) const
Definition: G4Fragment.cc:202
void SetAngularMomentum(const G4ThreeVector &)
Definition: G4Fragment.cc:281
G4double spin
Definition: G4Fragment.hh:224
void SetParticleDefinition(const G4ParticleDefinition *p)
Definition: G4Fragment.hh:459
void SetZandA_asInt(G4int Znew, G4int Anew)
Definition: G4Fragment.hh:281
G4LorentzVector theMomentum
Definition: G4Fragment.hh:204
G4NuclearPolarization * NuclearPolarization()
Definition: G4Fragment.hh:474
void SetSpin(G4double value)
Definition: G4Fragment.hh:438
void CalculateGroundStateMass()
Definition: G4Fragment.hh:265
G4int xLevel
Definition: G4Fragment.hh:220
G4int GetA_asInt() const
Definition: G4Fragment.hh:271
static const G4double minFragExcitation
Definition: G4Fragment.hh:227
static G4double GetNuclearMass(G4int A, G4int Z, G4int L)
static G4double GetNuclearMass(const G4double A, const G4double Z)
static constexpr double proton_mass_c2
static constexpr double neutron_mass_c2
T max(const T t1, const T t2)
brief Return the largest of the two arguments
int G4lrint(double ad)
Definition: templates.hh:134