Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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);
84
85 // A,Z,numberOfLambdas and 4-momentum
86 G4Fragment(G4int A, G4int Z, G4int numberOfLambdas,
87 const G4LorentzVector& aMomentum);
88
89 // 4-momentum and pointer to G4particleDefinition (for gammas, e-)
90 G4Fragment(const G4LorentzVector& aMomentum,
91 const G4ParticleDefinition* aParticleDefinition);
92
93 // ============= OPERATORS ==================
94
95 G4Fragment & operator=(const G4Fragment &right);
96 G4bool operator==(const G4Fragment &right) const;
97 G4bool operator!=(const G4Fragment &right) const;
98
99 friend std::ostream& operator<<(std::ostream&, const G4Fragment&);
100
101 // new/delete operators are overloded to use G4Allocator
102 inline void *operator new(size_t);
103 inline void operator delete(void *aFragment);
104
105 // ============= GENERAL METHODS ==================
106
107 inline G4int GetZ_asInt() const;
108 inline G4int GetA_asInt() const;
109
110 // update number of nucleons without check on input
111 // ground state mass is not recomputed
112 inline void SetZandA_asInt(G4int Znew, G4int Anew, G4int Lnew=0);
113
114 // non-negative number of lambdas/anti-lambdas
115 // ground state mass is not recomputed
116 void SetNumberOfLambdas(G4int numberOfLambdas);
117
118 inline G4int GetNumberOfLambdas() const;
119
120 inline G4double GetExcitationEnergy() const;
121
122 inline G4double GetGroundStateMass() const;
123
124 inline const G4LorentzVector& GetMomentum() const;
125
126 // ground state mass and excitation energy are recomputed
128
129 // update main fragment parameters full check on input
130 // ground state mass and excitation energy are recomputed
131 inline void SetMomentum(const G4LorentzVector& value);
132 inline void SetZAandMomentum(const G4LorentzVector&,
133 G4int Z, G4int A,
134 G4int nLambdas = 0);
135
136 // ground state mass is not recomputed
139
140 // computation of mass for any imput Z, A and number of Lambdas
141 // no check on input values
143 G4int nLambdas = 0) const;
144
145 // extra methods
146 inline G4double GetSpin() const;
147 inline void SetSpin(G4double value);
148
149 inline G4int GetCreatorModelID() const;
150 inline void SetCreatorModelID(G4int value);
151
152 inline G4bool IsLongLived() const;
153 inline void SetLongLived(G4bool value);
154
155 // obsolete methods
156 inline G4double GetZ() const;
157 inline G4double GetA() const;
158 inline void SetZ(G4double value);
159 inline void SetA(G4double value);
160
161 // ============= METHODS FOR PRE-COMPOUND MODEL ===============
162
163 inline G4int GetNumberOfExcitons() const;
164
165 inline G4int GetNumberOfParticles() const;
166 inline G4int GetNumberOfCharged() const;
167 inline void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP);
168
169 inline G4int GetNumberOfHoles() const;
170 inline G4int GetNumberOfChargedHoles() const;
171 inline void SetNumberOfHoles(G4int valueTot, G4int valueP=0);
172
173 // these methods will be removed in future
174 inline void SetNumberOfParticles(G4int value);
175 inline void SetNumberOfCharged(G4int value);
176
177 // ============= METHODS FOR PHOTON EVAPORATION ===============
178
179 inline G4int GetNumberOfElectrons() const;
180 inline void SetNumberOfElectrons(G4int value);
181
182 inline G4int GetFloatingLevelNumber() const;
183 inline void SetFloatingLevelNumber(G4int value);
184
185 inline const G4ParticleDefinition * GetParticleDefinition() const;
186 inline void SetParticleDefinition(const G4ParticleDefinition * p);
187
188 inline G4double GetCreationTime() const;
189 inline void SetCreationTime(G4double time);
190
191 // G4Fragment class is not responsible for creation and delition of
192 // G4NuclearPolarization object
196
199
200 // ============= PRIVATE METHODS ==============================
201
202private:
203
204 void CalculateMassAndExcitationEnergy();
205
206 void ExcitationEnergyWarning();
207
208 void NumberOfExitationWarning(const G4String&);
209
210 // ============= DATA MEMBERS ==================
211
212 G4int theA;
213
214 G4int theZ;
215
216 // Non-negative number of lambdas/anti-lambdas inside the nucleus/anti-nucleus
217 G4int theL;
218
219 G4double theExcitationEnergy;
220
221 G4double theGroundStateMass;
222
223 G4LorentzVector theMomentum;
224
225 // Nuclear polarisation by default is nullptr
226 G4NuclearPolarization* thePolarization;
227
228 // creator model type
229 G4int creatorModel;
230
231 // Exciton model data members
232 G4int numberOfParticles;
233 G4int numberOfCharged;
234 G4int numberOfHoles;
235 G4int numberOfChargedHoles;
236
237 // Gamma evaporation data members
238 G4int numberOfShellElectrons;
239 G4int xLevel;
240
241 const G4ParticleDefinition* theParticleDefinition;
242
243 G4double spin;
244 G4double theCreationTime;
245
246 G4bool isLongLived = false;
247};
248
249// ============= INLINE METHOD IMPLEMENTATIONS ===================
250
251#if defined G4HADRONIC_ALLOC_EXPORT
253#else
255#endif
256
257inline void * G4Fragment::operator new(size_t)
258{
259 if (!pFragmentAllocator()) {
261 }
262 return (void*) pFragmentAllocator()->MallocSingle();
263}
264
265inline void G4Fragment::operator delete(void * aFragment)
266{
267 pFragmentAllocator()->FreeSingle((G4Fragment *) aFragment);
268}
269
270inline G4double
272{
273 return ( nLambdas <= 0 )
276}
277
279{
280 CalculateMassAndExcitationEnergy();
281 return theGroundStateMass;
282}
283
285{
286 return theA;
287}
288
290{
291 return theZ;
292}
293
294inline void
296{
297 theZ = Znew;
298 theA = Anew;
299 theL = std::max(Lnew, 0);
300}
301
303{
304 theL = std::max(Lnew, 0);
305}
306
308{
309 return theL;
310}
311
313{
314 return theExcitationEnergy;
315}
316
318{
319 return theGroundStateMass;
320}
321
323{
324 return theMomentum;
325}
326
328{
329 theMomentum = value;
330 CalculateMassAndExcitationEnergy();
331}
332
333inline void
335 G4int Z, G4int A, G4int nLambdas)
336{
337 SetZandA_asInt(Z, A, nLambdas);
338 SetMomentum(v);
339}
340
342{
343 return static_cast<G4double>(theZ);
344}
345
347{
348 return static_cast<G4double>(theA);
349}
350
351inline void G4Fragment::SetZ(const G4double value)
352{
353 theZ = G4lrint(value);
354}
355
356inline void G4Fragment::SetA(const G4double value)
357{
358 theA = G4lrint(value);
359}
360
362{
363 return numberOfParticles + numberOfHoles;
364}
365
367{
368 return numberOfParticles;
369}
370
372{
373 return numberOfCharged;
374}
375
376inline
378{
379 numberOfParticles = valueTot;
380 numberOfCharged = valueP;
381 if(valueTot < valueP) {
382 NumberOfExitationWarning("SetNumberOfExcitedParticle");
383 }
384}
385
387{
388 return numberOfHoles;
389}
390
392{
393 return numberOfChargedHoles;
394}
395
396inline void G4Fragment::SetNumberOfHoles(G4int valueTot, G4int valueP)
397{
398 numberOfHoles = valueTot;
399 numberOfChargedHoles = valueP;
400 if(valueTot < valueP) {
401 NumberOfExitationWarning("SetNumberOfHoles");
402 }
403}
404
406{
407 numberOfParticles = value;
408}
409
411{
412 numberOfCharged = value;
413 if(value > numberOfParticles) {
414 NumberOfExitationWarning("SetNumberOfCharged");
415 }
416}
417
419{
420 return numberOfShellElectrons;
421}
422
424{
425 numberOfShellElectrons = value;
426}
427
429{
430 return creatorModel;
431}
432
434{
435 creatorModel = value;
436}
437
439{
440 return spin;
441}
442
444{
445 spin = value;
446}
447
449{
450 return isLongLived;
451}
452
454{
455 isLongLived = value;
456}
457
459{
460 return xLevel;
461}
462
464{
465 xLevel = value;
466}
467
468inline
470{
471 return theParticleDefinition;
472}
473
475{
476 theParticleDefinition = p;
477}
478
480{
481 return theCreationTime;
482}
483
485{
486 theCreationTime = time;
487}
488
490{
491 return thePolarization;
492}
493
495{
496 return thePolarization;
497}
498
500{
501 thePolarization = p;
502}
503
504#endif
G4DLLIMPORT G4Allocator< G4Fragment > *& pFragmentAllocator()
Definition: G4Fragment.cc:47
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
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]
void SetZ(G4double value)
Definition: G4Fragment.hh:351
void SetFloatingLevelNumber(G4int value)
Definition: G4Fragment.hh:463
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:317
void SetNumberOfLambdas(G4int numberOfLambdas)
Definition: G4Fragment.hh:302
G4int GetNumberOfParticles() const
Definition: G4Fragment.hh:366
G4int GetCreatorModelID() const
Definition: G4Fragment.hh:428
G4int GetNumberOfHoles() const
Definition: G4Fragment.hh:386
G4NuclearPolarization * GetNuclearPolarization() const
Definition: G4Fragment.hh:494
void SetZandA_asInt(G4int Znew, G4int Anew, G4int Lnew=0)
Definition: G4Fragment.hh:295
void SetNumberOfCharged(G4int value)
Definition: G4Fragment.hh:410
friend std::ostream & operator<<(std::ostream &, const G4Fragment &)
Definition: G4Fragment.cc:254
G4int GetNumberOfChargedHoles() const
Definition: G4Fragment.hh:391
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:312
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:322
void SetCreatorModelID(G4int value)
Definition: G4Fragment.hh:433
G4double GetCreationTime() const
Definition: G4Fragment.hh:479
void SetNuclearPolarization(G4NuclearPolarization *)
Definition: G4Fragment.hh:499
G4double GetSpin() const
Definition: G4Fragment.hh:438
G4double GetBindingEnergy() const
Definition: G4Fragment.cc:211
G4double GetZ() const
Definition: G4Fragment.hh:341
G4bool IsLongLived() const
Definition: G4Fragment.hh:448
G4bool operator!=(const G4Fragment &right) const
Definition: G4Fragment.cc:249
G4int GetZ_asInt() const
Definition: G4Fragment.hh:289
G4double GetA() const
Definition: G4Fragment.hh:346
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:484
G4int GetFloatingLevelNumber() const
Definition: G4Fragment.hh:458
void SetNumberOfElectrons(G4int value)
Definition: G4Fragment.hh:423
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:396
G4int GetNumberOfExcitons() const
Definition: G4Fragment.hh:361
G4int GetNumberOfLambdas() const
Definition: G4Fragment.hh:307
G4double ComputeGroundStateMass(G4int Z, G4int A, G4int nLambdas=0) const
Definition: G4Fragment.hh:271
G4ThreeVector GetAngularMomentum() const
Definition: G4Fragment.cc:330
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:469
void SetA(G4double value)
Definition: G4Fragment.hh:356
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:327
void SetExcEnergyAndMomentum(G4double eexc, const G4LorentzVector &)
Definition: G4Fragment.cc:203
G4int GetNumberOfElectrons() const
Definition: G4Fragment.hh:418
void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
Definition: G4Fragment.hh:377
G4Fragment & operator=(const G4Fragment &right)
Definition: G4Fragment.cc:219
void SetLongLived(G4bool value)
Definition: G4Fragment.hh:453
G4double RecomputeGroundStateMass()
Definition: G4Fragment.hh:278
void SetNumberOfParticles(G4int value)
Definition: G4Fragment.hh:405
G4int GetNumberOfCharged() const
Definition: G4Fragment.hh:371
G4bool operator==(const G4Fragment &right) const
Definition: G4Fragment.cc:244
void SetAngularMomentum(const G4ThreeVector &)
Definition: G4Fragment.cc:325
void SetParticleDefinition(const G4ParticleDefinition *p)
Definition: G4Fragment.hh:474
G4NuclearPolarization * NuclearPolarization()
Definition: G4Fragment.hh:489
void SetSpin(G4double value)
Definition: G4Fragment.hh:443
G4int GetA_asInt() const
Definition: G4Fragment.hh:284
void SetZAandMomentum(const G4LorentzVector &, G4int Z, G4int A, G4int nLambdas=0)
Definition: G4Fragment.hh:334
static G4double GetNuclearMass(G4int A, G4int Z, G4int L)
static G4double GetNuclearMass(const G4double A, const G4double Z)
int G4lrint(double ad)
Definition: templates.hh:134