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