Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 // Constructors
63
64 virtual ~G4VDecayChannel();
65 // Destructor
66
67 G4bool operator==(const G4VDecayChannel& r) const { return (this == &r); }
68 G4bool operator!=(const G4VDecayChannel& r) const { return (this != &r); }
69 // Equality operators
70
71 inline G4bool operator<(const G4VDecayChannel& right) const;
72 // Less-than operator is defined for G4DecayTable
73
74 virtual G4DecayProducts* DecayIt(G4double parentMass = -1.0) = 0;
75
76 inline const G4String& GetKinematicsName() const;
77 // Get kinematics name
78 inline G4double GetBR() const;
79 // Get branching ratio
80 inline G4int GetNumberOfDaughters() const;
81 // Get number of daughter particles
82
84 // Get the pointer to the parent particle
85 inline G4ParticleDefinition* GetDaughter(G4int anIndex);
86 // Get the pointer to a daughter particle
87
89 // Get the angular momentum of the decay
90 inline const G4String& GetParentName() const;
91 // Get the name of the parent particle
92 inline const G4String& GetDaughterName(G4int anIndex) const;
93 // Get the name of a daughter particle
94
95 inline G4double GetParentMass() const;
96 // Get mass of parent
97 inline G4double GetDaughterMass(G4int anIndex) const;
98 // Get mass of daughter
99
100 void SetParent(const G4ParticleDefinition * particle_type);
101 inline void SetParent(const G4String &particle_name);
102 // Set the parent particle (by name or by pointer)
103 void SetBR(G4double value);
104 // Set branching ratio
105 void SetNumberOfDaughters(G4int value);
106 // Set number of daughter particles
107 void SetDaughter(G4int anIndex, const G4ParticleDefinition* particle_type);
108 void SetDaughter(G4int anIndex, const G4String& particle_name);
109 // Set a daughter particle (by name or by pointer)
110
111 inline void SetVerboseLevel(G4int value);
112 inline G4int GetVerboseLevel() const;
113 void DumpInfo();
114
115 inline G4double GetRangeMass() const;
116 inline void SetRangeMass(G4double val);
117 virtual G4bool IsOKWithParentMass(G4double parentMass);
118
119 void SetPolarization(const G4ThreeVector&);
120 inline const G4ThreeVector& GetPolarization() const;
121
122 protected:
123
124 void ClearDaughtersName();
125 // Clear daughters array
126
127 inline void CheckAndFillDaughters();
128 inline void CheckAndFillParent();
129
131 G4double maxDev = 1.0) const;
132
134 // Default constructor
135
138 // Copy constructor and assignment operator
139
140 private:
141
142 void FillDaughters();
143 // Fill daughters array
144 void FillParent();
145 // Fill parent
146
147 const G4String& GetNoName() const;
148
149 protected:
150
152 // Kinematics name
154 // Branching ratio [0.0 - 1.0]
156 // Parent particle
158 // Daughter particles
159
161 // Range of mass allowed in decay
162
164 // Polarization of the parent particle
165
167 // Pointer to particle table
168
169 static const G4String noName;
170
178
180 // Number of daughters
181
183 // Control flag for output message
184 // 0: Silent
185 // 1: Warning message
186 // 2: More
187};
188
189// ------------------------------------------------------------
190// Inline methods
191// ------------------------------------------------------------
192
193inline
195{
196 return (this->rbranch < right.rbranch);
197}
198
199inline
201{
202 // pointers to daughter particles are filled, if they are not set yet
204
205 // get the pointer to a daughter particle
206 if ( (anIndex>=0) && (anIndex<numberOfDaughters) )
207 {
208 return G4MT_daughters[anIndex];
209 }
210 else
211 {
212 if (verboseLevel>0)
213 G4cout << "G4VDecayChannel::GetDaughter index out of range "
214 << anIndex << G4endl;
215 return nullptr;
216 }
217}
218
219inline
221{
222 if ( (anIndex>=0) && (anIndex<numberOfDaughters) )
223 {
224 return *daughters_name[anIndex];
225 }
226 else
227 {
228 if (verboseLevel>0)
229 {
230 G4cout << "G4VDecayChannel::GetDaughterName ";
231 G4cout << "index out of range " << anIndex << G4endl;
232 }
233 return GetNoName();
234 }
235}
236
237inline
239{
240 if ( (anIndex>=0) && (anIndex<numberOfDaughters) )
241 {
242 return G4MT_daughters_mass[anIndex];
243 }
244 else
245 {
246 if (verboseLevel>0)
247 {
248 G4cout << "G4VDecayChannel::GetDaughterMass ";
249 G4cout << "index out of range " << anIndex << G4endl;
250 }
251 return 0.0;
252 }
253}
254
255inline
257{
258 // the pointer to the parent particle is filled, if it is not set yet
260 // get the pointer to the parent particle
261 return G4MT_parent;
262}
263
264inline
266{
267 return *parent_name;
268}
269
270inline
272{
273 return G4MT_parent_mass;
274}
275
276inline
277void G4VDecayChannel::SetParent(const G4String& particle_name)
278{
279 if (parent_name != nullptr) delete parent_name;
280 parent_name = new G4String(particle_name);
281 G4MT_parent = nullptr;
282}
283
284inline
286{
287 return numberOfDaughters;
288}
289
290inline
292{
293 return kinematics_name;
294}
295
296inline
298{
299 return rbranch;
300}
301
302inline
304{
305 verboseLevel = value;
306}
307
308inline
310{
311 return verboseLevel;
312}
313
314inline
316{
317 return rangeMass;
318}
319
320inline
322
323inline
325{
326 parent_polarization = polar;
327}
328
329inline
331{
332 return parent_polarization;
333}
334
335inline
337{
339 if ( G4MT_daughters == nullptr )
340 {
341 l.unlock();
342 FillDaughters();
343 }
344}
345
346inline
348{
350 if ( G4MT_parent == nullptr )
351 {
352 l.unlock();
353 FillParent();
354 }
355}
356
357#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)