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 // GEANT4 Class header file 29 // 30 // 31 // File name: G4VEmProcess 32 // 33 // Author: Vladimir Ivanchenko 34 // 35 // Creation date: 01.10.2003 36 // 37 // Modifications: Vladimir Ivanchenko 38 // 39 // Class Description: 40 // 41 // It is the base class - EM discrete and rest/discrete process 42 43 // ------------------------------------------------------------------- 44 // 45 46 #ifndef G4VEmProcess_h 47 #define G4VEmProcess_h 1 48 49 #include <CLHEP/Units/SystemOfUnits.h> 50 51 #include "G4VDiscreteProcess.hh" 52 #include "globals.hh" 53 #include "G4Material.hh" 54 #include "G4MaterialCutsCouple.hh" 55 #include "G4Track.hh" 56 #include "G4UnitsTable.hh" 57 #include "G4ParticleDefinition.hh" 58 #include "G4ParticleChangeForGamma.hh" 59 #include "G4EmParameters.hh" 60 #include "G4EmDataHandler.hh" 61 #include "G4EmTableType.hh" 62 #include "G4EmModelManager.hh" 63 #include "G4EmSecondaryParticleType.hh" 64 65 class G4Step; 66 class G4VEmModel; 67 class G4DataVector; 68 class G4VParticleChange; 69 class G4PhysicsTable; 70 class G4PhysicsVector; 71 class G4EmBiasingManager; 72 class G4LossTableManager; 73 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 75 76 class G4VEmProcess : public G4VDiscreteProcess 77 { 78 public: 79 80 G4VEmProcess(const G4String& name, G4ProcessType type = fElectromagnetic); 81 82 ~G4VEmProcess() override; 83 84 //------------------------------------------------------------------------ 85 // Virtual methods to be implemented in concrete processes 86 //------------------------------------------------------------------------ 87 88 void ProcessDescription(std::ostream& outFile) const override; 89 90 protected: 91 92 virtual void StreamProcessInfo(std::ostream&) const {}; 93 94 virtual void InitialiseProcess(const G4ParticleDefinition*) = 0; 95 96 //------------------------------------------------------------------------ 97 // Implementation of virtual methods common to all Discrete processes 98 //------------------------------------------------------------------------ 99 100 public: 101 102 // Initialise for build of tables 103 void PreparePhysicsTable(const G4ParticleDefinition&) override; 104 105 // Build physics table during initialisation 106 void BuildPhysicsTable(const G4ParticleDefinition&) override; 107 108 // Called before tracking of each new G4Track 109 void StartTracking(G4Track*) override; 110 111 // implementation of virtual method, specific for G4VEmProcess 112 G4double PostStepGetPhysicalInteractionLength( 113 const G4Track& track, 114 G4double previousStepSize, 115 G4ForceCondition* condition) override; 116 117 // implementation of virtual method, specific for G4VEmProcess 118 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override; 119 120 // Store PhysicsTable in a file. 121 // Return false in case of failure at I/O 122 G4bool StorePhysicsTable(const G4ParticleDefinition*, 123 const G4String& directory, 124 G4bool ascii = false) override; 125 126 // Retrieve Physics from a file. 127 // (return true if the Physics Table can be build by using file) 128 // (return false if the process has no functionality or in case of failure) 129 // File name should is constructed as processName+particleName and the 130 // should be placed under the directory specified by the argument. 131 G4bool RetrievePhysicsTable(const G4ParticleDefinition*, 132 const G4String& directory, 133 G4bool ascii) override; 134 135 // allowing check process name 136 virtual G4VEmProcess* GetEmProcess(const G4String& name); 137 138 //------------------------------------------------------------------------ 139 // Specific methods for Discrete EM post step simulation 140 //------------------------------------------------------------------------ 141 142 // The main method to access cross section per volume 143 inline G4double GetLambda(G4double kinEnergy, 144 const G4MaterialCutsCouple* couple, 145 G4double logKinEnergy); 146 147 // It returns the cross section per volume for energy/material 148 G4double GetCrossSection(const G4double kinEnergy, 149 const G4MaterialCutsCouple* couple) override; 150 151 // It returns the cross section of the process per atom 152 G4double ComputeCrossSectionPerAtom(G4double kineticEnergy, 153 G4double Z, G4double A=0., 154 G4double cut=0.0); 155 156 inline G4double MeanFreePath(const G4Track& track); 157 158 //------------------------------------------------------------------------ 159 // Specific methods to build and access Physics Tables 160 //------------------------------------------------------------------------ 161 162 // Binning for lambda table 163 void SetLambdaBinning(G4int nbins); 164 165 // Min kinetic energy for tables 166 void SetMinKinEnergy(G4double e); 167 168 // Min kinetic energy for high energy table 169 void SetMinKinEnergyPrim(G4double e); 170 171 // Max kinetic energy for tables 172 void SetMaxKinEnergy(G4double e); 173 174 // Cross section table pointers 175 inline G4PhysicsTable* LambdaTable() const; 176 inline G4PhysicsTable* LambdaTablePrim() const; 177 inline void SetLambdaTable(G4PhysicsTable*); 178 inline void SetLambdaTablePrim(G4PhysicsTable*); 179 180 // Integral method type and peak positions 181 inline std::vector<G4double>* EnergyOfCrossSectionMax() const; 182 inline void SetEnergyOfCrossSectionMax(std::vector<G4double>*); 183 inline G4CrossSectionType CrossSectionType() const; 184 inline void SetCrossSectionType(G4CrossSectionType val); 185 186 //------------------------------------------------------------------------ 187 // Define and access particle type 188 //------------------------------------------------------------------------ 189 190 inline const G4ParticleDefinition* Particle() const; 191 inline const G4ParticleDefinition* SecondaryParticle() const; 192 193 protected: 194 195 //------------------------------------------------------------------------ 196 // Specific methods to set, access, modify models and basic parameters 197 //------------------------------------------------------------------------ 198 199 // Select model in run time 200 inline G4VEmModel* SelectModel(G4double kinEnergy, std::size_t); 201 202 public: 203 204 // Select model by energy and couple index 205 inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy, 206 std::size_t idxCouple) const; 207 208 // Add model for region, smaller value of order defines which 209 // model will be selected for a given energy interval 210 void AddEmModel(G4int, G4VEmModel*, const G4Region* region = nullptr); 211 212 // Assign a model to a process local list, to enable the list in run time 213 // the derived process should execute AddEmModel(..) for all such models 214 void SetEmModel(G4VEmModel*, G4int index = 0); 215 216 inline G4int NumberOfModels() const; 217 218 // return a model from the local list 219 inline G4VEmModel* EmModel(std::size_t index = 0) const; 220 221 // Access to active model 222 inline const G4VEmModel* GetCurrentModel() const; 223 224 // Access to models 225 inline G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const; 226 227 // Access to the current G4Element 228 const G4Element* GetCurrentElement() const; 229 230 // Biasing parameters 231 void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true); 232 inline G4double CrossSectionBiasingFactor() const; 233 234 // Activate forced interaction 235 void ActivateForcedInteraction(G4double length = 0.0, 236 const G4String& r = "", 237 G4bool flag = true); 238 239 void ActivateSecondaryBiasing(const G4String& region, G4double factor, 240 G4double energyLimit); 241 242 inline void SetEmMasterProcess(const G4VEmProcess*); 243 244 inline void SetBuildTableFlag(G4bool val); 245 246 inline void CurrentSetup(const G4MaterialCutsCouple*, G4double energy); 247 248 inline G4bool UseBaseMaterial() const; 249 250 void BuildLambdaTable(); 251 252 void StreamInfo(std::ostream& outFile, const G4ParticleDefinition&, 253 G4bool rst=false) const; 254 255 // hide copy constructor and assignment operator 256 G4VEmProcess(G4VEmProcess &) = delete; 257 G4VEmProcess & operator=(const G4VEmProcess &right) = delete; 258 259 //------------------------------------------------------------------------ 260 // Other generic methods 261 //------------------------------------------------------------------------ 262 263 protected: 264 265 G4double GetMeanFreePath(const G4Track& track, 266 G4double previousStepSize, 267 G4ForceCondition* condition) override; 268 269 G4PhysicsVector* LambdaPhysicsVector(const G4MaterialCutsCouple*); 270 271 inline void DefineMaterial(const G4MaterialCutsCouple* couple); 272 273 inline G4int LambdaBinning() const; 274 275 inline G4double MinKinEnergy() const; 276 277 inline G4double MaxKinEnergy() const; 278 279 // Single scattering parameters 280 inline G4double PolarAngleLimit() const; 281 282 inline G4ParticleChangeForGamma* GetParticleChange(); 283 284 inline void SetParticle(const G4ParticleDefinition* p); 285 286 inline void SetSecondaryParticle(const G4ParticleDefinition* p); 287 288 inline std::size_t CurrentMaterialCutsCoupleIndex() const; 289 290 inline const G4MaterialCutsCouple* MaterialCutsCouple() const; 291 292 inline G4bool ApplyCuts() const; 293 294 inline G4double GetGammaEnergyCut(); 295 296 inline G4double GetElectronEnergyCut(); 297 298 inline void SetStartFromNullFlag(G4bool val); 299 300 inline void SetSplineFlag(G4bool val); 301 302 const G4Element* GetTargetElement() const; 303 304 const G4Isotope* GetTargetIsotope() const; 305 306 // these two methods assume that vectors are initilized 307 // and idx is within vector length 308 inline G4int DensityIndex(G4int idx) const; 309 inline G4double DensityFactor(G4int idx) const; 310 311 private: 312 313 void PrintWarning(G4String tit, G4double val); 314 315 void ComputeIntegralLambda(G4double kinEnergy, const G4Track&); 316 317 inline G4double LogEkin(const G4Track&); 318 319 inline G4double GetLambdaFromTable(G4double kinEnergy); 320 321 inline G4double GetLambdaFromTable(G4double kinEnergy, G4double logKinEnergy); 322 323 inline G4double GetLambdaFromTablePrim(G4double kinEnergy); 324 325 inline G4double GetLambdaFromTablePrim(G4double kinEnergy, G4double logKinEnergy); 326 327 inline G4double GetCurrentLambda(G4double kinEnergy); 328 329 inline G4double GetCurrentLambda(G4double kinEnergy, G4double logKinEnergy); 330 331 inline G4double ComputeCurrentLambda(G4double kinEnergy); 332 333 // ======== pointers ========= 334 335 G4EmModelManager* modelManager = nullptr; 336 const G4ParticleDefinition* particle = nullptr; 337 const G4ParticleDefinition* currentParticle = nullptr; 338 const G4ParticleDefinition* theGamma = nullptr; 339 const G4ParticleDefinition* theElectron = nullptr; 340 const G4ParticleDefinition* thePositron = nullptr; 341 const G4ParticleDefinition* secondaryParticle = nullptr; 342 const G4VEmProcess* masterProc = nullptr; 343 G4EmDataHandler* theData = nullptr; 344 G4VEmModel* currentModel = nullptr; 345 G4LossTableManager* lManager = nullptr; 346 G4EmParameters* theParameters = nullptr; 347 const G4Material* baseMaterial = nullptr; 348 349 // ======== tables and vectors ======== 350 G4PhysicsTable* theLambdaTable = nullptr; 351 G4PhysicsTable* theLambdaTablePrim = nullptr; 352 353 const std::vector<G4double>* theCuts = nullptr; 354 const std::vector<G4double>* theCutsGamma = nullptr; 355 const std::vector<G4double>* theCutsElectron = nullptr; 356 const std::vector<G4double>* theCutsPositron = nullptr; 357 358 protected: 359 360 // ======== pointers ========= 361 362 const G4MaterialCutsCouple* currentCouple = nullptr; 363 const G4Material* currentMaterial = nullptr; 364 G4EmBiasingManager* biasManager = nullptr; 365 std::vector<G4double>* theEnergyOfCrossSectionMax = nullptr; 366 367 private: 368 369 const std::vector<G4double>* theDensityFactor = nullptr; 370 const std::vector<G4int>* theDensityIdx = nullptr; 371 372 // ======== parameters ========= 373 G4double minKinEnergy; 374 G4double maxKinEnergy; 375 G4double minKinEnergyPrim = DBL_MAX; 376 G4double lambdaFactor = 0.8; 377 G4double invLambdaFactor; 378 G4double biasFactor = 1.0; 379 G4double massRatio = 1.0; 380 G4double fFactor = 1.0; 381 G4double fLambda = 0.0; 382 G4double fLambdaEnergy = 0.0; 383 384 protected: 385 386 G4double mfpKinEnergy = DBL_MAX; 387 G4double preStepKinEnergy = 0.0; 388 G4double preStepLambda = 0.0; 389 390 private: 391 392 G4CrossSectionType fXSType = fEmNoIntegral; 393 394 G4int numberOfModels = 0; 395 G4int nLambdaBins = 84; 396 397 protected: 398 399 G4int mainSecondaries = 1; 400 G4int secID = _EM; 401 G4int fluoID = _Fluorescence; 402 G4int augerID = _AugerElectron; 403 G4int biasID = _EM; 404 G4int tripletID = _TripletElectron; 405 std::size_t currentCoupleIndex = 0; 406 std::size_t basedCoupleIndex = 0; 407 std::size_t coupleIdxLambda = 0; 408 std::size_t idxLambda = 0; 409 410 G4bool isTheMaster = false; 411 G4bool baseMat = false; 412 413 private: 414 415 G4bool buildLambdaTable = true; 416 G4bool applyCuts = false; 417 G4bool startFromNull = false; 418 G4bool splineFlag = true; 419 G4bool actMinKinEnergy = false; 420 G4bool actMaxKinEnergy = false; 421 G4bool actBinning = false; 422 G4bool isIon = false; 423 G4bool biasFlag = false; 424 G4bool weightFlag = false; 425 426 protected: 427 428 // ======== particle change ========= 429 std::vector<G4DynamicParticle*> secParticles; 430 G4ParticleChangeForGamma fParticleChange; 431 432 private: 433 434 // ======== local vectors ========= 435 std::vector<G4VEmModel*> emModels; 436 437 }; 438 439 // ======== Run time inline methods ================ 440 441 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 442 443 inline std::size_t G4VEmProcess::CurrentMaterialCutsCoupleIndex() const 444 { 445 return currentCoupleIndex; 446 } 447 448 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 449 450 inline const G4MaterialCutsCouple* G4VEmProcess::MaterialCutsCouple() const 451 { 452 return currentCouple; 453 } 454 455 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 456 457 inline G4double G4VEmProcess::GetGammaEnergyCut() 458 { 459 return (*theCutsGamma)[currentCoupleIndex]; 460 } 461 462 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 463 464 inline G4double G4VEmProcess::GetElectronEnergyCut() 465 { 466 return (*theCutsElectron)[currentCoupleIndex]; 467 } 468 469 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 470 471 inline void G4VEmProcess::DefineMaterial(const G4MaterialCutsCouple* couple) 472 { 473 if (couple != currentCouple) { 474 currentCouple = couple; 475 baseMaterial = currentMaterial = couple->GetMaterial(); 476 basedCoupleIndex = currentCoupleIndex = couple->GetIndex(); 477 fFactor = biasFactor; 478 mfpKinEnergy = DBL_MAX; 479 if (baseMat) { 480 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex]; 481 if (nullptr != currentMaterial->GetBaseMaterial()) 482 baseMaterial = currentMaterial->GetBaseMaterial(); 483 fFactor *= (*theDensityFactor)[currentCoupleIndex]; 484 } 485 } 486 } 487 488 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 489 490 inline 491 G4VEmModel* G4VEmProcess::SelectModel(G4double kinEnergy, std::size_t) 492 { 493 if(1 < numberOfModels) { 494 currentModel = modelManager->SelectModel(kinEnergy, currentCoupleIndex); 495 } 496 currentModel->SetCurrentCouple(currentCouple); 497 return currentModel; 498 } 499 500 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 501 502 inline 503 G4VEmModel* G4VEmProcess::SelectModelForMaterial(G4double kinEnergy, 504 std::size_t idxCouple) const 505 { 506 return modelManager->SelectModel(kinEnergy, idxCouple); 507 } 508 509 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 510 511 inline G4double G4VEmProcess::GetLambdaFromTable(G4double e) 512 { 513 return ((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda); 514 } 515 516 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 517 518 inline G4double G4VEmProcess::LogEkin(const G4Track& track) 519 { 520 return track.GetDynamicParticle()->GetLogKineticEnergy(); 521 } 522 523 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 524 525 inline G4double G4VEmProcess::GetLambdaFromTable(G4double e, G4double loge) 526 { 527 return ((*theLambdaTable)[basedCoupleIndex])->LogVectorValue(e, loge); 528 } 529 530 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 531 532 inline G4double G4VEmProcess::GetLambdaFromTablePrim(G4double e) 533 { 534 return ((*theLambdaTablePrim)[basedCoupleIndex])->Value(e, idxLambda)/e; 535 } 536 537 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 538 539 inline G4double G4VEmProcess::GetLambdaFromTablePrim(G4double e, G4double loge) 540 { 541 return ((*theLambdaTablePrim)[basedCoupleIndex])->LogVectorValue(e, loge)/e; 542 } 543 544 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 545 546 inline G4double G4VEmProcess::ComputeCurrentLambda(G4double e) 547 { 548 return currentModel->CrossSectionPerVolume(baseMaterial, currentParticle, e); 549 } 550 551 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 552 553 inline G4double G4VEmProcess::GetCurrentLambda(G4double e) 554 { 555 if(currentCoupleIndex != coupleIdxLambda || fLambdaEnergy != e) { 556 coupleIdxLambda = currentCoupleIndex; 557 fLambdaEnergy = e; 558 if(e >= minKinEnergyPrim) { fLambda = GetLambdaFromTablePrim(e); } 559 else if(nullptr != theLambdaTable) { fLambda = GetLambdaFromTable(e); } 560 else { fLambda = ComputeCurrentLambda(e); } 561 fLambda *= fFactor; 562 } 563 return fLambda; 564 } 565 566 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 567 568 inline G4double G4VEmProcess::GetCurrentLambda(G4double e, G4double loge) 569 { 570 if(currentCoupleIndex != coupleIdxLambda || fLambdaEnergy != e) { 571 coupleIdxLambda = currentCoupleIndex; 572 fLambdaEnergy = e; 573 if(e >= minKinEnergyPrim) { fLambda = GetLambdaFromTablePrim(e, loge); } 574 else if(nullptr != theLambdaTable) { fLambda = GetLambdaFromTable(e, loge); } 575 else { fLambda = ComputeCurrentLambda(e); } 576 fLambda *= fFactor; 577 } 578 return fLambda; 579 } 580 581 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 582 583 inline void 584 G4VEmProcess::CurrentSetup(const G4MaterialCutsCouple* couple, G4double energy) 585 { 586 DefineMaterial(couple); 587 SelectModel(energy*massRatio, currentCoupleIndex); 588 } 589 590 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 591 592 inline G4double 593 G4VEmProcess::GetLambda(G4double kinEnergy, const G4MaterialCutsCouple* couple, 594 G4double logKinEnergy) 595 { 596 CurrentSetup(couple, kinEnergy); 597 return GetCurrentLambda(kinEnergy, logKinEnergy); 598 } 599 600 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 601 602 G4double G4VEmProcess::MeanFreePath(const G4Track& track) 603 { 604 const G4double kinEnergy = track.GetKineticEnergy(); 605 CurrentSetup(track.GetMaterialCutsCouple(), kinEnergy); 606 const G4double xs = GetCurrentLambda(kinEnergy, 607 track.GetDynamicParticle()->GetLogKineticEnergy()); 608 return (0.0 < xs) ? 1.0/xs : DBL_MAX; 609 } 610 611 // ======== Get/Set inline methods used at initialisation ================ 612 613 inline G4bool G4VEmProcess::ApplyCuts() const 614 { 615 return applyCuts; 616 } 617 618 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 619 620 inline G4int G4VEmProcess::LambdaBinning() const 621 { 622 return nLambdaBins; 623 } 624 625 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 626 627 inline G4double G4VEmProcess::MinKinEnergy() const 628 { 629 return minKinEnergy; 630 } 631 632 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 633 634 inline G4double G4VEmProcess::MaxKinEnergy() const 635 { 636 return maxKinEnergy; 637 } 638 639 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 640 641 inline G4double G4VEmProcess::CrossSectionBiasingFactor() const 642 { 643 return biasFactor; 644 } 645 646 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 647 648 inline G4PhysicsTable* G4VEmProcess::LambdaTable() const 649 { 650 return theLambdaTable; 651 } 652 653 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 654 655 inline G4PhysicsTable* G4VEmProcess::LambdaTablePrim() const 656 { 657 return theLambdaTablePrim; 658 } 659 660 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 661 662 inline void G4VEmProcess::SetLambdaTable(G4PhysicsTable* ptr) 663 { 664 theLambdaTable = ptr; 665 } 666 667 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 668 669 inline void G4VEmProcess::SetLambdaTablePrim(G4PhysicsTable* ptr) 670 { 671 theLambdaTablePrim = ptr; 672 } 673 674 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 675 676 inline std::vector<G4double>* G4VEmProcess::EnergyOfCrossSectionMax() const 677 { 678 return theEnergyOfCrossSectionMax; 679 } 680 681 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 682 683 inline void 684 G4VEmProcess::SetEnergyOfCrossSectionMax(std::vector<G4double>* ptr) 685 { 686 theEnergyOfCrossSectionMax = ptr; 687 } 688 689 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 690 691 inline const G4ParticleDefinition* G4VEmProcess::Particle() const 692 { 693 return particle; 694 } 695 696 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 697 698 inline const G4ParticleDefinition* G4VEmProcess::SecondaryParticle() const 699 { 700 return secondaryParticle; 701 } 702 703 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 704 705 inline void G4VEmProcess::SetCrossSectionType(G4CrossSectionType val) 706 { 707 fXSType = val; 708 } 709 710 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 711 712 inline G4CrossSectionType G4VEmProcess::CrossSectionType() const 713 { 714 return fXSType; 715 } 716 717 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 718 719 inline void G4VEmProcess::SetBuildTableFlag(G4bool val) 720 { 721 buildLambdaTable = val; 722 } 723 724 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 725 726 inline G4ParticleChangeForGamma* G4VEmProcess::GetParticleChange() 727 { 728 return &fParticleChange; 729 } 730 731 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 732 733 inline void G4VEmProcess::SetParticle(const G4ParticleDefinition* p) 734 { 735 particle = p; 736 currentParticle = p; 737 } 738 739 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 740 741 inline void G4VEmProcess::SetSecondaryParticle(const G4ParticleDefinition* p) 742 { 743 secondaryParticle = p; 744 } 745 746 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 747 748 inline void G4VEmProcess::SetStartFromNullFlag(G4bool val) 749 { 750 startFromNull = val; 751 } 752 753 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 754 755 inline void G4VEmProcess::SetSplineFlag(G4bool val) 756 { 757 splineFlag = val; 758 } 759 760 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 761 762 inline G4int G4VEmProcess::DensityIndex(G4int idx) const 763 { 764 return (*theDensityIdx)[idx]; 765 } 766 767 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 768 769 inline G4double G4VEmProcess::DensityFactor(G4int idx) const 770 { 771 return (*theDensityFactor)[idx]; 772 } 773 774 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 775 776 inline G4bool G4VEmProcess::UseBaseMaterial() const 777 { 778 return baseMat; 779 } 780 781 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 782 783 inline const G4VEmModel* G4VEmProcess::GetCurrentModel() const 784 { 785 return currentModel; 786 } 787 788 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 789 790 inline void G4VEmProcess::SetEmMasterProcess(const G4VEmProcess* ptr) 791 { 792 masterProc = ptr; 793 } 794 795 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 796 797 inline G4int G4VEmProcess::NumberOfModels() const 798 { 799 return numberOfModels; 800 } 801 802 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 803 804 inline G4VEmModel* G4VEmProcess::EmModel(std::size_t index) const 805 { 806 return (index < emModels.size()) ? emModels[index] : nullptr; 807 } 808 809 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 810 811 inline G4VEmModel* G4VEmProcess::GetModelByIndex(G4int idx, G4bool ver) const 812 { 813 return modelManager->GetModel(idx, ver); 814 } 815 816 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 817 818 #endif 819