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 // 27 // 28 // 29 // ------------------------------------------------------------ 30 // GEANT 4 class header file 31 // 32 // History: first implementation, based on object model of 33 // 7 July 1996 H.Kurashige 34 // ------------------------------------------------------------ 35 // New Physics scheme 18 Jan. 1997 H.Kurahige 36 // ------------------------------------------------------------ 37 // modified 4 Feb. 1997 H.Kurahige 38 // modified 8 Sep. 1997 H.Kurahige 39 // remove BuildPhysicsTable() 27 Nov. 1997 H.Kurashige 40 // modified for new ParticleChange 12 Mar. 1998 H.Kurashige 41 // added aPhysicsTable 2 Aug. 1998 H.Kurashige 42 // PreAssignedDecayTime 18 Jan. 2001 H.Kurashige 43 // Add External Decayer 23 Feb. 2001 H.Kurashige 44 // Remove PhysicsTable 12 Feb. 2002 H.Kurashige 45 // Fixed bug in PostStepGPIL 46 // in case of stopping during AlongStepDoIt 12 Mar. 2004 H.Kurashige 47 // Add GetRemainderLifeTime 10 Aug/2004 H.Kurashige 48 // Add DaughterPolarization 23 July 2008 H.Kurashige 49 50 51 #ifndef G4Decay_h 52 #define G4Decay_h 1 53 54 #include "G4ios.hh" 55 #include "globals.hh" 56 #include "G4VRestDiscreteProcess.hh" 57 #include "G4ParticleChangeForDecay.hh" 58 #include "G4DecayProcessType.hh" 59 60 class G4VExtDecayer; 61 62 class G4Decay : public G4VRestDiscreteProcess 63 { 64 // Class Description 65 // This class is a decay process 66 67 public: 68 // Constructors 69 G4Decay(const G4String& processName ="Decay"); 70 71 // Destructor 72 virtual ~G4Decay(); 73 74 private: 75 // copy constructor 76 G4Decay(const G4Decay &right); 77 78 // Assignment Operation (generated) 79 G4Decay & operator=(const G4Decay &right); 80 81 public: //With Description 82 // G4Decay Process has both 83 // PostStepDoIt (for decay in flight) 84 // and 85 // AtRestDoIt (for decay at rest) 86 87 virtual G4VParticleChange *PostStepDoIt( 88 const G4Track& aTrack, 89 const G4Step& aStep 90 ) override; 91 92 virtual G4VParticleChange* AtRestDoIt( 93 const G4Track& aTrack, 94 const G4Step& aStep 95 ) override; 96 97 virtual void BuildPhysicsTable(const G4ParticleDefinition&) override; 98 // In G4Decay, thePhysicsTable stores values of 99 // beta * std::sqrt( 1 - beta*beta) 100 // as a function of normalized kinetic enregy (=Ekin/mass), 101 // becasuse this table is universal for all particle types, 102 103 104 virtual G4bool IsApplicable(const G4ParticleDefinition&) override; 105 // returns "true" if the decay process can be applied to 106 // the particle type. 107 108 protected: // With Description 109 virtual G4VParticleChange* DecayIt( 110 const G4Track& aTrack, 111 const G4Step& aStep 112 ); 113 // The DecayIt() method returns by pointer a particle-change object, 114 // which has information of daughter particles. 115 116 // Set daughter polarization 117 // NO OPERATION in the base class of G4Decay 118 virtual void DaughterPolarization(const G4Track& aTrack, 119 G4DecayProducts* products); 120 121 public: 122 virtual G4double AtRestGetPhysicalInteractionLength( 123 const G4Track& track, 124 G4ForceCondition* condition 125 ) override; 126 127 virtual G4double PostStepGetPhysicalInteractionLength( 128 const G4Track& track, 129 G4double previousStepSize, 130 G4ForceCondition* condition 131 ) override; 132 133 protected: // With Description 134 // GetMeanFreePath returns ctau*beta*gamma for decay in flight 135 // GetMeanLifeTime returns ctau for decay at rest 136 virtual G4double GetMeanFreePath(const G4Track& aTrack, 137 G4double previousStepSize, 138 G4ForceCondition* condition 139 ) override; 140 141 virtual G4double GetMeanLifeTime(const G4Track& aTrack, 142 G4ForceCondition* condition 143 ) override; 144 145 public: //With Description 146 virtual void StartTracking(G4Track*) override; 147 virtual void EndTracking() override; 148 // inform Start/End of tracking for each track to the physics process 149 150 public: //With Description 151 void SetExtDecayer(G4VExtDecayer*); 152 const G4VExtDecayer* GetExtDecayer() const; 153 // Set/Get External Decayer 154 155 G4double GetRemainderLifeTime() const; 156 //Get Remainder of life time at rest decay 157 158 virtual void ProcessDescription(std::ostream& outFile) const override; 159 // 160 161 protected: 162 G4int verboseLevel; 163 // controle flag for output message 164 // 0: Silent 165 // 1: Warning message 166 // 2: More 167 168 protected: 169 // HighestValue. 170 const G4double HighestValue; 171 172 // Remainder of life time at rest 173 G4double fRemainderLifeTime; 174 175 // ParticleChange for decay process 176 G4ParticleChangeForDecay fParticleChangeForDecay; 177 178 // External Decayer 179 G4VExtDecayer* pExtDecayer; 180 }; 181 182 inline 183 G4VParticleChange* G4Decay::AtRestDoIt( 184 const G4Track& aTrack, 185 const G4Step& aStep 186 ) 187 { 188 return DecayIt(aTrack, aStep); 189 } 190 191 192 inline 193 const G4VExtDecayer* G4Decay::GetExtDecayer() const 194 { 195 return pExtDecayer; 196 } 197 198 inline 199 G4double G4Decay::GetRemainderLifeTime() const 200 { 201 return fRemainderLifeTime; 202 } 203 204 #endif 205 206 207 208 209 210 211 212 213 214 215