[819] | 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 | // $Id: G4VEnergyLossProcess.hh,v 1.76 2007/11/07 18:38:49 vnivanch Exp $ |
---|
| 27 | // GEANT4 tag $Name: |
---|
| 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // GEANT4 Class header file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4VEnergyLossProcess |
---|
| 35 | // |
---|
| 36 | // Author: Vladimir Ivanchenko on base of Laszlo Urban code |
---|
| 37 | // |
---|
| 38 | // Creation date: 03.01.2002 |
---|
| 39 | // |
---|
| 40 | // Modifications: |
---|
| 41 | // |
---|
| 42 | // 26-12-02 Secondary production moved to derived classes (V.Ivanchenko) |
---|
| 43 | // 20-01-03 Migrade to cut per region (V.Ivanchenko) |
---|
| 44 | // 24-01-03 Make models region aware (V.Ivanchenko) |
---|
| 45 | // 05-02-03 Fix compilation warnings (V.Ivanchenko) |
---|
| 46 | // 13-02-03 SubCutoffProcessors defined for regions (V.Ivanchenko) |
---|
| 47 | // 17-02-03 Fix problem of store/restore tables (V.Ivanchenko) |
---|
| 48 | // 26-02-03 Region dependent step limit (V.Ivanchenko) |
---|
| 49 | // 26-03-03 Add GetDEDXDispersion (V.Ivanchenko) |
---|
| 50 | // 09-04-03 Fix problem of negative range limit for non integral (V.Ivanchenko) |
---|
| 51 | // 13-05-03 Add calculation of precise range (V.Ivanchenko) |
---|
| 52 | // 21-07-03 Add UpdateEmModel method (V.Ivanchenko) |
---|
| 53 | // 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko) |
---|
| 54 | // 14-01-04 Activate precise range calculation (V.Ivanchenko) |
---|
| 55 | // 10-03-04 Fix problem of step limit calculation (V.Ivanchenko) |
---|
| 56 | // 30-06-04 make destructor virtual (V.Ivanchenko) |
---|
| 57 | // 05-07-04 fix problem of GenericIons seen at small cuts (V.Ivanchenko) |
---|
| 58 | // 03-08-04 Add DEDX table to all processes for control on integral range(VI) |
---|
| 59 | // 06-08-04 Clear up names of member functions (V.Ivanchenko) |
---|
| 60 | // 27-08-04 Add NeedBuildTables method (V.Ivanchneko) |
---|
| 61 | // 09-09-04 Bug fix for the integral mode with 2 peaks (V.Ivanchneko) |
---|
| 62 | // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko) |
---|
| 63 | // 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko) |
---|
| 64 | // 11-04-05 Use MaxSecondaryEnergy from a model (V.Ivanchenko) |
---|
| 65 | // 10-01-05 Remove SetStepLimits (V.Ivanchenko) |
---|
| 66 | // 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko) |
---|
| 67 | // 13-01-06 Remove AddSubCutSecondaries and cleanup (V.Ivantchenko) |
---|
| 68 | // 20-01-06 Introduce G4EmTableType and reducing number of methods (VI) |
---|
| 69 | // 26-01-06 Add public method GetCSDARange (V.Ivanchenko) |
---|
| 70 | // 22-03-06 Add SetDynamicMassCharge (V.Ivanchenko) |
---|
| 71 | // 23-03-06 Use isIonisation flag (V.Ivanchenko) |
---|
| 72 | // 13-05-06 Add method to access model by index (V.Ivanchenko) |
---|
| 73 | // 14-01-07 add SetEmModel(index) and SetFluctModel() (mma) |
---|
| 74 | // 15-01-07 Add separate ionisation tables and reorganise get/set methods for |
---|
| 75 | // dedx tables (V.Ivanchenko) |
---|
| 76 | // 13-03-07 use SafetyHelper instead of navigator (V.Ivanchenko) |
---|
| 77 | // 27-07-07 use stl vector for emModels instead of C-array (V.Ivanchenko) |
---|
| 78 | // 25-09-07 More accurate handling zero xsect in |
---|
| 79 | // PostStepGetPhysicalInteractionLength (V.Ivanchenko) |
---|
| 80 | // 27-10-07 Virtual functions moved to source (V.Ivanchenko) |
---|
| 81 | // |
---|
| 82 | // Class Description: |
---|
| 83 | // |
---|
| 84 | // It is the unified energy loss process it calculates the continuous |
---|
| 85 | // energy loss for charged particles using a set of Energy Loss |
---|
| 86 | // models valid for different energy regions. There are a possibility |
---|
| 87 | // to create and access to dE/dx and range tables, or to calculate |
---|
| 88 | // that information on fly. |
---|
| 89 | |
---|
| 90 | // ------------------------------------------------------------------- |
---|
| 91 | // |
---|
| 92 | |
---|
| 93 | #ifndef G4VEnergyLossProcess_h |
---|
| 94 | #define G4VEnergyLossProcess_h 1 |
---|
| 95 | |
---|
| 96 | #include "G4VContinuousDiscreteProcess.hh" |
---|
| 97 | #include "globals.hh" |
---|
| 98 | #include "G4Material.hh" |
---|
| 99 | #include "G4MaterialCutsCouple.hh" |
---|
| 100 | #include "G4Track.hh" |
---|
| 101 | #include "G4EmModelManager.hh" |
---|
| 102 | #include "G4UnitsTable.hh" |
---|
| 103 | #include "G4ParticleChangeForLoss.hh" |
---|
| 104 | #include "G4EmTableType.hh" |
---|
| 105 | #include "G4PhysicsTable.hh" |
---|
| 106 | #include "G4PhysicsVector.hh" |
---|
| 107 | |
---|
| 108 | class G4Step; |
---|
| 109 | class G4ParticleDefinition; |
---|
| 110 | class G4VEmModel; |
---|
| 111 | class G4VEmFluctuationModel; |
---|
| 112 | class G4DataVector; |
---|
| 113 | class G4Region; |
---|
| 114 | class G4SafetyHelper; |
---|
| 115 | |
---|
| 116 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 117 | |
---|
| 118 | class G4VEnergyLossProcess : public G4VContinuousDiscreteProcess |
---|
| 119 | { |
---|
| 120 | public: |
---|
| 121 | |
---|
| 122 | G4VEnergyLossProcess(const G4String& name = "EnergyLoss", |
---|
| 123 | G4ProcessType type = fElectromagnetic); |
---|
| 124 | |
---|
| 125 | virtual ~G4VEnergyLossProcess(); |
---|
| 126 | |
---|
| 127 | //------------------------------------------------------------------------ |
---|
| 128 | // Virtual methods to be implemented in concrete processes |
---|
| 129 | //------------------------------------------------------------------------ |
---|
| 130 | |
---|
| 131 | virtual G4bool IsApplicable(const G4ParticleDefinition& p) = 0; |
---|
| 132 | |
---|
| 133 | virtual void PrintInfo() = 0; |
---|
| 134 | |
---|
| 135 | protected: |
---|
| 136 | |
---|
| 137 | virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition*, |
---|
| 138 | const G4ParticleDefinition*) = 0; |
---|
| 139 | |
---|
| 140 | //------------------------------------------------------------------------ |
---|
| 141 | // Methods with standard implementation; may be overwritten if needed |
---|
| 142 | //------------------------------------------------------------------------ |
---|
| 143 | protected: |
---|
| 144 | |
---|
| 145 | virtual G4double MinPrimaryEnergy(const G4ParticleDefinition*, |
---|
| 146 | const G4Material*, G4double cut); |
---|
| 147 | |
---|
| 148 | virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*, |
---|
| 149 | const G4DynamicParticle*, |
---|
| 150 | G4double& eloss, |
---|
| 151 | G4double& length); |
---|
| 152 | |
---|
| 153 | //------------------------------------------------------------------------ |
---|
| 154 | // Generic methods common to all ContinuousDiscrete processes |
---|
| 155 | //------------------------------------------------------------------------ |
---|
| 156 | public: |
---|
| 157 | |
---|
| 158 | void PrintInfoDefinition(); |
---|
| 159 | |
---|
| 160 | void PreparePhysicsTable(const G4ParticleDefinition&); |
---|
| 161 | |
---|
| 162 | void BuildPhysicsTable(const G4ParticleDefinition&); |
---|
| 163 | |
---|
| 164 | G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&); |
---|
| 165 | |
---|
| 166 | G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&); |
---|
| 167 | |
---|
| 168 | // Store PhysicsTable in a file. |
---|
| 169 | // Return false in case of failure at I/O |
---|
| 170 | G4bool StorePhysicsTable(const G4ParticleDefinition*, |
---|
| 171 | const G4String& directory, |
---|
| 172 | G4bool ascii = false); |
---|
| 173 | |
---|
| 174 | // Retrieve Physics from a file. |
---|
| 175 | // (return true if the Physics Table can be build by using file) |
---|
| 176 | // (return false if the process has no functionality or in case of failure) |
---|
| 177 | // File name should is constructed as processName+particleName and the |
---|
| 178 | // should be placed under the directory specifed by the argument. |
---|
| 179 | G4bool RetrievePhysicsTable(const G4ParticleDefinition*, |
---|
| 180 | const G4String& directory, |
---|
| 181 | G4bool ascii); |
---|
| 182 | |
---|
| 183 | protected: |
---|
| 184 | |
---|
| 185 | G4double GetMeanFreePath(const G4Track& track, |
---|
| 186 | G4double previousStepSize, |
---|
| 187 | G4ForceCondition* condition); |
---|
| 188 | |
---|
| 189 | G4double GetContinuousStepLimit(const G4Track& track, |
---|
| 190 | G4double previousStepSize, |
---|
| 191 | G4double currentMinimumStep, |
---|
| 192 | G4double& currentSafety); |
---|
| 193 | |
---|
| 194 | //------------------------------------------------------------------------ |
---|
| 195 | // Specific methods for along/post step EM processes |
---|
| 196 | //------------------------------------------------------------------------ |
---|
| 197 | |
---|
| 198 | public: |
---|
| 199 | |
---|
| 200 | void AddCollaborativeProcess(G4VEnergyLossProcess*); |
---|
| 201 | |
---|
| 202 | void SampleSubCutSecondaries(std::vector<G4Track*>&, const G4Step&, |
---|
| 203 | G4VEmModel* model, G4int matIdx, |
---|
| 204 | G4double& extraEdep); |
---|
| 205 | |
---|
| 206 | G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, |
---|
| 207 | const G4DynamicParticle* dp, |
---|
| 208 | G4double length); |
---|
| 209 | |
---|
| 210 | |
---|
| 211 | virtual G4double AlongStepGetPhysicalInteractionLength( |
---|
| 212 | const G4Track&, |
---|
| 213 | G4double previousStepSize, |
---|
| 214 | G4double currentMinimumStep, |
---|
| 215 | G4double& currentSafety, |
---|
| 216 | G4GPILSelection* selection |
---|
| 217 | ); |
---|
| 218 | |
---|
| 219 | virtual G4double PostStepGetPhysicalInteractionLength( |
---|
| 220 | const G4Track& track, |
---|
| 221 | G4double previousStepSize, |
---|
| 222 | G4ForceCondition* condition |
---|
| 223 | ); |
---|
| 224 | |
---|
| 225 | //------------------------------------------------------------------------ |
---|
| 226 | // Specific methods to build and access Physics Tables |
---|
| 227 | //------------------------------------------------------------------------ |
---|
| 228 | |
---|
| 229 | G4double MicroscopicCrossSection(G4double kineticEnergy, |
---|
| 230 | const G4MaterialCutsCouple* couple); |
---|
| 231 | |
---|
| 232 | G4PhysicsTable* BuildDEDXTable(G4EmTableType tType = fRestricted); |
---|
| 233 | |
---|
| 234 | G4PhysicsTable* BuildLambdaTable(G4EmTableType tType = fRestricted); |
---|
| 235 | |
---|
| 236 | void SetDEDXTable(G4PhysicsTable* p, G4EmTableType tType); |
---|
| 237 | void SetCSDARangeTable(G4PhysicsTable* pRange); |
---|
| 238 | void SetRangeTableForLoss(G4PhysicsTable* p); |
---|
| 239 | void SetInverseRangeTable(G4PhysicsTable* p); |
---|
| 240 | void SetSecondaryRangeTable(G4PhysicsTable* p); |
---|
| 241 | |
---|
| 242 | void SetLambdaTable(G4PhysicsTable* p); |
---|
| 243 | void SetSubLambdaTable(G4PhysicsTable* p); |
---|
| 244 | |
---|
| 245 | // Binning for dEdx, range, inverse range and labda tables |
---|
| 246 | inline void SetDEDXBinning(G4int nbins); |
---|
| 247 | inline void SetLambdaBinning(G4int nbins); |
---|
| 248 | |
---|
| 249 | // Binning for dEdx, range, and inverse range tables |
---|
| 250 | inline void SetDEDXBinningForCSDARange(G4int nbins); |
---|
| 251 | |
---|
| 252 | // Min kinetic energy for tables |
---|
| 253 | inline void SetMinKinEnergy(G4double e); |
---|
| 254 | inline G4double MinKinEnergy() const; |
---|
| 255 | |
---|
| 256 | // Max kinetic energy for tables |
---|
| 257 | inline void SetMaxKinEnergy(G4double e); |
---|
| 258 | inline G4double MaxKinEnergy() const; |
---|
| 259 | |
---|
| 260 | // Max kinetic energy for tables |
---|
| 261 | inline void SetMaxKinEnergyForCSDARange(G4double e); |
---|
| 262 | |
---|
| 263 | // Access to specific tables |
---|
| 264 | inline G4PhysicsTable* DEDXTable() const; |
---|
| 265 | inline G4PhysicsTable* DEDXTableForSubsec() const; |
---|
| 266 | inline G4PhysicsTable* DEDXunRestrictedTable() const; |
---|
| 267 | inline G4PhysicsTable* IonisationTable() const; |
---|
| 268 | inline G4PhysicsTable* IonisationTableForSubsec() const; |
---|
| 269 | inline G4PhysicsTable* CSDARangeTable() const; |
---|
| 270 | inline G4PhysicsTable* RangeTableForLoss() const; |
---|
| 271 | inline G4PhysicsTable* InverseRangeTable() const; |
---|
| 272 | inline G4PhysicsTable* LambdaTable(); |
---|
| 273 | inline G4PhysicsTable* SubLambdaTable(); |
---|
| 274 | |
---|
| 275 | // Return values for given G4MaterialCutsCouple |
---|
| 276 | inline G4double GetDEDX(G4double& kineticEnergy, const G4MaterialCutsCouple*); |
---|
| 277 | inline G4double GetDEDXForSubsec(G4double& kineticEnergy, |
---|
| 278 | const G4MaterialCutsCouple*); |
---|
| 279 | inline G4double GetRange(G4double& kineticEnergy, const G4MaterialCutsCouple*); |
---|
| 280 | inline G4double GetCSDARange(G4double& kineticEnergy, const G4MaterialCutsCouple*); |
---|
| 281 | inline G4double GetRangeForLoss(G4double& kineticEnergy, const G4MaterialCutsCouple*); |
---|
| 282 | inline G4double GetKineticEnergy(G4double& range, const G4MaterialCutsCouple*); |
---|
| 283 | inline G4double GetLambda(G4double& kineticEnergy, const G4MaterialCutsCouple*); |
---|
| 284 | |
---|
| 285 | inline G4bool TablesAreBuilt() const; |
---|
| 286 | |
---|
| 287 | //------------------------------------------------------------------------ |
---|
| 288 | // Define and access particle type |
---|
| 289 | //------------------------------------------------------------------------ |
---|
| 290 | |
---|
| 291 | inline void SetBaseParticle(const G4ParticleDefinition* p); |
---|
| 292 | inline const G4ParticleDefinition* Particle() const; |
---|
| 293 | inline const G4ParticleDefinition* BaseParticle() const; |
---|
| 294 | inline const G4ParticleDefinition* SecondaryParticle() const; |
---|
| 295 | |
---|
| 296 | //------------------------------------------------------------------------ |
---|
| 297 | // Specific methods to set, access, modify models |
---|
| 298 | //------------------------------------------------------------------------ |
---|
| 299 | |
---|
| 300 | // Add EM model coupled with fluctuation model for the region |
---|
| 301 | inline void AddEmModel(G4int, G4VEmModel*, G4VEmFluctuationModel* fluc = 0, |
---|
| 302 | const G4Region* region = 0); |
---|
| 303 | |
---|
| 304 | // Assign a model to a process |
---|
| 305 | inline void SetEmModel(G4VEmModel*, G4int index=1); |
---|
| 306 | |
---|
| 307 | // return the assigned model |
---|
| 308 | inline G4VEmModel* EmModel(G4int index=1); |
---|
| 309 | |
---|
| 310 | // Assign a fluctuation model to a process |
---|
| 311 | inline void SetFluctModel(G4VEmFluctuationModel*); |
---|
| 312 | |
---|
| 313 | // return the assigned fluctuation model |
---|
| 314 | inline G4VEmFluctuationModel* FluctModel(); |
---|
| 315 | |
---|
| 316 | // Define new energy range for the model identified by the name |
---|
| 317 | inline void UpdateEmModel(const G4String&, G4double, G4double); |
---|
| 318 | |
---|
| 319 | // Access to models |
---|
| 320 | inline G4VEmModel* GetModelByIndex(G4int idx = 0); |
---|
| 321 | |
---|
| 322 | inline G4int NumberOfModels(); |
---|
| 323 | |
---|
| 324 | //------------------------------------------------------------------------ |
---|
| 325 | // Get/set parameters used for simulation of energy loss |
---|
| 326 | //------------------------------------------------------------------------ |
---|
| 327 | |
---|
| 328 | inline void SetLossFluctuations(G4bool val); |
---|
| 329 | inline void SetRandomStep(G4bool val); |
---|
| 330 | inline void SetIntegral(G4bool val); |
---|
| 331 | inline G4bool IsIntegral() const; |
---|
| 332 | |
---|
| 333 | // Set/Get flag "isIonisation" |
---|
| 334 | inline void SetIonisation(G4bool val); |
---|
| 335 | inline G4bool IsIonisationProcess() const; |
---|
| 336 | |
---|
| 337 | // Redefine parameteters for stepping control |
---|
| 338 | // |
---|
| 339 | inline void SetLinearLossLimit(G4double val); |
---|
| 340 | inline void SetMinSubRange(G4double val); |
---|
| 341 | inline void SetStepFunction(G4double v1, G4double v2); |
---|
| 342 | inline void SetLambdaFactor(G4double val); |
---|
| 343 | |
---|
| 344 | |
---|
| 345 | // Add subcutoff option for the region |
---|
| 346 | void ActivateSubCutoff(G4bool val, const G4Region* region = 0); |
---|
| 347 | |
---|
| 348 | inline G4int NumberOfSubCutoffRegions() const; |
---|
| 349 | |
---|
| 350 | // Activate deexcitation code |
---|
| 351 | virtual void ActivateDeexcitation(G4bool, const G4Region* region = 0); |
---|
| 352 | |
---|
| 353 | //------------------------------------------------------------------------ |
---|
| 354 | // Run time method for simulation of ionisation |
---|
| 355 | //------------------------------------------------------------------------ |
---|
| 356 | |
---|
| 357 | inline G4double SampleRange(); |
---|
| 358 | |
---|
| 359 | inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy, size_t& idx) const; |
---|
| 360 | |
---|
| 361 | |
---|
| 362 | // Set scaling parameters |
---|
| 363 | inline void SetDynamicMassCharge(G4double massratio, G4double charge2ratio); |
---|
| 364 | |
---|
| 365 | // Helper functions |
---|
| 366 | inline G4double MeanFreePath(const G4Track& track); |
---|
| 367 | |
---|
| 368 | inline G4double ContinuousStepLimit(const G4Track& track, |
---|
| 369 | G4double previousStepSize, |
---|
| 370 | G4double currentMinimumStep, |
---|
| 371 | G4double& currentSafety); |
---|
| 372 | |
---|
| 373 | protected: |
---|
| 374 | |
---|
| 375 | G4PhysicsVector* LambdaPhysicsVector(const G4MaterialCutsCouple*, |
---|
| 376 | G4double cut); |
---|
| 377 | |
---|
| 378 | inline virtual void InitialiseMassCharge(const G4Track&); |
---|
| 379 | |
---|
| 380 | inline void SetParticle(const G4ParticleDefinition* p); |
---|
| 381 | |
---|
| 382 | inline void SetSecondaryParticle(const G4ParticleDefinition* p); |
---|
| 383 | |
---|
| 384 | inline G4VEmModel* SelectModel(G4double kinEnergy); |
---|
| 385 | |
---|
| 386 | inline size_t CurrentMaterialCutsCoupleIndex() const; |
---|
| 387 | |
---|
| 388 | inline G4double GetCurrentRange() const; |
---|
| 389 | |
---|
| 390 | private: |
---|
| 391 | |
---|
| 392 | // Clear tables |
---|
| 393 | void Clear(); |
---|
| 394 | |
---|
| 395 | inline void InitialiseStep(const G4Track&); |
---|
| 396 | |
---|
| 397 | inline void DefineMaterial(const G4MaterialCutsCouple* couple); |
---|
| 398 | |
---|
| 399 | // Returnd values for scaled energy and base particles mass |
---|
| 400 | // |
---|
| 401 | inline G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy); |
---|
| 402 | inline G4double GetSubDEDXForScaledEnergy(G4double scaledKinEnergy); |
---|
| 403 | inline G4double GetIonisationForScaledEnergy(G4double scaledKinEnergy); |
---|
| 404 | inline G4double GetSubIonisationForScaledEnergy(G4double scaledKinEnergy); |
---|
| 405 | inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy); |
---|
| 406 | inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinEnergy); |
---|
| 407 | inline G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy); |
---|
| 408 | inline G4double ScaledKinEnergyForLoss(G4double range); |
---|
| 409 | inline void ComputeLambdaForScaledEnergy(G4double scaledKinEnergy); |
---|
| 410 | |
---|
| 411 | // hide assignment operator |
---|
| 412 | |
---|
| 413 | G4VEnergyLossProcess(G4VEnergyLossProcess &); |
---|
| 414 | G4VEnergyLossProcess & operator=(const G4VEnergyLossProcess &right); |
---|
| 415 | |
---|
| 416 | // ===================================================================== |
---|
| 417 | |
---|
| 418 | protected: |
---|
| 419 | |
---|
| 420 | G4ParticleChangeForLoss fParticleChange; |
---|
| 421 | |
---|
| 422 | private: |
---|
| 423 | |
---|
| 424 | G4EmModelManager* modelManager; |
---|
| 425 | std::vector<G4VEmModel*> emModels; |
---|
| 426 | G4VEmFluctuationModel* fluctModel; |
---|
| 427 | std::vector<const G4Region*> scoffRegions; |
---|
| 428 | G4int nSCoffRegions; |
---|
| 429 | G4int* idxSCoffRegions; |
---|
| 430 | std::vector<G4DynamicParticle*> secParticles; |
---|
| 431 | std::vector<G4Track*> scTracks; |
---|
| 432 | std::vector<G4VEnergyLossProcess*> scProcesses; |
---|
| 433 | G4int nProcesses; |
---|
| 434 | |
---|
| 435 | // tables and vectors |
---|
| 436 | G4PhysicsTable* theDEDXTable; |
---|
| 437 | G4PhysicsTable* theDEDXSubTable; |
---|
| 438 | G4PhysicsTable* theDEDXunRestrictedTable; |
---|
| 439 | G4PhysicsTable* theIonisationTable; |
---|
| 440 | G4PhysicsTable* theIonisationSubTable; |
---|
| 441 | G4PhysicsTable* theRangeTableForLoss; |
---|
| 442 | G4PhysicsTable* theCSDARangeTable; |
---|
| 443 | G4PhysicsTable* theSecondaryRangeTable; |
---|
| 444 | G4PhysicsTable* theInverseRangeTable; |
---|
| 445 | G4PhysicsTable* theLambdaTable; |
---|
| 446 | G4PhysicsTable* theSubLambdaTable; |
---|
| 447 | G4double* theDEDXAtMaxEnergy; |
---|
| 448 | G4double* theRangeAtMaxEnergy; |
---|
| 449 | G4double* theEnergyOfCrossSectionMax; |
---|
| 450 | G4double* theCrossSectionMax; |
---|
| 451 | |
---|
| 452 | const G4DataVector* theCuts; |
---|
| 453 | const G4DataVector* theSubCuts; |
---|
| 454 | |
---|
| 455 | G4SafetyHelper* safetyHelper; |
---|
| 456 | |
---|
| 457 | const G4ParticleDefinition* particle; |
---|
| 458 | const G4ParticleDefinition* baseParticle; |
---|
| 459 | const G4ParticleDefinition* secondaryParticle; |
---|
| 460 | const G4ParticleDefinition* theElectron; |
---|
| 461 | const G4ParticleDefinition* thePositron; |
---|
| 462 | |
---|
| 463 | G4PhysicsVector* vstrag; |
---|
| 464 | |
---|
| 465 | // cash |
---|
| 466 | const G4Material* currentMaterial; |
---|
| 467 | const G4MaterialCutsCouple* currentCouple; |
---|
| 468 | size_t currentMaterialIndex; |
---|
| 469 | |
---|
| 470 | G4int nBins; |
---|
| 471 | G4int nBinsCSDA; |
---|
| 472 | G4int nWarnings; |
---|
| 473 | |
---|
| 474 | G4double lowestKinEnergy; |
---|
| 475 | G4double minKinEnergy; |
---|
| 476 | G4double maxKinEnergy; |
---|
| 477 | G4double maxKinEnergyCSDA; |
---|
| 478 | |
---|
| 479 | G4double massRatio; |
---|
| 480 | G4double reduceFactor; |
---|
| 481 | G4double chargeSquare; |
---|
| 482 | G4double chargeSqRatio; |
---|
| 483 | |
---|
| 484 | G4double preStepLambda; |
---|
| 485 | G4double fRange; |
---|
| 486 | G4double preStepKinEnergy; |
---|
| 487 | G4double preStepScaledEnergy; |
---|
| 488 | G4double linLossLimit; |
---|
| 489 | G4double minSubRange; |
---|
| 490 | G4double dRoverRange; |
---|
| 491 | G4double finalRange; |
---|
| 492 | G4double lambdaFactor; |
---|
| 493 | G4double mfpKinEnergy; |
---|
| 494 | |
---|
| 495 | G4GPILSelection aGPILSelection; |
---|
| 496 | |
---|
| 497 | G4bool lossFluctuationFlag; |
---|
| 498 | G4bool rndmStepFlag; |
---|
| 499 | G4bool tablesAreBuilt; |
---|
| 500 | G4bool integral; |
---|
| 501 | G4bool isIonisation; |
---|
| 502 | G4bool useSubCutoff; |
---|
| 503 | }; |
---|
| 504 | |
---|
| 505 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 506 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 507 | |
---|
| 508 | inline void G4VEnergyLossProcess::DefineMaterial( |
---|
| 509 | const G4MaterialCutsCouple* couple) |
---|
| 510 | { |
---|
| 511 | if(couple != currentCouple) { |
---|
| 512 | currentCouple = couple; |
---|
| 513 | currentMaterial = couple->GetMaterial(); |
---|
| 514 | currentMaterialIndex = couple->GetIndex(); |
---|
| 515 | mfpKinEnergy = DBL_MAX; |
---|
| 516 | } |
---|
| 517 | } |
---|
| 518 | |
---|
| 519 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 520 | |
---|
| 521 | inline void G4VEnergyLossProcess::InitialiseStep(const G4Track& track) |
---|
| 522 | { |
---|
| 523 | InitialiseMassCharge(track); |
---|
| 524 | preStepKinEnergy = track.GetKineticEnergy(); |
---|
| 525 | preStepScaledEnergy = preStepKinEnergy*massRatio; |
---|
| 526 | DefineMaterial(track.GetMaterialCutsCouple()); |
---|
| 527 | if (theNumberOfInteractionLengthLeft < 0.0) mfpKinEnergy = DBL_MAX; |
---|
| 528 | } |
---|
| 529 | |
---|
| 530 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 531 | |
---|
| 532 | inline void G4VEnergyLossProcess::InitialiseMassCharge(const G4Track&) |
---|
| 533 | {} |
---|
| 534 | |
---|
| 535 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 536 | |
---|
| 537 | inline G4double G4VEnergyLossProcess::GetDEDX(G4double& kineticEnergy, |
---|
| 538 | const G4MaterialCutsCouple* couple) |
---|
| 539 | { |
---|
| 540 | DefineMaterial(couple); |
---|
| 541 | return GetDEDXForScaledEnergy(kineticEnergy*massRatio); |
---|
| 542 | } |
---|
| 543 | |
---|
| 544 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 545 | |
---|
| 546 | inline G4double G4VEnergyLossProcess::GetDEDXForSubsec(G4double& kineticEnergy, |
---|
| 547 | const G4MaterialCutsCouple* couple) |
---|
| 548 | { |
---|
| 549 | DefineMaterial(couple); |
---|
| 550 | return GetSubDEDXForScaledEnergy(kineticEnergy*massRatio); |
---|
| 551 | } |
---|
| 552 | |
---|
| 553 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 554 | |
---|
| 555 | inline G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e) |
---|
| 556 | { |
---|
| 557 | G4bool b; |
---|
| 558 | G4double x = |
---|
| 559 | ((*theDEDXTable)[currentMaterialIndex]->GetValue(e, b))*chargeSqRatio; |
---|
| 560 | if(e < minKinEnergy) x *= std::sqrt(e/minKinEnergy); |
---|
| 561 | return x; |
---|
| 562 | } |
---|
| 563 | |
---|
| 564 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 565 | |
---|
| 566 | inline G4double G4VEnergyLossProcess::GetSubDEDXForScaledEnergy(G4double e) |
---|
| 567 | { |
---|
| 568 | G4bool b; |
---|
| 569 | G4double x = |
---|
| 570 | ((*theDEDXSubTable)[currentMaterialIndex]->GetValue(e, b))*chargeSqRatio; |
---|
| 571 | if(e < minKinEnergy) x *= std::sqrt(e/minKinEnergy); |
---|
| 572 | return x; |
---|
| 573 | } |
---|
| 574 | |
---|
| 575 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 576 | |
---|
| 577 | inline G4double G4VEnergyLossProcess::GetIonisationForScaledEnergy(G4double e) |
---|
| 578 | { |
---|
| 579 | G4bool b; |
---|
| 580 | G4double x = 0.0; |
---|
| 581 | // if(theIonisationTable) { |
---|
| 582 | x = ((*theIonisationTable)[currentMaterialIndex]->GetValue(e, b)) |
---|
| 583 | *chargeSqRatio; |
---|
| 584 | if(e < minKinEnergy) x *= std::sqrt(e/minKinEnergy); |
---|
| 585 | //} |
---|
| 586 | return x; |
---|
| 587 | } |
---|
| 588 | |
---|
| 589 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 590 | |
---|
| 591 | inline |
---|
| 592 | G4double G4VEnergyLossProcess::GetSubIonisationForScaledEnergy(G4double e) |
---|
| 593 | { |
---|
| 594 | G4bool b; |
---|
| 595 | G4double x = 0.0; |
---|
| 596 | //if(theIonisationSubTable) { |
---|
| 597 | x = ((*theIonisationSubTable)[currentMaterialIndex]->GetValue(e, b)) |
---|
| 598 | *chargeSqRatio; |
---|
| 599 | if(e < minKinEnergy) x *= std::sqrt(e/minKinEnergy); |
---|
| 600 | //} |
---|
| 601 | return x; |
---|
| 602 | } |
---|
| 603 | |
---|
| 604 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 605 | |
---|
| 606 | inline G4double G4VEnergyLossProcess::GetRange(G4double& kineticEnergy, |
---|
| 607 | const G4MaterialCutsCouple* couple) |
---|
| 608 | { |
---|
| 609 | G4double x = fRange; |
---|
| 610 | if(kineticEnergy != preStepKinEnergy || couple != currentCouple) { |
---|
| 611 | DefineMaterial(couple); |
---|
| 612 | if(theCSDARangeTable) |
---|
| 613 | x = GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio) |
---|
| 614 | * reduceFactor; |
---|
| 615 | else if(theRangeTableForLoss) |
---|
| 616 | x = GetScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor; |
---|
| 617 | } |
---|
| 618 | return x; |
---|
| 619 | } |
---|
| 620 | |
---|
| 621 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 622 | |
---|
| 623 | inline G4double G4VEnergyLossProcess::GetCSDARange( |
---|
| 624 | G4double& kineticEnergy, const G4MaterialCutsCouple* couple) |
---|
| 625 | { |
---|
| 626 | DefineMaterial(couple); |
---|
| 627 | G4double x = DBL_MAX; |
---|
| 628 | if(theCSDARangeTable) |
---|
| 629 | x = GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio) |
---|
| 630 | * reduceFactor; |
---|
| 631 | return x; |
---|
| 632 | } |
---|
| 633 | |
---|
| 634 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 635 | |
---|
| 636 | inline G4double G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy( |
---|
| 637 | G4double e) |
---|
| 638 | { |
---|
| 639 | G4bool b; |
---|
| 640 | G4double x; |
---|
| 641 | |
---|
| 642 | if (e < maxKinEnergyCSDA) { |
---|
| 643 | x = ((*theCSDARangeTable)[currentMaterialIndex])->GetValue(e, b); |
---|
| 644 | if(e < minKinEnergy) x *= std::sqrt(e/minKinEnergy); |
---|
| 645 | } else { |
---|
| 646 | x = theRangeAtMaxEnergy[currentMaterialIndex] + |
---|
| 647 | (e - maxKinEnergyCSDA)/theDEDXAtMaxEnergy[currentMaterialIndex]; |
---|
| 648 | } |
---|
| 649 | return x; |
---|
| 650 | } |
---|
| 651 | |
---|
| 652 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 653 | |
---|
| 654 | inline G4double G4VEnergyLossProcess::GetRangeForLoss( |
---|
| 655 | G4double& kineticEnergy, |
---|
| 656 | const G4MaterialCutsCouple* couple) |
---|
| 657 | { |
---|
| 658 | DefineMaterial(couple); |
---|
| 659 | G4double x = DBL_MAX; |
---|
| 660 | if(theRangeTableForLoss) |
---|
| 661 | x = GetScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor; |
---|
| 662 | // G4cout << "Range from " << GetProcessName() |
---|
| 663 | // << " e= " << kineticEnergy << " r= " << x << G4endl; |
---|
| 664 | return x; |
---|
| 665 | } |
---|
| 666 | |
---|
| 667 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 668 | |
---|
| 669 | inline G4double G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e) |
---|
| 670 | { |
---|
| 671 | G4bool b; |
---|
| 672 | G4double x = ((*theRangeTableForLoss)[currentMaterialIndex])->GetValue(e, b); |
---|
| 673 | if(e < minKinEnergy) x *= std::sqrt(e/minKinEnergy); |
---|
| 674 | return x; |
---|
| 675 | } |
---|
| 676 | |
---|
| 677 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 678 | |
---|
| 679 | inline G4double G4VEnergyLossProcess::GetKineticEnergy( |
---|
| 680 | G4double& range, |
---|
| 681 | const G4MaterialCutsCouple* couple) |
---|
| 682 | { |
---|
| 683 | DefineMaterial(couple); |
---|
| 684 | G4double r = range/reduceFactor; |
---|
| 685 | G4double e = ScaledKinEnergyForLoss(r)/massRatio; |
---|
| 686 | return e; |
---|
| 687 | } |
---|
| 688 | |
---|
| 689 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 690 | |
---|
| 691 | inline G4double G4VEnergyLossProcess::ScaledKinEnergyForLoss(G4double r) |
---|
| 692 | { |
---|
| 693 | G4PhysicsVector* v = (*theInverseRangeTable)[currentMaterialIndex]; |
---|
| 694 | G4double rmin = v->GetLowEdgeEnergy(0); |
---|
| 695 | G4double e = 0.0; |
---|
| 696 | if(r >= rmin) { |
---|
| 697 | G4bool b; |
---|
| 698 | e = v->GetValue(r, b); |
---|
| 699 | } else if(r > 0.0) { |
---|
| 700 | G4double x = r/rmin; |
---|
| 701 | e = minKinEnergy*x*x; |
---|
| 702 | } |
---|
| 703 | return e; |
---|
| 704 | } |
---|
| 705 | |
---|
| 706 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 707 | |
---|
| 708 | inline G4double G4VEnergyLossProcess::GetLambda(G4double& kineticEnergy, |
---|
| 709 | const G4MaterialCutsCouple* couple) |
---|
| 710 | { |
---|
| 711 | DefineMaterial(couple); |
---|
| 712 | G4double x = 0.0; |
---|
| 713 | if(theLambdaTable) x = GetLambdaForScaledEnergy(kineticEnergy*massRatio); |
---|
| 714 | return x; |
---|
| 715 | } |
---|
| 716 | |
---|
| 717 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 718 | |
---|
| 719 | inline G4double G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e) |
---|
| 720 | { |
---|
| 721 | G4bool b; |
---|
| 722 | return |
---|
| 723 | chargeSqRatio*(((*theLambdaTable)[currentMaterialIndex])->GetValue(e, b)); |
---|
| 724 | } |
---|
| 725 | |
---|
| 726 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 727 | |
---|
| 728 | inline void G4VEnergyLossProcess::ComputeLambdaForScaledEnergy(G4double e) |
---|
| 729 | { |
---|
| 730 | mfpKinEnergy = theEnergyOfCrossSectionMax[currentMaterialIndex]; |
---|
| 731 | if (e <= mfpKinEnergy) { |
---|
| 732 | preStepLambda = GetLambdaForScaledEnergy(e); |
---|
| 733 | |
---|
| 734 | } else { |
---|
| 735 | G4double e1 = e*lambdaFactor; |
---|
| 736 | if(e1 > mfpKinEnergy) { |
---|
| 737 | preStepLambda = GetLambdaForScaledEnergy(e); |
---|
| 738 | G4double preStepLambda1 = GetLambdaForScaledEnergy(e1); |
---|
| 739 | if(preStepLambda1 > preStepLambda) { |
---|
| 740 | mfpKinEnergy = e1; |
---|
| 741 | preStepLambda = preStepLambda1; |
---|
| 742 | } |
---|
| 743 | } else { |
---|
| 744 | preStepLambda = chargeSqRatio*theCrossSectionMax[currentMaterialIndex]; |
---|
| 745 | } |
---|
| 746 | } |
---|
| 747 | } |
---|
| 748 | |
---|
| 749 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 750 | |
---|
| 751 | inline G4double G4VEnergyLossProcess::ContinuousStepLimit( |
---|
| 752 | const G4Track& track, G4double x, G4double y, G4double& z) |
---|
| 753 | { |
---|
| 754 | G4GPILSelection sel; |
---|
| 755 | return AlongStepGetPhysicalInteractionLength(track, x, y, z, &sel); |
---|
| 756 | } |
---|
| 757 | |
---|
| 758 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 759 | |
---|
| 760 | inline G4double G4VEnergyLossProcess::SampleRange() |
---|
| 761 | { |
---|
| 762 | G4double e = amu_c2*preStepKinEnergy/particle->GetPDGMass(); |
---|
| 763 | G4bool b; |
---|
| 764 | G4double s = fRange*std::pow(10.,vstrag->GetValue(e,b)); |
---|
| 765 | G4double x = fRange + G4RandGauss::shoot(0.0,s); |
---|
| 766 | if(x > 0.0) fRange = x; |
---|
| 767 | return fRange; |
---|
| 768 | } |
---|
| 769 | |
---|
| 770 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 771 | |
---|
| 772 | inline G4double G4VEnergyLossProcess::MeanFreePath(const G4Track& track) |
---|
| 773 | { |
---|
| 774 | DefineMaterial(track.GetMaterialCutsCouple()); |
---|
| 775 | preStepLambda = GetLambdaForScaledEnergy(track.GetKineticEnergy()*massRatio); |
---|
| 776 | G4double x = DBL_MAX; |
---|
| 777 | if(DBL_MIN < preStepLambda) x = 1.0/preStepLambda; |
---|
| 778 | return x; |
---|
| 779 | } |
---|
| 780 | |
---|
| 781 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 782 | |
---|
| 783 | inline G4double G4VEnergyLossProcess::MinPrimaryEnergy( |
---|
| 784 | const G4ParticleDefinition*, const G4Material*, G4double cut) |
---|
| 785 | { |
---|
| 786 | return cut; |
---|
| 787 | } |
---|
| 788 | |
---|
| 789 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 790 | |
---|
| 791 | inline G4VEmModel* G4VEnergyLossProcess::SelectModel(G4double kinEnergy) |
---|
| 792 | { |
---|
| 793 | return modelManager->SelectModel(kinEnergy, currentMaterialIndex); |
---|
| 794 | } |
---|
| 795 | |
---|
| 796 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 797 | |
---|
| 798 | inline G4VEmModel* G4VEnergyLossProcess::SelectModelForMaterial( |
---|
| 799 | G4double kinEnergy, size_t& idx) const |
---|
| 800 | { |
---|
| 801 | return modelManager->SelectModel(kinEnergy, idx); |
---|
| 802 | } |
---|
| 803 | |
---|
| 804 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 805 | |
---|
| 806 | inline const G4ParticleDefinition* G4VEnergyLossProcess::Particle() const |
---|
| 807 | { |
---|
| 808 | return particle; |
---|
| 809 | } |
---|
| 810 | |
---|
| 811 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 812 | |
---|
| 813 | inline const G4ParticleDefinition* G4VEnergyLossProcess::BaseParticle() const |
---|
| 814 | { |
---|
| 815 | return baseParticle; |
---|
| 816 | } |
---|
| 817 | |
---|
| 818 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 819 | |
---|
| 820 | inline const G4ParticleDefinition* G4VEnergyLossProcess::SecondaryParticle() const |
---|
| 821 | { |
---|
| 822 | return secondaryParticle; |
---|
| 823 | } |
---|
| 824 | |
---|
| 825 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 826 | |
---|
| 827 | inline void G4VEnergyLossProcess::CorrectionsAlongStep( |
---|
| 828 | const G4MaterialCutsCouple*, |
---|
| 829 | const G4DynamicParticle*, |
---|
| 830 | G4double&, |
---|
| 831 | G4double&) |
---|
| 832 | {} |
---|
| 833 | |
---|
| 834 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 835 | |
---|
| 836 | inline G4PhysicsTable* G4VEnergyLossProcess::DEDXTable() const |
---|
| 837 | { |
---|
| 838 | return theDEDXTable; |
---|
| 839 | } |
---|
| 840 | |
---|
| 841 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 842 | |
---|
| 843 | inline G4PhysicsTable* G4VEnergyLossProcess::DEDXTableForSubsec() const |
---|
| 844 | { |
---|
| 845 | return theDEDXSubTable; |
---|
| 846 | } |
---|
| 847 | |
---|
| 848 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 849 | |
---|
| 850 | inline G4PhysicsTable* G4VEnergyLossProcess::DEDXunRestrictedTable() const |
---|
| 851 | { |
---|
| 852 | return theDEDXunRestrictedTable; |
---|
| 853 | } |
---|
| 854 | |
---|
| 855 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 856 | |
---|
| 857 | inline G4PhysicsTable* G4VEnergyLossProcess::IonisationTable() const |
---|
| 858 | { |
---|
| 859 | G4PhysicsTable* t = theDEDXTable; |
---|
| 860 | if(theIonisationTable) t = theIonisationTable; |
---|
| 861 | return t; |
---|
| 862 | } |
---|
| 863 | |
---|
| 864 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 865 | |
---|
| 866 | inline G4PhysicsTable* G4VEnergyLossProcess::IonisationTableForSubsec() const |
---|
| 867 | { |
---|
| 868 | G4PhysicsTable* t = theDEDXSubTable; |
---|
| 869 | if(theIonisationSubTable) t = theIonisationSubTable; |
---|
| 870 | return t; |
---|
| 871 | } |
---|
| 872 | |
---|
| 873 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 874 | |
---|
| 875 | inline G4PhysicsTable* G4VEnergyLossProcess::CSDARangeTable() const |
---|
| 876 | { |
---|
| 877 | return theCSDARangeTable; |
---|
| 878 | } |
---|
| 879 | |
---|
| 880 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 881 | |
---|
| 882 | inline G4PhysicsTable* G4VEnergyLossProcess::RangeTableForLoss() const |
---|
| 883 | { |
---|
| 884 | return theRangeTableForLoss; |
---|
| 885 | } |
---|
| 886 | |
---|
| 887 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 888 | |
---|
| 889 | inline G4PhysicsTable* G4VEnergyLossProcess::InverseRangeTable() const |
---|
| 890 | { |
---|
| 891 | return theInverseRangeTable; |
---|
| 892 | } |
---|
| 893 | |
---|
| 894 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 895 | |
---|
| 896 | inline G4PhysicsTable* G4VEnergyLossProcess::LambdaTable() |
---|
| 897 | { |
---|
| 898 | return theLambdaTable; |
---|
| 899 | } |
---|
| 900 | |
---|
| 901 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 902 | |
---|
| 903 | inline G4PhysicsTable* G4VEnergyLossProcess::SubLambdaTable() |
---|
| 904 | { |
---|
| 905 | return theSubLambdaTable; |
---|
| 906 | } |
---|
| 907 | |
---|
| 908 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 909 | |
---|
| 910 | inline G4bool G4VEnergyLossProcess::IsIntegral() const |
---|
| 911 | { |
---|
| 912 | return integral; |
---|
| 913 | } |
---|
| 914 | |
---|
| 915 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 916 | |
---|
| 917 | inline size_t G4VEnergyLossProcess::CurrentMaterialCutsCoupleIndex() const |
---|
| 918 | { |
---|
| 919 | return currentMaterialIndex; |
---|
| 920 | } |
---|
| 921 | |
---|
| 922 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 923 | |
---|
| 924 | inline void G4VEnergyLossProcess::SetDynamicMassCharge(G4double massratio, |
---|
| 925 | G4double charge2ratio) |
---|
| 926 | { |
---|
| 927 | massRatio = massratio; |
---|
| 928 | chargeSqRatio = charge2ratio; |
---|
| 929 | chargeSquare = charge2ratio*eplus*eplus; |
---|
| 930 | if(chargeSqRatio > 0.0) reduceFactor = 1.0/(chargeSqRatio*massRatio); |
---|
| 931 | } |
---|
| 932 | |
---|
| 933 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 934 | |
---|
| 935 | inline G4double G4VEnergyLossProcess::GetCurrentRange() const |
---|
| 936 | { |
---|
| 937 | return fRange; |
---|
| 938 | } |
---|
| 939 | |
---|
| 940 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 941 | |
---|
| 942 | void G4VEnergyLossProcess::AddEmModel(G4int order, G4VEmModel* p, |
---|
| 943 | G4VEmFluctuationModel* fluc, |
---|
| 944 | const G4Region* region) |
---|
| 945 | { |
---|
| 946 | modelManager->AddEmModel(order, p, fluc, region); |
---|
| 947 | if(p) p->SetParticleChange(pParticleChange, fluc); |
---|
| 948 | } |
---|
| 949 | |
---|
| 950 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 951 | |
---|
| 952 | inline G4VEmModel* G4VEnergyLossProcess::GetModelByIndex(G4int idx) |
---|
| 953 | { |
---|
| 954 | return modelManager->GetModel(idx); |
---|
| 955 | } |
---|
| 956 | |
---|
| 957 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 958 | |
---|
| 959 | inline G4int G4VEnergyLossProcess::NumberOfModels() |
---|
| 960 | { |
---|
| 961 | return modelManager->NumberOfModels(); |
---|
| 962 | } |
---|
| 963 | |
---|
| 964 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 965 | |
---|
| 966 | inline void G4VEnergyLossProcess::SetEmModel(G4VEmModel* p, G4int index) |
---|
| 967 | { |
---|
| 968 | G4int n = emModels.size(); |
---|
| 969 | if(index >= n) for(G4int i=n; i<index+1; i++) {emModels.push_back(0);} |
---|
| 970 | emModels[index] = p; |
---|
| 971 | } |
---|
| 972 | |
---|
| 973 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 974 | |
---|
| 975 | inline G4VEmModel* G4VEnergyLossProcess::EmModel(G4int index) |
---|
| 976 | { |
---|
| 977 | G4VEmModel* p = 0; |
---|
| 978 | if(index >= 0 && index < G4int(emModels.size())) p = emModels[index]; |
---|
| 979 | return p; |
---|
| 980 | } |
---|
| 981 | |
---|
| 982 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 983 | |
---|
| 984 | inline void G4VEnergyLossProcess::SetFluctModel(G4VEmFluctuationModel* p) |
---|
| 985 | { |
---|
| 986 | fluctModel = p; |
---|
| 987 | } |
---|
| 988 | |
---|
| 989 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 990 | |
---|
| 991 | inline G4VEmFluctuationModel* G4VEnergyLossProcess::FluctModel() |
---|
| 992 | { |
---|
| 993 | return fluctModel; |
---|
| 994 | } |
---|
| 995 | |
---|
| 996 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 997 | |
---|
| 998 | inline void G4VEnergyLossProcess::UpdateEmModel(const G4String& nam, |
---|
| 999 | G4double emin, G4double emax) |
---|
| 1000 | { |
---|
| 1001 | modelManager->UpdateEmModel(nam, emin, emax); |
---|
| 1002 | } |
---|
| 1003 | |
---|
| 1004 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1005 | |
---|
| 1006 | inline void G4VEnergyLossProcess::SetIntegral(G4bool val) |
---|
| 1007 | { |
---|
| 1008 | integral = val; |
---|
| 1009 | } |
---|
| 1010 | |
---|
| 1011 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1012 | |
---|
| 1013 | inline void G4VEnergyLossProcess::SetParticle(const G4ParticleDefinition* p) |
---|
| 1014 | { |
---|
| 1015 | particle = p; |
---|
| 1016 | } |
---|
| 1017 | |
---|
| 1018 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1019 | |
---|
| 1020 | inline void G4VEnergyLossProcess::SetBaseParticle(const G4ParticleDefinition* p) |
---|
| 1021 | { |
---|
| 1022 | baseParticle = p; |
---|
| 1023 | } |
---|
| 1024 | |
---|
| 1025 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1026 | |
---|
| 1027 | inline void G4VEnergyLossProcess::SetSecondaryParticle(const G4ParticleDefinition* p) |
---|
| 1028 | { |
---|
| 1029 | secondaryParticle = p; |
---|
| 1030 | } |
---|
| 1031 | |
---|
| 1032 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1033 | |
---|
| 1034 | inline void G4VEnergyLossProcess::SetLinearLossLimit(G4double val) |
---|
| 1035 | { |
---|
| 1036 | linLossLimit = val; |
---|
| 1037 | } |
---|
| 1038 | |
---|
| 1039 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1040 | |
---|
| 1041 | inline void G4VEnergyLossProcess::SetLossFluctuations(G4bool val) |
---|
| 1042 | { |
---|
| 1043 | lossFluctuationFlag = val; |
---|
| 1044 | } |
---|
| 1045 | |
---|
| 1046 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1047 | |
---|
| 1048 | inline void G4VEnergyLossProcess::SetRandomStep(G4bool val) |
---|
| 1049 | { |
---|
| 1050 | rndmStepFlag = val; |
---|
| 1051 | } |
---|
| 1052 | |
---|
| 1053 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1054 | |
---|
| 1055 | inline void G4VEnergyLossProcess::SetMinSubRange(G4double val) |
---|
| 1056 | { |
---|
| 1057 | minSubRange = val; |
---|
| 1058 | } |
---|
| 1059 | |
---|
| 1060 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1061 | |
---|
| 1062 | inline G4bool G4VEnergyLossProcess::TablesAreBuilt() const |
---|
| 1063 | { |
---|
| 1064 | return tablesAreBuilt; |
---|
| 1065 | } |
---|
| 1066 | |
---|
| 1067 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1068 | |
---|
| 1069 | inline G4int G4VEnergyLossProcess::NumberOfSubCutoffRegions() const |
---|
| 1070 | { |
---|
| 1071 | return nSCoffRegions; |
---|
| 1072 | } |
---|
| 1073 | |
---|
| 1074 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1075 | |
---|
| 1076 | inline void G4VEnergyLossProcess::SetDEDXBinning(G4int nbins) |
---|
| 1077 | { |
---|
| 1078 | nBins = nbins; |
---|
| 1079 | } |
---|
| 1080 | |
---|
| 1081 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1082 | |
---|
| 1083 | inline void G4VEnergyLossProcess::SetLambdaBinning(G4int nbins) |
---|
| 1084 | { |
---|
| 1085 | nBins = nbins; |
---|
| 1086 | } |
---|
| 1087 | |
---|
| 1088 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1089 | |
---|
| 1090 | inline void G4VEnergyLossProcess::SetDEDXBinningForCSDARange(G4int nbins) |
---|
| 1091 | { |
---|
| 1092 | nBinsCSDA = nbins; |
---|
| 1093 | } |
---|
| 1094 | |
---|
| 1095 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1096 | |
---|
| 1097 | inline G4double G4VEnergyLossProcess::MinKinEnergy() const |
---|
| 1098 | { |
---|
| 1099 | return minKinEnergy; |
---|
| 1100 | } |
---|
| 1101 | |
---|
| 1102 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1103 | |
---|
| 1104 | inline void G4VEnergyLossProcess::SetMinKinEnergy(G4double e) |
---|
| 1105 | { |
---|
| 1106 | minKinEnergy = e; |
---|
| 1107 | } |
---|
| 1108 | |
---|
| 1109 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1110 | |
---|
| 1111 | inline void G4VEnergyLossProcess::SetMaxKinEnergy(G4double e) |
---|
| 1112 | { |
---|
| 1113 | maxKinEnergy = e; |
---|
| 1114 | if(e < maxKinEnergyCSDA) maxKinEnergyCSDA = e; |
---|
| 1115 | } |
---|
| 1116 | |
---|
| 1117 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1118 | |
---|
| 1119 | inline void G4VEnergyLossProcess::SetMaxKinEnergyForCSDARange(G4double e) |
---|
| 1120 | { |
---|
| 1121 | maxKinEnergyCSDA = e; |
---|
| 1122 | } |
---|
| 1123 | |
---|
| 1124 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1125 | |
---|
| 1126 | inline G4double G4VEnergyLossProcess::MaxKinEnergy() const |
---|
| 1127 | { |
---|
| 1128 | return maxKinEnergy; |
---|
| 1129 | } |
---|
| 1130 | |
---|
| 1131 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1132 | |
---|
| 1133 | inline void G4VEnergyLossProcess::SetLambdaFactor(G4double val) |
---|
| 1134 | { |
---|
| 1135 | if(val > 0.0 && val <= 1.0) lambdaFactor = val; |
---|
| 1136 | } |
---|
| 1137 | |
---|
| 1138 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1139 | |
---|
| 1140 | inline void G4VEnergyLossProcess::SetIonisation(G4bool val) |
---|
| 1141 | { |
---|
| 1142 | isIonisation = val; |
---|
| 1143 | if(val) aGPILSelection = CandidateForSelection; |
---|
| 1144 | else aGPILSelection = NotCandidateForSelection; |
---|
| 1145 | } |
---|
| 1146 | |
---|
| 1147 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1148 | |
---|
| 1149 | inline G4bool G4VEnergyLossProcess::IsIonisationProcess() const |
---|
| 1150 | { |
---|
| 1151 | return isIonisation; |
---|
| 1152 | } |
---|
| 1153 | |
---|
| 1154 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1155 | |
---|
| 1156 | void G4VEnergyLossProcess::SetStepFunction(G4double v1, G4double v2) |
---|
| 1157 | { |
---|
| 1158 | dRoverRange = v1; |
---|
| 1159 | finalRange = v2; |
---|
| 1160 | if (dRoverRange > 0.999) dRoverRange = 1.0; |
---|
| 1161 | currentCouple = 0; |
---|
| 1162 | mfpKinEnergy = DBL_MAX; |
---|
| 1163 | } |
---|
| 1164 | |
---|
| 1165 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 1166 | |
---|
| 1167 | #endif |
---|