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
G4VDecayChannel.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// G4VDecayChannel
27//
28// Class description:
29//
30// Abstract class to describe decay kinematics
31
32// Author: H.Kurashige, 27 July 1996
33// --------------------------------------------------------------------
34#ifndef G4VDecayChannel_hh
35#define G4VDecayChannel_hh 1
36
37#include <cmath>
38
39#include "G4ios.hh"
40#include "globals.hh"
41#include "G4ThreeVector.hh"
42#include "G4Threading.hh"
43#include "G4AutoLock.hh"
44
46class G4DecayProducts;
47class G4ParticleTable;
48
50{
51 public:
52
53 G4VDecayChannel(const G4String& aName, G4int Verbose = 1);
54 G4VDecayChannel(const G4String& aName,
55 const G4String& theParentName,
56 G4double theBR,
57 G4int theNumberOfDaughters,
58 const G4String& theDaughterName1,
59 const G4String& theDaughterName2 = "",
60 const G4String& theDaughterName3 = "",
61 const G4String& theDaughterName4 = "",
62 const G4String& theDaughterName5 = "" );
63 // Constructors
64
65 virtual ~G4VDecayChannel();
66 // Destructor
67
68 G4bool operator==(const G4VDecayChannel& r) const { return (this == &r); }
69 G4bool operator!=(const G4VDecayChannel& r) const { return (this != &r); }
70 // Equality operators
71
72 inline G4bool operator<(const G4VDecayChannel& right) const;
73 // Less-than operator is defined for G4DecayTable
74
75 virtual G4DecayProducts* DecayIt(G4double parentMass = -1.0) = 0;
76
77 inline const G4String& GetKinematicsName() const;
78 // Get kinematics name
79 inline G4double GetBR() const;
80 // Get branching ratio
81 inline G4int GetNumberOfDaughters() const;
82 // Get number of daughter particles
83
85 // Get the pointer to the parent particle
86 inline G4ParticleDefinition* GetDaughter(G4int anIndex);
87 // Get the pointer to a daughter particle
88
90 // Get the angular momentum of the decay
91 inline const G4String& GetParentName() const;
92 // Get the name of the parent particle
93 inline const G4String& GetDaughterName(G4int anIndex) const;
94 // Get the name of a daughter particle
95
96 inline G4double GetParentMass() const;
97 // Get mass of parent
98 inline G4double GetDaughterMass(G4int anIndex) const;
99 // Get mass of daughter
100
101 void SetParent(const G4ParticleDefinition * particle_type);
102 inline void SetParent(const G4String &particle_name);
103 // Set the parent particle (by name or by pointer)
104 void SetBR(G4double value);
105 // Set branching ratio
106 void SetNumberOfDaughters(G4int value);
107 // Set number of daughter particles
108 void SetDaughter(G4int anIndex, const G4ParticleDefinition* particle_type);
109 void SetDaughter(G4int anIndex, const G4String& particle_name);
110 // Set a daughter particle (by name or by pointer)
111
112 inline void SetVerboseLevel(G4int value);
113 inline G4int GetVerboseLevel() const;
114 void DumpInfo();
115
116 inline G4double GetRangeMass() const;
117 inline void SetRangeMass(G4double val);
118 virtual G4bool IsOKWithParentMass(G4double parentMass);
119
120 void SetPolarization(const G4ThreeVector&);
121 inline const G4ThreeVector& GetPolarization() const;
122
123 protected:
124
125 void ClearDaughtersName();
126 // Clear daughters array
127
128 inline void CheckAndFillDaughters();
129 inline void CheckAndFillParent();
130
132 G4double maxDev = 1.0) const;
133
135 // Default constructor
136
139 // Copy constructor and assignment operator
140
141 private:
142
143 void FillDaughters();
144 // Fill daughters array
145 void FillParent();
146 // Fill parent
147
148 const G4String& GetNoName() const;
149
150 protected:
151
153 // Kinematics name
155 // Branching ratio [0.0 - 1.0]
157 // Parent particle
159 // Daughter particles
160
162 // Range of mass allowed in decay
163
165 // Polarization of the parent particle
166
168 // Pointer to particle table
169
170 static const G4String noName;
171
179
181 // Number of daughters
182
184 // Control flag for output message
185 // 0: Silent
186 // 1: Warning message
187 // 2: More
188};
189
190// ------------------------------------------------------------
191// Inline methods
192// ------------------------------------------------------------
193
194inline
196{
197 return (this->rbranch < right.rbranch);
198}
199
200inline
202{
203 // pointers to daughter particles are filled, if they are not set yet
205
206 // get the pointer to a daughter particle
207 if ( (anIndex>=0) && (anIndex<numberOfDaughters) )
208 {
209 return G4MT_daughters[anIndex];
210 }
211 else
212 {
213 if (verboseLevel>0)
214 G4cout << "G4VDecayChannel::GetDaughter index out of range "
215 << anIndex << G4endl;
216 return nullptr;
217 }
218}
219
220inline
222{
223 if ( (anIndex>=0) && (anIndex<numberOfDaughters) )
224 {
225 return *daughters_name[anIndex];
226 }
227 else
228 {
229 if (verboseLevel>0)
230 {
231 G4cout << "G4VDecayChannel::GetDaughterName ";
232 G4cout << "index out of range " << anIndex << G4endl;
233 }
234 return GetNoName();
235 }
236}
237
238inline
240{
241 if ( (anIndex>=0) && (anIndex<numberOfDaughters) )
242 {
243 return G4MT_daughters_mass[anIndex];
244 }
245 else
246 {
247 if (verboseLevel>0)
248 {
249 G4cout << "G4VDecayChannel::GetDaughterMass ";
250 G4cout << "index out of range " << anIndex << G4endl;
251 }
252 return 0.0;
253 }
254}
255
256inline
258{
259 // the pointer to the parent particle is filled, if it is not set yet
261 // get the pointer to the parent particle
262 return G4MT_parent;
263}
264
265inline
267{
268 return *parent_name;
269}
270
271inline
273{
274 return G4MT_parent_mass;
275}
276
277inline
278void G4VDecayChannel::SetParent(const G4String& particle_name)
279{
280 if (parent_name != nullptr) delete parent_name;
281 parent_name = new G4String(particle_name);
282 G4MT_parent = nullptr;
283}
284
285inline
287{
288 return numberOfDaughters;
289}
290
291inline
293{
294 return kinematics_name;
295}
296
297inline
299{
300 return rbranch;
301}
302
303inline
305{
306 verboseLevel = value;
307}
308
309inline
311{
312 return verboseLevel;
313}
314
315inline
317{
318 return rangeMass;
319}
320
321inline
323
324inline
326{
327 parent_polarization = polar;
328}
329
330inline
332{
333 return parent_polarization;
334}
335
336inline
338{
340 if ( G4MT_daughters == nullptr )
341 {
342 l.unlock();
343 FillDaughters();
344 }
345}
346
347inline
349{
351 if ( G4MT_parent == nullptr )
352 {
353 l.unlock();
354 FillParent();
355 }
356}
357
358#endif
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition ** G4MT_daughters
virtual ~G4VDecayChannel()
G4double GetDaughterMass(G4int anIndex) const
G4String * parent_name
void SetPolarization(const G4ThreeVector &)
G4double GetBR() const
const G4String & GetParentName() const
G4bool operator<(const G4VDecayChannel &right) const
void SetBR(G4double value)
const G4ThreeVector & GetPolarization() const
void SetVerboseLevel(G4int value)
G4String ** daughters_name
G4double DynamicalMass(G4double massPDG, G4double width, G4double maxDev=1.0) const
G4int GetVerboseLevel() const
G4double G4MT_parent_mass
G4double GetRangeMass() const
void SetNumberOfDaughters(G4int value)
const G4String & GetKinematicsName() const
G4int GetNumberOfDaughters() const
static const G4String noName
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
G4ParticleDefinition * G4MT_parent
G4bool operator!=(const G4VDecayChannel &r) const
virtual G4bool IsOKWithParentMass(G4double parentMass)
void CheckAndFillDaughters()
G4ParticleDefinition * GetDaughter(G4int anIndex)
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
G4double * G4MT_daughters_mass
G4VDecayChannel & operator=(const G4VDecayChannel &)
G4int GetAngularMomentum()
G4double GetParentMass() const
G4ParticleDefinition * GetParent()
void SetRangeMass(G4double val)
const G4String & GetDaughterName(G4int anIndex) const
G4ParticleTable * particletable
G4bool operator==(const G4VDecayChannel &r) const
G4ThreeVector parent_polarization
G4String kinematics_name
G4double * G4MT_daughters_width
void SetParent(const G4ParticleDefinition *particle_type)