Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 // 27 // ------------------------------------------------------------------- 28 // 29 // GEANT4 Class header file 30 // 31 // 32 // File name: G4VEnergyLossProcess 33 // 34 // Author: Vladimir Ivanchenko on base of Laszlo Urban code 35 // 36 // Creation date: 03.01.2002 37 // 38 // Modifications: Vladimir Ivanchenko 39 // 40 // Class Description: 41 // 42 // It is the unified energy loss process it calculates the continuous 43 // energy loss for charged particles using a set of Energy Loss 44 // models valid for different energy regions. There are a possibility 45 // to create and access to dE/dx and range tables, or to calculate 46 // that information on fly. 47 48 // ------------------------------------------------------------------- 49 // 50 51 #ifndef G4VEnergyLossProcess_h 52 #define G4VEnergyLossProcess_h 1 53 54 #include "G4VContinuousDiscreteProcess.hh" 55 #include "globals.hh" 56 #include "G4Material.hh" 57 #include "G4MaterialCutsCouple.hh" 58 #include "G4Track.hh" 59 #include "G4EmModelManager.hh" 60 #include "G4ParticleChangeForLoss.hh" 61 #include "G4EmTableType.hh" 62 #include "G4EmSecondaryParticleType.hh" 63 #include "G4PhysicsTable.hh" 64 #include "G4PhysicsVector.hh" 65 66 class G4Step; 67 class G4ParticleDefinition; 68 class G4EmParameters; 69 class G4VEmModel; 70 class G4VEmFluctuationModel; 71 class G4DataVector; 72 class G4Region; 73 class G4SafetyHelper; 74 class G4VAtomDeexcitation; 75 class G4VSubCutProducer; 76 class G4EmBiasingManager; 77 class G4LossTableManager; 78 class G4EmDataHandler; 79 80 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 81 82 class G4VEnergyLossProcess : public G4VContinuousDiscreteProcess 83 { 84 public: 85 86 G4VEnergyLossProcess(const G4String& name = "EnergyLoss", 87 G4ProcessType type = fElectromagnetic); 88 89 ~G4VEnergyLossProcess() override; 90 91 //------------------------------------------------------------------------ 92 // Virtual methods to be implemented in concrete processes 93 //------------------------------------------------------------------------ 94 95 protected: 96 97 // description of specific process parameters 98 virtual void StreamProcessInfo(std::ostream&) const {}; 99 100 virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition*, 101 const G4ParticleDefinition*) = 0; 102 103 public: 104 105 // used as low energy limit LambdaTable 106 virtual G4double MinPrimaryEnergy(const G4ParticleDefinition*, 107 const G4Material*, G4double cut); 108 109 // print documentation in html format 110 void ProcessDescription(std::ostream& outFile) const override; 111 112 // prepare all tables 113 void PreparePhysicsTable(const G4ParticleDefinition&) override; 114 115 // build all tables 116 void BuildPhysicsTable(const G4ParticleDefinition&) override; 117 118 // build a table 119 G4PhysicsTable* BuildDEDXTable(G4EmTableType tType = fRestricted); 120 121 // build a table 122 G4PhysicsTable* BuildLambdaTable(G4EmTableType tType = fRestricted); 123 124 // Called before tracking of each new G4Track 125 void StartTracking(G4Track*) override; 126 127 // Step limit from AlongStep 128 G4double AlongStepGetPhysicalInteractionLength( 129 const G4Track&, 130 G4double previousStepSize, 131 G4double currentMinimumStep, 132 G4double& currentSafety, 133 G4GPILSelection* selection) override; 134 135 // Step limit from cross section 136 G4double PostStepGetPhysicalInteractionLength( 137 const G4Track& track, 138 G4double previousStepSize, 139 G4ForceCondition* condition) override; 140 141 // AlongStep computations 142 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&) override; 143 144 // PostStep sampling of secondaries 145 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override; 146 147 // Store all PhysicsTable in files. 148 // Return false in case of any fatal failure at I/O 149 G4bool StorePhysicsTable(const G4ParticleDefinition*, 150 const G4String& directory, 151 G4bool ascii = false) override; 152 153 // Retrieve all Physics from a files. 154 // Return true if all the Physics Table are built. 155 // Return false if any fatal failure. 156 G4bool RetrievePhysicsTable(const G4ParticleDefinition*, 157 const G4String& directory, 158 G4bool ascii) override; 159 160 private: 161 162 // summary printout after initialisation 163 void StreamInfo(std::ostream& out, const G4ParticleDefinition& part, 164 G4bool rst=false) const; 165 166 //------------------------------------------------------------------------ 167 // Public interface to cross section, mfp and sampling of fluctuations 168 // These methods are not used in run time 169 //------------------------------------------------------------------------ 170 171 public: 172 173 // access to dispersion of restricted energy loss 174 G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, 175 const G4DynamicParticle* dp, 176 G4double length); 177 178 // Access to cross section table 179 G4double CrossSectionPerVolume(G4double kineticEnergy, 180 const G4MaterialCutsCouple* couple); 181 G4double CrossSectionPerVolume(G4double kineticEnergy, 182 const G4MaterialCutsCouple* couple, 183 G4double logKineticEnergy); 184 185 // access to cross section 186 G4double MeanFreePath(const G4Track& track); 187 188 // access to step limit 189 G4double ContinuousStepLimit(const G4Track& track, 190 G4double previousStepSize, 191 G4double currentMinimumStep, 192 G4double& currentSafety); 193 194 protected: 195 196 // implementation of the pure virtual method 197 G4double GetMeanFreePath(const G4Track& track, 198 G4double previousStepSize, 199 G4ForceCondition* condition) override; 200 201 // implementation of the pure virtual method 202 G4double GetContinuousStepLimit(const G4Track& track, 203 G4double previousStepSize, 204 G4double currentMinimumStep, 205 G4double& currentSafety) override; 206 207 // creation of an empty vector for cross sections for derived processes 208 G4PhysicsVector* LambdaPhysicsVector(const G4MaterialCutsCouple*, 209 G4double cut); 210 211 inline std::size_t CurrentMaterialCutsCoupleIndex() const; 212 213 //------------------------------------------------------------------------ 214 // Specific methods to set, access, modify models 215 //------------------------------------------------------------------------ 216 217 // Select model in run time 218 inline void SelectModel(G4double kinEnergy); 219 220 public: 221 // Select model by energy and couple index 222 // Not for run time processing 223 inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy, 224 std::size_t& idxCouple) const; 225 226 // Add EM model coupled with fluctuation model for region, smaller value 227 // of order defines which pair of models will be selected for a given 228 // energy interval 229 void AddEmModel(G4int, G4VEmModel*, 230 G4VEmFluctuationModel* fluc = nullptr, 231 const G4Region* region = nullptr); 232 233 // Assign a model to a process local list, to enable the list in run time 234 // the derived process should execute AddEmModel(..) for all such models 235 void SetEmModel(G4VEmModel*, G4int index=0); 236 237 // Access to models 238 inline std::size_t NumberOfModels() const; 239 240 // Return a model from the local list 241 inline G4VEmModel* EmModel(std::size_t index=0) const; 242 243 // Access to models from G4EmModelManager list 244 inline G4VEmModel* GetModelByIndex(std::size_t idx = 0, G4bool ver = false) const; 245 246 // Assign a fluctuation model to a process 247 inline void SetFluctModel(G4VEmFluctuationModel*); 248 249 // Return the assigned fluctuation model 250 inline G4VEmFluctuationModel* FluctModel() const; 251 252 //------------------------------------------------------------------------ 253 // Define and access particle type 254 //------------------------------------------------------------------------ 255 256 protected: 257 inline void SetParticle(const G4ParticleDefinition* p); 258 inline void SetSecondaryParticle(const G4ParticleDefinition* p); 259 260 public: 261 inline void SetBaseParticle(const G4ParticleDefinition* p); 262 inline const G4ParticleDefinition* Particle() const; 263 inline const G4ParticleDefinition* BaseParticle() const; 264 inline const G4ParticleDefinition* SecondaryParticle() const; 265 266 // hide assignment operator 267 G4VEnergyLossProcess(G4VEnergyLossProcess &) = delete; 268 G4VEnergyLossProcess & operator=(const G4VEnergyLossProcess &right) = delete; 269 270 //------------------------------------------------------------------------ 271 // Get/set parameters to configure the process at initialisation time 272 //------------------------------------------------------------------------ 273 274 // Add subcut processor for the region 275 void ActivateSubCutoff(const G4Region* region); 276 277 // Activate biasing 278 void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true); 279 280 void ActivateForcedInteraction(G4double length, 281 const G4String& region, 282 G4bool flag = true); 283 284 void ActivateSecondaryBiasing(const G4String& region, G4double factor, 285 G4double energyLimit); 286 287 inline void SetLossFluctuations(G4bool val); 288 289 inline void SetSpline(G4bool val); 290 inline void SetCrossSectionType(G4CrossSectionType val); 291 inline G4CrossSectionType CrossSectionType() const; 292 293 // Set/Get flag "isIonisation" 294 void SetIonisation(G4bool val); 295 inline G4bool IsIonisationProcess() const; 296 297 // Redefine parameteters for stepping control 298 void SetLinearLossLimit(G4double val); 299 void SetStepFunction(G4double v1, G4double v2); 300 void SetLowestEnergyLimit(G4double); 301 302 inline G4int NumberOfSubCutoffRegions() const; 303 304 //------------------------------------------------------------------------ 305 // Specific methods to path Physics Tables to the process 306 //------------------------------------------------------------------------ 307 308 void SetDEDXTable(G4PhysicsTable* p, G4EmTableType tType); 309 void SetCSDARangeTable(G4PhysicsTable* pRange); 310 void SetRangeTableForLoss(G4PhysicsTable* p); 311 void SetInverseRangeTable(G4PhysicsTable* p); 312 void SetLambdaTable(G4PhysicsTable* p); 313 314 void SetTwoPeaksXS(std::vector<G4TwoPeaksXS*>*); 315 void SetEnergyOfCrossSectionMax(std::vector<G4double>*); 316 317 //------------------------------------------------------------------------ 318 // Specific methods to define custom Physics Tables to the process 319 //------------------------------------------------------------------------ 320 321 // Binning for dEdx, range, inverse range and lambda tables 322 void SetDEDXBinning(G4int nbins); 323 324 // Min kinetic energy for tables 325 void SetMinKinEnergy(G4double e); 326 inline G4double MinKinEnergy() const; 327 328 // Max kinetic energy for tables 329 void SetMaxKinEnergy(G4double e); 330 inline G4double MaxKinEnergy() const; 331 332 // Biasing parameters 333 inline G4double CrossSectionBiasingFactor() const; 334 335 // Return values for given G4MaterialCutsCouple 336 inline G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple*); 337 inline G4double GetCSDADEDX(G4double kineticEnergy, 338 const G4MaterialCutsCouple*); 339 inline G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple*, 340 G4double logKineticEnergy); 341 inline G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple*); 342 inline G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple*, 343 G4double logKineticEnergy); 344 inline G4double GetCSDARange(G4double kineticEnergy, 345 const G4MaterialCutsCouple*); 346 inline G4double GetKineticEnergy(G4double range, 347 const G4MaterialCutsCouple*); 348 inline G4double GetLambda(G4double kineticEnergy,const G4MaterialCutsCouple*); 349 inline G4double GetLambda(G4double kineticEnergy,const G4MaterialCutsCouple*, 350 G4double logKineticEnergy); 351 352 inline G4bool TablesAreBuilt() const; 353 354 // Access to specific tables 355 inline G4PhysicsTable* DEDXTable() const; 356 inline G4PhysicsTable* DEDXunRestrictedTable() const; 357 inline G4PhysicsTable* IonisationTable() const; 358 inline G4PhysicsTable* CSDARangeTable() const; 359 inline G4PhysicsTable* RangeTableForLoss() const; 360 inline G4PhysicsTable* InverseRangeTable() const; 361 inline G4PhysicsTable* LambdaTable() const; 362 inline std::vector<G4TwoPeaksXS*>* TwoPeaksXS() const; 363 inline std::vector<G4double>* EnergyOfCrossSectionMax() const; 364 365 inline G4bool UseBaseMaterial() const; 366 367 //------------------------------------------------------------------------ 368 // Run time method for simulation of ionisation 369 //------------------------------------------------------------------------ 370 371 // access atom on which interaction happens 372 const G4Element* GetCurrentElement() const; 373 374 // Set scaling parameters for ions is needed to G4EmCalculator 375 void SetDynamicMassCharge(G4double massratio, G4double charge2ratio); 376 377 private: 378 379 void FillSecondariesAlongStep(G4double weight); 380 381 void PrintWarning(const G4String&, G4double val) const; 382 383 // define material and indexes 384 inline void DefineMaterial(const G4MaterialCutsCouple* couple); 385 386 //------------------------------------------------------------------------ 387 // Compute values using scaling relation, mass and charge of based particle 388 //------------------------------------------------------------------------ 389 inline G4double GetDEDXForScaledEnergy(G4double scaledKinE); 390 inline G4double GetDEDXForScaledEnergy(G4double scaledKinE, 391 G4double logScaledKinE); 392 inline G4double GetIonisationForScaledEnergy(G4double scaledKinE); 393 inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinE); 394 inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinE, 395 G4double logScaledKinE); 396 397 inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinE); 398 inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinE, 399 G4double logScaledKinE); 400 401 inline G4double ScaledKinEnergyForLoss(G4double range); 402 inline G4double GetLambdaForScaledEnergy(G4double scaledKinE); 403 inline G4double GetLambdaForScaledEnergy(G4double scaledKinE, 404 G4double logScaledKinE); 405 406 inline G4double LogScaledEkin(const G4Track& aTrack); 407 408 void ComputeLambdaForScaledEnergy(G4double scaledKinE, 409 const G4Track& aTrack); 410 411 G4bool IsRegionForCubcutProcessor(const G4Track& aTrack); 412 413 protected: 414 415 G4ParticleChangeForLoss fParticleChange; 416 const G4Material* currentMaterial = nullptr; 417 const G4MaterialCutsCouple* currentCouple = nullptr; 418 419 private: 420 421 G4LossTableManager* lManager; 422 G4EmModelManager* modelManager; 423 G4VEmModel* currentModel = nullptr; 424 G4EmBiasingManager* biasManager = nullptr; 425 G4SafetyHelper* safetyHelper; 426 G4EmParameters* theParameters; 427 G4VEmFluctuationModel* fluctModel = nullptr; 428 G4VAtomDeexcitation* atomDeexcitation = nullptr; 429 G4VSubCutProducer* subcutProducer = nullptr; 430 431 const G4ParticleDefinition* particle = nullptr; 432 const G4ParticleDefinition* baseParticle = nullptr; 433 const G4ParticleDefinition* secondaryParticle = nullptr; 434 G4EmDataHandler* theData = nullptr; 435 436 G4PhysicsTable* theDEDXTable = nullptr; 437 G4PhysicsTable* theDEDXunRestrictedTable = nullptr; 438 G4PhysicsTable* theIonisationTable = nullptr; 439 G4PhysicsTable* theRangeTableForLoss = nullptr; 440 G4PhysicsTable* theCSDARangeTable = nullptr; 441 G4PhysicsTable* theInverseRangeTable = nullptr; 442 G4PhysicsTable* theLambdaTable = nullptr; 443 444 std::vector<const G4Region*>* scoffRegions = nullptr; 445 std::vector<G4VEmModel*>* emModels = nullptr; 446 const std::vector<G4int>* theDensityIdx = nullptr; 447 const std::vector<G4double>* theDensityFactor = nullptr; 448 const G4DataVector* theCuts = nullptr; 449 450 std::vector<G4double>* theEnergyOfCrossSectionMax = nullptr; 451 std::vector<G4TwoPeaksXS*>* fXSpeaks = nullptr; 452 453 G4double lowestKinEnergy; 454 G4double minKinEnergy; 455 G4double maxKinEnergy; 456 G4double maxKinEnergyCSDA; 457 458 G4double linLossLimit = 0.01; 459 G4double dRoverRange = 0.2; 460 G4double finalRange; 461 G4double lambdaFactor = 0.8; 462 G4double invLambdaFactor; 463 G4double biasFactor = 1.0; 464 465 G4double massRatio = 1.0; 466 G4double logMassRatio = 0.0; 467 G4double fFactor = 1.0; 468 G4double reduceFactor = 1.0; 469 G4double chargeSqRatio = 1.0; 470 G4double fRange = 0.0; 471 G4double fRangeEnergy = 0.0; 472 473 protected: 474 475 G4double preStepLambda = 0.0; 476 G4double preStepKinEnergy = 0.0; 477 G4double preStepScaledEnergy = 0.0; 478 G4double mfpKinEnergy = 0.0; 479 480 std::size_t currentCoupleIndex = 0; 481 482 private: 483 484 G4int nBins; 485 G4int nBinsCSDA; 486 G4int numberOfModels = 0; 487 G4int nSCoffRegions = 0; 488 G4int secID = _DeltaElectron; 489 G4int tripletID = _TripletElectron; 490 G4int biasID = _DeltaEBelowCut; 491 G4int epixeID = _ePIXE; 492 G4int gpixeID = _GammaPIXE; 493 G4int mainSecondaries = 1; 494 495 std::size_t basedCoupleIndex = 0; 496 std::size_t coupleIdxRange = 0; 497 std::size_t idxDEDX = 0; 498 std::size_t idxDEDXunRestricted = 0; 499 std::size_t idxIonisation = 0; 500 std::size_t idxRange = 0; 501 std::size_t idxCSDA = 0; 502 std::size_t idxSecRange = 0; 503 std::size_t idxInverseRange = 0; 504 std::size_t idxLambda = 0; 505 506 G4GPILSelection aGPILSelection; 507 G4CrossSectionType fXSType = fEmOnePeak; 508 509 G4bool lossFluctuationFlag = true; 510 G4bool useCutAsFinalRange = false; 511 G4bool tablesAreBuilt = false; 512 G4bool spline = true; 513 G4bool isIon = false; 514 G4bool isIonisation = false; 515 G4bool useDeexcitation = false; 516 G4bool biasFlag = false; 517 G4bool weightFlag = false; 518 G4bool isMaster = false; 519 G4bool baseMat = false; 520 G4bool actLinLossLimit = false; 521 G4bool actLossFluc = false; 522 G4bool actBinning = false; 523 G4bool actMinKinEnergy = false; 524 G4bool actMaxKinEnergy = false; 525 526 std::vector<G4DynamicParticle*> secParticles; 527 std::vector<G4Track*> scTracks; 528 }; 529 530 // ======== Run time inline methods ================ 531 532 inline std::size_t G4VEnergyLossProcess::CurrentMaterialCutsCoupleIndex() const 533 { 534 return currentCoupleIndex; 535 } 536 537 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 538 539 inline void G4VEnergyLossProcess::SelectModel(G4double kinEnergy) 540 { 541 currentModel = modelManager->SelectModel(kinEnergy, currentCoupleIndex); 542 currentModel->SetCurrentCouple(currentCouple); 543 } 544 545 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 546 547 inline G4VEmModel* G4VEnergyLossProcess::SelectModelForMaterial( 548 G4double kinEnergy, std::size_t& idx) const 549 { 550 return modelManager->SelectModel(kinEnergy, idx); 551 } 552 553 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 554 555 inline void 556 G4VEnergyLossProcess::DefineMaterial(const G4MaterialCutsCouple* couple) 557 { 558 if(couple != currentCouple) { 559 currentCouple = couple; 560 currentMaterial = couple->GetMaterial(); 561 basedCoupleIndex = currentCoupleIndex = couple->GetIndex(); 562 fFactor = chargeSqRatio*biasFactor; 563 mfpKinEnergy = DBL_MAX; 564 idxLambda = 0; 565 if(baseMat) { 566 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex]; 567 fFactor *= (*theDensityFactor)[currentCoupleIndex]; 568 } 569 reduceFactor = 1.0/(fFactor*massRatio); 570 } 571 } 572 573 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 574 575 inline G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e) 576 { 577 /* 578 G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= " 579 << basedCoupleIndex << " E(MeV)= " << e 580 << " Emin= " << minKinEnergy << " Factor= " << fFactor 581 << " " << theDEDXTable << G4endl; */ 582 G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->Value(e, idxDEDX); 583 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); } 584 return x; 585 } 586 587 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 588 589 inline 590 G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e, G4double loge) 591 { 592 /* 593 G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= " 594 << basedCoupleIndex << " E(MeV)= " << e 595 << " Emin= " << minKinEnergy << " Factor= " << fFactor 596 << " " << theDEDXTable << G4endl; */ 597 G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->LogVectorValue(e,loge); 598 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); } 599 return x; 600 } 601 602 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 603 604 inline G4double G4VEnergyLossProcess::GetIonisationForScaledEnergy(G4double e) 605 { 606 G4double x = 607 fFactor*(*theIonisationTable)[basedCoupleIndex]->Value(e, idxIonisation); 608 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); } 609 return x; 610 } 611 612 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 613 614 inline G4double G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e) 615 { 616 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= " 617 // << basedCoupleIndex << " E(MeV)= " << e 618 // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl; 619 if(currentCoupleIndex != coupleIdxRange || fRangeEnergy != e) { 620 coupleIdxRange = currentCoupleIndex; 621 fRangeEnergy = e; 622 fRange = reduceFactor*((*theRangeTableForLoss)[basedCoupleIndex])->Value(e, idxRange); 623 if (fRange < 0.0) { fRange = 0.0; } 624 else if (e < minKinEnergy) { fRange *= std::sqrt(e/minKinEnergy); } 625 } 626 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= " 627 // << basedCoupleIndex << " E(MeV)= " << e 628 // << " R= " << computedRange << " " << theRangeTableForLoss << G4endl; 629 return fRange; 630 } 631 632 inline G4double 633 G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e, G4double loge) 634 { 635 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= " 636 // << basedCoupleIndex << " E(MeV)= " << e 637 // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl; 638 if(currentCoupleIndex != coupleIdxRange || fRangeEnergy != e) { 639 coupleIdxRange = currentCoupleIndex; 640 fRangeEnergy = e; 641 fRange = reduceFactor*((*theRangeTableForLoss)[basedCoupleIndex])->LogVectorValue(e, loge); 642 if (fRange < 0.0) { fRange = 0.0; } 643 else if (e < minKinEnergy) { fRange *= std::sqrt(e/minKinEnergy); } 644 } 645 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= " 646 // << basedCoupleIndex << " E(MeV)= " << e 647 // << " R= " << fRange << " " << theRangeTableForLoss << G4endl; 648 return fRange; 649 } 650 651 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 652 653 inline G4double 654 G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e) 655 { 656 G4double x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e, idxCSDA); 657 if (x < 0.0) { x = 0.0; } 658 else if (e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); } 659 return x; 660 } 661 662 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 663 664 inline G4double 665 G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e, 666 G4double loge) 667 { 668 G4double x = ((*theCSDARangeTable)[basedCoupleIndex])->LogVectorValue(e, loge); 669 if (x < 0.0) { x = 0.0; } 670 else if (e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); } 671 return x; 672 } 673 674 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 675 676 inline G4double G4VEnergyLossProcess::ScaledKinEnergyForLoss(G4double r) 677 { 678 //G4cout << "G4VEnergyLossProcess::GetEnergy: Idx= " 679 // << basedCoupleIndex << " R(mm)= " << r << " " 680 // << theInverseRangeTable << G4endl; 681 G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex]; 682 G4double rmin = v->Energy(0); 683 G4double e = 0.0; 684 if(r >= rmin) { e = v->Value(r, idxInverseRange); } 685 else if(r > 0.0) { 686 G4double x = r/rmin; 687 e = minKinEnergy*x*x; 688 } 689 return e; 690 } 691 692 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 693 694 inline G4double G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e) 695 { 696 return fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda); 697 } 698 699 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 700 701 inline G4double 702 G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e, G4double loge) 703 { 704 return fFactor*((*theLambdaTable)[basedCoupleIndex])->LogVectorValue(e, loge); 705 } 706 707 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 708 709 inline G4double G4VEnergyLossProcess::LogScaledEkin(const G4Track& track) 710 { 711 return track.GetDynamicParticle()->GetLogKineticEnergy() + logMassRatio; 712 } 713 714 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 715 716 inline G4double 717 G4VEnergyLossProcess::GetDEDX(G4double kinEnergy, 718 const G4MaterialCutsCouple* couple) 719 { 720 DefineMaterial(couple); 721 return GetDEDXForScaledEnergy(kinEnergy*massRatio); 722 } 723 724 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 725 726 inline G4double 727 G4VEnergyLossProcess::GetDEDX(G4double kinEnergy, 728 const G4MaterialCutsCouple* couple, 729 G4double logKinEnergy) 730 { 731 DefineMaterial(couple); 732 return GetDEDXForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio); 733 } 734 735 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 736 737 inline G4double 738 G4VEnergyLossProcess::GetRange(G4double kinEnergy, 739 const G4MaterialCutsCouple* couple) 740 { 741 DefineMaterial(couple); 742 return GetScaledRangeForScaledEnergy(kinEnergy*massRatio); 743 } 744 745 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 746 747 inline G4double 748 G4VEnergyLossProcess::GetRange(G4double kinEnergy, 749 const G4MaterialCutsCouple* couple, 750 G4double logKinEnergy) 751 { 752 DefineMaterial(couple); 753 return GetScaledRangeForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio); 754 } 755 756 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 757 758 inline G4double 759 G4VEnergyLossProcess::GetCSDARange(G4double kineticEnergy, 760 const G4MaterialCutsCouple* couple) 761 { 762 DefineMaterial(couple); 763 return (nullptr == theCSDARangeTable) ? DBL_MAX : 764 GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor; 765 } 766 767 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 768 769 inline G4double 770 G4VEnergyLossProcess::GetKineticEnergy(G4double range, 771 const G4MaterialCutsCouple* couple) 772 { 773 DefineMaterial(couple); 774 return ScaledKinEnergyForLoss(range/reduceFactor)/massRatio; 775 } 776 777 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 778 779 inline G4double 780 G4VEnergyLossProcess::GetLambda(G4double kinEnergy, 781 const G4MaterialCutsCouple* couple) 782 { 783 DefineMaterial(couple); 784 return (nullptr != theLambdaTable) ? 785 GetLambdaForScaledEnergy(kinEnergy*massRatio) : 0.0; 786 } 787 788 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 789 790 inline G4double 791 G4VEnergyLossProcess::GetLambda(G4double kinEnergy, 792 const G4MaterialCutsCouple* couple, 793 G4double logKinEnergy) 794 { 795 DefineMaterial(couple); 796 return (nullptr != theLambdaTable) ? 797 GetLambdaForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio) 798 : 0.0; 799 } 800 801 // ======== Get/Set inline methods used at initialisation ================ 802 803 inline void G4VEnergyLossProcess::SetFluctModel(G4VEmFluctuationModel* p) 804 { 805 fluctModel = p; 806 } 807 808 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 809 810 inline G4VEmFluctuationModel* G4VEnergyLossProcess::FluctModel() const 811 { 812 return fluctModel; 813 } 814 815 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 816 817 inline void G4VEnergyLossProcess::SetParticle(const G4ParticleDefinition* p) 818 { 819 particle = p; 820 } 821 822 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 823 824 inline void 825 G4VEnergyLossProcess::SetSecondaryParticle(const G4ParticleDefinition* p) 826 { 827 secondaryParticle = p; 828 } 829 830 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 831 832 inline void 833 G4VEnergyLossProcess::SetBaseParticle(const G4ParticleDefinition* p) 834 { 835 baseParticle = p; 836 } 837 838 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 839 840 inline const G4ParticleDefinition* G4VEnergyLossProcess::Particle() const 841 { 842 return particle; 843 } 844 845 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 846 847 inline const G4ParticleDefinition* G4VEnergyLossProcess::BaseParticle() const 848 { 849 return baseParticle; 850 } 851 852 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 853 854 inline const G4ParticleDefinition* 855 G4VEnergyLossProcess::SecondaryParticle() const 856 { 857 return secondaryParticle; 858 } 859 860 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 861 862 inline void G4VEnergyLossProcess::SetLossFluctuations(G4bool val) 863 { 864 lossFluctuationFlag = val; 865 actLossFluc = true; 866 } 867 868 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 869 870 inline void G4VEnergyLossProcess::SetSpline(G4bool val) 871 { 872 spline = val; 873 } 874 875 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 876 877 inline void G4VEnergyLossProcess::SetCrossSectionType(G4CrossSectionType val) 878 { 879 fXSType = val; 880 } 881 882 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 883 884 inline G4CrossSectionType G4VEnergyLossProcess::CrossSectionType() const 885 { 886 return fXSType; 887 } 888 889 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 890 891 inline G4bool G4VEnergyLossProcess::IsIonisationProcess() const 892 { 893 return isIonisation; 894 } 895 896 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 897 898 inline G4int G4VEnergyLossProcess::NumberOfSubCutoffRegions() const 899 { 900 return nSCoffRegions; 901 } 902 903 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 904 905 inline G4double G4VEnergyLossProcess::MinKinEnergy() const 906 { 907 return minKinEnergy; 908 } 909 910 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 911 912 inline G4double G4VEnergyLossProcess::MaxKinEnergy() const 913 { 914 return maxKinEnergy; 915 } 916 917 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 918 919 inline G4double G4VEnergyLossProcess::CrossSectionBiasingFactor() const 920 { 921 return biasFactor; 922 } 923 924 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 925 926 inline G4bool G4VEnergyLossProcess::TablesAreBuilt() const 927 { 928 return tablesAreBuilt; 929 } 930 931 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 932 933 inline G4PhysicsTable* G4VEnergyLossProcess::DEDXTable() const 934 { 935 return theDEDXTable; 936 } 937 938 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 939 940 inline G4PhysicsTable* G4VEnergyLossProcess::DEDXunRestrictedTable() const 941 { 942 return theDEDXunRestrictedTable; 943 } 944 945 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 946 947 inline G4PhysicsTable* G4VEnergyLossProcess::IonisationTable() const 948 { 949 return theIonisationTable; 950 } 951 952 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 953 954 inline G4PhysicsTable* G4VEnergyLossProcess::CSDARangeTable() const 955 { 956 return theCSDARangeTable; 957 } 958 959 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 960 961 inline G4PhysicsTable* G4VEnergyLossProcess::RangeTableForLoss() const 962 { 963 return theRangeTableForLoss; 964 } 965 966 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 967 968 inline G4PhysicsTable* G4VEnergyLossProcess::InverseRangeTable() const 969 { 970 return theInverseRangeTable; 971 } 972 973 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 974 975 inline G4PhysicsTable* G4VEnergyLossProcess::LambdaTable() const 976 { 977 return theLambdaTable; 978 } 979 980 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 981 982 inline G4bool G4VEnergyLossProcess::UseBaseMaterial() const 983 { 984 return baseMat; 985 } 986 987 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 988 989 inline std::vector<G4double>* 990 G4VEnergyLossProcess::EnergyOfCrossSectionMax() const 991 { 992 return theEnergyOfCrossSectionMax; 993 } 994 995 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 996 997 inline std::vector<G4TwoPeaksXS*>* G4VEnergyLossProcess::TwoPeaksXS() const 998 { 999 return fXSpeaks; 1000 } 1001 1002 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 1003 1004 inline std::size_t G4VEnergyLossProcess::NumberOfModels() const 1005 { 1006 return numberOfModels; 1007 } 1008 1009 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 1010 1011 inline G4VEmModel* G4VEnergyLossProcess::EmModel(std::size_t index) const 1012 { 1013 return (index < emModels->size()) ? (*emModels)[index] : nullptr; 1014 } 1015 1016 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 1017 1018 inline G4VEmModel* 1019 G4VEnergyLossProcess::GetModelByIndex(std::size_t idx, G4bool ver) const 1020 { 1021 return modelManager->GetModel((G4int)idx, ver); 1022 } 1023 1024 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 1025 1026 #endif 1027