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