Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 //
27 // $Id$
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 "G4ios.hh"
42 #include "globals.hh"
43 #include <cmath>
44 
46 class G4DecayProducts;
47 class G4ParticleTable;
48 
50 {
51  // Class Description
52  // This class is a abstract class to describe decay kinematics
53  //
54 
55  public:
56  //Constructors
57  G4VDecayChannel(const G4String &aName, G4int Verbose = 1);
58  G4VDecayChannel(const G4String &aName,
59  const G4String& theParentName,
60  G4double theBR,
61  G4int theNumberOfDaughters,
62  const G4String& theDaughterName1,
63  const G4String& theDaughterName2 = "",
64  const G4String& theDaughterName3 = "",
65  const G4String& theDaughterName4 = "" );
66 
67  // Destructor
68  virtual ~G4VDecayChannel();
69 
70  protected:
71  // default constructor
73  // copy constructor and assignment operatotr
76 
77  public:
78  // equality operators
79  G4int operator==(const G4VDecayChannel &right) const {return (this == &right);};
80  G4int operator!=(const G4VDecayChannel &right) const {return (this != &right);};
81 
82  // less-than operator is defined for G4DecayTable
83  G4int operator<(const G4VDecayChannel &right) const;
84 
85  public: // With Description
86  virtual G4DecayProducts* DecayIt(G4double parentMass = -1.0) = 0;
87 
88  public: // With Description
89  //get kinematics name
90  const G4String& GetKinematicsName() const;
91  //get branching ratio
92  G4double GetBR() const;
93  //get number of daughter particles
94  G4int GetNumberOfDaughters() const;
95 
96  //get the pointer to the parent particle
98  //get the pointer to a daughter particle
100 
101  //get the angular momentum of the decay
103  //get the name of the parent particle
104  const G4String& GetParentName() const;
105  //get the name of a daughter particle
106  const G4String& GetDaughterName(G4int anIndex) const;
107 
108  // get mass of parent
109  G4double GetParentMass() const;
110  G4double GetDaughterMass(G4int anIndex) const;
111 
112  //set the parent particle (by name or by pointer)
113  void SetParent(const G4ParticleDefinition * particle_type);
114  void SetParent(const G4String &particle_name);
115  //set branching ratio
116  void SetBR(G4double value);
117  //set number of daughter particles
119  //set a daughter particle (by name or by pointer)
120  void SetDaughter(G4int anIndex,
121  const G4ParticleDefinition * particle_type);
122  void SetDaughter(G4int anIndex,
123  const G4String &particle_name);
124 
125  protected:
126  // kinematics name
128  // branching ratio [0.0 - 1.0]
130  // number of daughters
132  // parent particle
134  //daughter particles
136 
137  protected: // With Description
138  // celar daughters array
139  void ClearDaughtersName();
140 
141  protected:
142  // pointer to particle table
144 
145  // temporary buffers of pointers to G4ParticleDefinition
148 
149  // parent mass
152 
153 
154  // fill daughters array
155  void FillDaughters();
156  // fill parent
157  void FillParent();
158 
159  public: // With Description
161  G4int GetVerboseLevel() const;
162  void DumpInfo();
163 
164  private:
165  const G4String& GetNoName() const;
166 
167  protected:
168  // controle flag for output message
170  // 0: Silent
171  // 1: Warning message
172  // 2: More
173 
174  static const G4String noName;
175 };
176 
177 inline
179 {
180  return (this->rbranch < right.rbranch);
181 }
182 
183 inline
185  {
186  //pointers to daughter particles are filled, if they are not set yet
187  if (daughters == 0) FillDaughters();
188 
189  //get the pointer to a daughter particle
190  if ( (anIndex>=0) && (anIndex<numberOfDaughters) ) {
191  return daughters[anIndex];
192  } else {
193  if (verboseLevel>0)
194  G4cout << "G4VDecayChannel::GetDaughter index out of range "<<anIndex<<G4endl;
195  return 0;
196  }
197 }
198 
199 inline
201 {
202  if ( (anIndex>=0) && (anIndex<numberOfDaughters) ) {
203  return *daughters_name[anIndex];
204  } else {
205  if (verboseLevel>0){
206  G4cout << "G4VDecayChannel::GetDaughterName ";
207  G4cout << "index out of range " << anIndex << G4endl;
208  }
209  return GetNoName();
210  }
211 }
212 
213 inline
215 {
216  if ( (anIndex>=0) && (anIndex<numberOfDaughters) ) {
217  return daughters_mass[anIndex];
218  } else {
219  if (verboseLevel>0){
220  G4cout << "G4VDecayChannel::GetDaughterMass ";
221  G4cout << "index out of range " << anIndex << G4endl;
222  }
223  return 0.0;
224  }
225 }
226 
227 inline
229 {
230  //the pointer to the parent particle is filled, if it is not set yet
231  if (parent == 0) FillParent();
232  //get the pointer to the parent particle
233  return parent;
234 }
235 
236 inline
238 {
239  return *parent_name;
240 }
241 
242 inline
244 {
245  return parent_mass;
246 }
247 
248 
249 inline
250  void G4VDecayChannel::SetParent(const G4String &particle_name)
251 {
252  if (parent_name != 0) delete parent_name;
253  parent_name = new G4String(particle_name);
254  parent = 0;
255 }
256 
257 inline
259 {
260  return numberOfDaughters;
261 }
262 
263 inline
265 
266 inline
268 
269 inline
271 
272 inline
274 
275 inline
277 
278 
279 
280 #endif
281 
282 
283 
284 
285 
286