Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4VDecayChannel.hh
이 파일의 문서화 페이지로 가기
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 // $Id: G4VDecayChannel.hh 105720 2017-08-16 12:38:10Z gcosmo $
28 //
29 //
30 // ------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 // History: first implementation, based on object model of
34 // 27 July 1996 H.Kurashige
35 // 30 May 1997 H.Kurashige
36 // 23 Mar. 2000 H.Weber : add GetAngularMomentum()
37 // ------------------------------------------------------------
38 #ifndef G4VDecayChannel_h
39 #define G4VDecayChannel_h 1
40 
41 #include <cmath>
42 
43 #include "G4ios.hh"
44 #include "globals.hh"
45 #include "G4ThreeVector.hh"
46 #include "G4Threading.hh"
47 #include "G4AutoLock.hh"
48 
50 class G4DecayProducts;
51 class G4ParticleTable;
52 
54 {
55  // Class Description:
56  //
57  // Abstract class to describe decay kinematics
58 
59  public: // with description
60 
61  // Constructors
62  G4VDecayChannel(const G4String &aName, G4int Verbose = 1);
63  G4VDecayChannel(const G4String &aName,
64  const G4String& theParentName,
65  G4double theBR,
66  G4int theNumberOfDaughters,
67  const G4String& theDaughterName1,
68  const G4String& theDaughterName2 = "",
69  const G4String& theDaughterName3 = "",
70  const G4String& theDaughterName4 = "" );
71 
72  // Destructor
73  virtual ~G4VDecayChannel();
74 
75  // equality operators
76  G4int operator==(const G4VDecayChannel &right) const {return (this == &right);}
77  G4int operator!=(const G4VDecayChannel &right) const {return (this != &right);}
78 
79  // less-than operator is defined for G4DecayTable
80  G4int operator<(const G4VDecayChannel &right) const;
81 
82  virtual G4DecayProducts* DecayIt(G4double parentMass = -1.0) = 0;
83 
84  // get kinematics name
85  const G4String& GetKinematicsName() const;
86  // get branching ratio
87  G4double GetBR() const;
88  // get number of daughter particles
89  G4int GetNumberOfDaughters() const;
90 
91  // get the pointer to the parent particle
93  // get the pointer to a daughter particle
95 
96  // get the angular momentum of the decay
98  // get the name of the parent particle
99  const G4String& GetParentName() const;
100  //get the name of a daughter particle
101  const G4String& GetDaughterName(G4int anIndex) const;
102 
103  // get mass of parent
104  G4double GetParentMass() const;
105  G4double GetDaughterMass(G4int anIndex) const;
106 
107  // set the parent particle (by name or by pointer)
108  void SetParent(const G4ParticleDefinition * particle_type);
109  void SetParent(const G4String &particle_name);
110  // set branching ratio
111  void SetBR(G4double value);
112  // set number of daughter particles
114  // set a daughter particle (by name or by pointer)
115  void SetDaughter(G4int anIndex,
116  const G4ParticleDefinition * particle_type);
117  void SetDaughter(G4int anIndex,
118  const G4String &particle_name);
119 
121  G4int GetVerboseLevel() const;
122  void DumpInfo();
123 
124  G4double GetRangeMass() const;
125  void SetRangeMass(G4double val);
126  virtual G4bool IsOKWithParentMass(G4double parentMass);
127 
128  void SetPolarization(const G4ThreeVector&);
129  const G4ThreeVector& GetPolarization() const;
130 
131  protected: // with description
132 
133  // clear daughters array
134  void ClearDaughtersName();
135 
136  protected:
137  inline void CheckAndFillDaughters();
138  inline void CheckAndFillParent();
139 
140  private:
141  // fill daughters array
142  void FillDaughters();
143  // fill parent
144  void FillParent();
145 
146  protected: // without description
147 
148  // default constructor
149  G4VDecayChannel();
150 
151  // copy constructor and assignment operator
154 
155  private:
156 
157  const G4String& GetNoName() const;
158 
159  protected:
160  G4double DynamicalMass(G4double massPDG, G4double width, G4double maxDev = +1.) const;
161 
162  protected:
163 
164  // kinematics name
166  // branching ratio [0.0 - 1.0]
168  // number of daughters
170  // parent particle
172  // daughter particles
174 
175  // range of mass allowed in decay
177 
178  // polarization of the parent particle
180 
181  // pointer to particle table
183 
184  // control flag for output message
186  // 0: Silent
187  // 1: Warning message
188  // 2: More
189 
190  static const G4String noName;
191 
199 };
200 
201 // ------------------------------------------------------------
202 // Inline methods
203 // ------------------------------------------------------------
204 
205 inline
207 {
208  return (this->rbranch < right.rbranch);
209 }
210 
211 inline
213  {
214  // pointers to daughter particles are filled, if they are not set yet
216 
217  // get the pointer to a daughter particle
218  if ( (anIndex>=0) && (anIndex<numberOfDaughters) ) {
219  return G4MT_daughters[anIndex];
220  } else {
221  if (verboseLevel>0)
222  G4cout << "G4VDecayChannel::GetDaughter index out of range "<<anIndex<<G4endl;
223  return 0;
224  }
225 }
226 
227 inline
229 {
230  if ( (anIndex>=0) && (anIndex<numberOfDaughters) ) {
231  return *daughters_name[anIndex];
232  } else {
233  if (verboseLevel>0){
234  G4cout << "G4VDecayChannel::GetDaughterName ";
235  G4cout << "index out of range " << anIndex << G4endl;
236  }
237  return GetNoName();
238  }
239 }
240 
241 inline
243 {
244  if ( (anIndex>=0) && (anIndex<numberOfDaughters) ) {
245  return G4MT_daughters_mass[anIndex];
246  } else {
247  if (verboseLevel>0){
248  G4cout << "G4VDecayChannel::GetDaughterMass ";
249  G4cout << "index out of range " << anIndex << G4endl;
250  }
251  return 0.0;
252  }
253 }
254 
255 inline
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 
264 inline
266 {
267  return *parent_name;
268 }
269 
270 inline
272 {
273  return G4MT_parent_mass;
274 }
275 
276 inline
277  void G4VDecayChannel::SetParent(const G4String &particle_name)
278 {
279  if (parent_name != 0) delete parent_name;
280  parent_name = new G4String(particle_name);
281  G4MT_parent = 0;
282 }
283 
284 inline
286 {
287  return numberOfDaughters;
288 }
289 
290 inline
292 
293 inline
295 
296 inline
298 
299 inline
301 
302 inline
304 
305 inline
306  void G4VDecayChannel::SetRangeMass(G4double val){ if(val>=0.) rangeMass=val; }
307 
308 inline
310 {
311  parent_polarization = polar;
312 }
313 
314 inline
316 {
317  return parent_polarization;
318 }
319 
320 inline
323  if ( G4MT_daughters ==0 ) {
324  l.unlock();
325  FillDaughters();
326  }
327 }
328 
331  if ( G4MT_parent == 0 ) {
332  l.unlock();
333  FillParent();
334  }
335 }
336 #endif
G4String ** daughters_name
G4int GetAngularMomentum()
G4ParticleTable * particletable
G4double GetBR() const
void SetRangeMass(G4double val)
void SetVerboseLevel(G4int value)
G4int GetNumberOfDaughters() const
#define G4endl
Definition: G4ios.hh:61
G4String * parent_name
G4double GetRangeMass() const
const G4ThreeVector & GetPolarization() const
virtual ~G4VDecayChannel()
const G4String & GetParentName() const
void SetBR(G4double value)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
#define width
const XML_Char int const XML_Char * value
Definition: expat.h:331
G4ParticleDefinition * G4MT_parent
static const G4String noName
G4ParticleDefinition ** G4MT_daughters
G4double GetParentMass() const
void SetNumberOfDaughters(G4int value)
G4VDecayChannel & operator=(const G4VDecayChannel &)
G4ThreeVector parent_polarization
G4String kinematics_name
G4ParticleDefinition * GetDaughter(G4int anIndex)
void SetParent(const G4ParticleDefinition *particle_type)
G4int GetVerboseLevel() const
const G4String & GetDaughterName(G4int anIndex) const
G4double DynamicalMass(G4double massPDG, G4double width, G4double maxDev=+1.) const
void SetPolarization(const G4ThreeVector &)
G4double * G4MT_daughters_width
int G4int
Definition: G4Types.hh:78
G4int operator!=(const G4VDecayChannel &right) const
G4int operator<(const G4VDecayChannel &right) const
const G4String & GetKinematicsName() const
const G4String & GetNoName() const
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition * GetParent()
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
virtual G4bool IsOKWithParentMass(G4double parentMass)
G4int operator==(const G4VDecayChannel &right) const
G4double GetDaughterMass(G4int anIndex) const
G4double G4MT_parent_mass
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
void CheckAndFillDaughters()
G4double * G4MT_daughters_mass
std::mutex G4Mutex
Definition: G4Threading.hh:84