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 // G4VParticleChange inline methods implementation 27 // 28 // Author: Hisaya Kurashige, 23 March 1998 29 // -------------------------------------------------------------------- 30 31 inline G4Track* G4VParticleChange::GetSecondary(G4int anIndex) const 32 { 33 return theListOfSecondaries[anIndex]; 34 } 35 36 inline G4int G4VParticleChange::GetNumberOfSecondaries() const 37 { 38 return theNumberOfSecondaries; 39 } 40 41 inline void G4VParticleChange::ProposeTrackStatus(G4TrackStatus aStatus) 42 { 43 theStatusChange = aStatus; 44 } 45 46 inline const G4Track* G4VParticleChange::GetCurrentTrack() const 47 { 48 return theCurrentTrack; 49 } 50 51 inline G4TrackStatus G4VParticleChange::GetTrackStatus() const 52 { 53 return theStatusChange; 54 } 55 56 inline G4SteppingControl G4VParticleChange::GetSteppingControl() const 57 { 58 return theSteppingControlFlag; 59 } 60 61 inline void 62 G4VParticleChange::ProposeSteppingControl(G4SteppingControl StepControlFlag) 63 { 64 theSteppingControlFlag = StepControlFlag; 65 } 66 67 inline G4bool G4VParticleChange::GetFirstStepInVolume() const 68 { 69 return theFirstStepInVolume; 70 } 71 72 inline G4bool G4VParticleChange::GetLastStepInVolume() const 73 { 74 return theLastStepInVolume; 75 } 76 77 inline void G4VParticleChange::ProposeFirstStepInVolume(G4bool flag) 78 { 79 theFirstStepInVolume = flag; 80 } 81 82 inline void G4VParticleChange::ProposeLastStepInVolume(G4bool flag) 83 { 84 theLastStepInVolume = flag; 85 } 86 87 inline G4double G4VParticleChange::GetLocalEnergyDeposit() const 88 { 89 return theLocalEnergyDeposit; 90 } 91 92 inline void G4VParticleChange::ProposeLocalEnergyDeposit(G4double anEnergyPart) 93 { 94 theLocalEnergyDeposit = anEnergyPart; 95 } 96 97 inline G4double G4VParticleChange::GetNonIonizingEnergyDeposit() const 98 { 99 return theNonIonizingEnergyDeposit; 100 } 101 102 inline void 103 G4VParticleChange::ProposeNonIonizingEnergyDeposit(G4double anEnergyPart) 104 { 105 theNonIonizingEnergyDeposit = anEnergyPart; 106 } 107 108 inline G4double G4VParticleChange::GetTrueStepLength() const 109 { 110 return theTrueStepLength; 111 } 112 113 inline void G4VParticleChange::ProposeTrueStepLength(G4double aLength) 114 { 115 theTrueStepLength = aLength; 116 } 117 118 inline void G4VParticleChange::SetVerboseLevel(G4int vLevel) 119 { 120 verboseLevel = vLevel; 121 } 122 123 inline G4int G4VParticleChange::GetVerboseLevel() const 124 { 125 return verboseLevel; 126 } 127 128 inline G4double G4VParticleChange::GetParentWeight() const 129 { 130 return theParentWeight; 131 } 132 133 inline G4double G4VParticleChange::GetWeight() const 134 { 135 return theParentWeight; 136 } 137 138 // ---------------------------------------------------------------- 139 // inline functions for Initialization 140 // 141 inline void G4VParticleChange::InitializeLocalEnergyDeposit() 142 { 143 // clear theLocalEnergyDeposited 144 theLocalEnergyDeposit = 0.0; 145 theNonIonizingEnergyDeposit = 0.0; 146 } 147 148 inline void G4VParticleChange::InitializeSteppingControl() 149 { 150 // SteppingControlFlag 151 theSteppingControlFlag = NormalCondition; 152 } 153 154 inline void G4VParticleChange::Clear() 155 { 156 theNumberOfSecondaries = 0; 157 theFirstStepInVolume = false; 158 theLastStepInVolume = false; 159 } 160 161 inline void G4VParticleChange::InitializeStatusChange(const G4Track& track) 162 { 163 // set TrackStatus equal to the parent track's one 164 theStatusChange = track.GetTrackStatus(); 165 theCurrentTrack = &track; 166 } 167 168 inline void G4VParticleChange::InitializeParentWeight(const G4Track& track) 169 { 170 // set the parent track's weight 171 theParentWeight = track.GetWeight(); 172 isParentWeightProposed = false; 173 } 174 175 inline void G4VParticleChange::InitializeFromStep(const G4Step* step) 176 { 177 // set the parent track's global time at the pre-step point 178 theParentGlobalTime = step->GetPreStepPoint()->GetGlobalTime(); 179 180 theTrueStepLength = step->GetStepLength(); 181 theFirstStepInVolume = step->IsFirstStepInVolume(); 182 theLastStepInVolume = step->IsLastStepInVolume(); 183 } 184 185 // ---------------------------------------------------------------- 186 // methods for handling secondaries 187 // 188 inline void G4VParticleChange::InitializeSecondaries() 189 { 190 theNumberOfSecondaries = 0; 191 } 192 193 inline void G4VParticleChange::SetNumberOfSecondaries(G4int totSecondaries) 194 { 195 if(totSecondaries > theSizeOftheListOfSecondaries) 196 { 197 theListOfSecondaries.resize(totSecondaries, nullptr); 198 theSizeOftheListOfSecondaries = totSecondaries; 199 } 200 theNumberOfSecondaries = 0; 201 } 202 203 // ---------------------------------------------------------------- 204 // total 205 // 206 inline void G4VParticleChange::Initialize(const G4Track& track) 207 { 208 InitializeStatusChange(track); 209 InitializeLocalEnergyDeposit(); 210 InitializeSteppingControl(); 211 InitializeSecondaries(); 212 InitializeParentWeight(track); 213 InitializeFromStep(track.GetStep()); 214 } 215 216 inline void G4VParticleChange::ClearDebugFlag() 217 { 218 debugFlag = false; 219 } 220 221 inline void G4VParticleChange::SetDebugFlag() 222 { 223 debugFlag = true; 224 } 225 226 inline G4bool G4VParticleChange::GetDebugFlag() const 227 { 228 return debugFlag; 229 } 230 231 inline void G4VParticleChange::SetSecondaryWeightByProcess(G4bool flag) 232 { 233 fSetSecondaryWeightByProcess = flag; 234 } 235 236 inline G4bool G4VParticleChange::IsSecondaryWeightSetByProcess() const 237 { 238 return fSetSecondaryWeightByProcess; 239 } 240 241 inline void G4VParticleChange::ProposeWeight(G4double w) 242 { 243 theParentWeight = w; 244 isParentWeightProposed = true; 245 } 246 247 inline void G4VParticleChange::ProposeParentWeight(G4double w) 248 { 249 ProposeWeight(w); 250 } 251 252 inline G4double G4VParticleChange::ComputeBeta(G4double ekin) 253 { 254 G4double mass = theCurrentTrack->GetDefinition()->GetPDGMass(); 255 return std::sqrt(ekin*(ekin + 2*mass))/(ekin + mass); 256 } 257