Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VEnergyLossProcess.hh
Go to the documentation of this file.
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
55#include "globals.hh"
56#include "G4Material.hh"
58#include "G4Track.hh"
59#include "G4EmModelManager.hh"
61#include "G4EmTableType.hh"
63#include "G4PhysicsTable.hh"
64#include "G4PhysicsVector.hh"
65
66class G4Step;
68class G4EmParameters;
69class G4VEmModel;
71class G4DataVector;
72class G4Region;
73class G4SafetyHelper;
78class G4EmDataHandler;
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81
83{
84public:
85
86 G4VEnergyLossProcess(const G4String& name = "EnergyLoss",
88
89 ~G4VEnergyLossProcess() override;
90
91 //------------------------------------------------------------------------
92 // Virtual methods to be implemented in concrete processes
93 //------------------------------------------------------------------------
94
95protected:
96
97 // description of specific process parameters
98 virtual void StreamProcessInfo(std::ostream&) const {};
99
101 const G4ParticleDefinition*) = 0;
102
103public:
104
105 // used as low energy limit LambdaTable
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
120
121 // build a table
123
124 // Called before tracking of each new G4Track
125 void StartTracking(G4Track*) override;
126
127 // Step limit from AlongStep
129 const G4Track&,
130 G4double previousStepSize,
131 G4double currentMinimumStep,
132 G4double& currentSafety,
133 G4GPILSelection* selection) override;
134
135 // Step limit from cross section
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
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.
157 const G4String& directory,
158 G4bool ascii) override;
159
160private:
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
171public:
172
173 // access to dispersion of restricted energy loss
175 const G4DynamicParticle* dp,
176 G4double length);
177
178 // Access to cross section table
180 const G4MaterialCutsCouple* couple);
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
190 G4double previousStepSize,
191 G4double currentMinimumStep,
192 G4double& currentSafety);
193
194protected:
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
203 G4double previousStepSize,
204 G4double currentMinimumStep,
205 G4double& currentSafety) override;
206
207 // creation of an empty vector for cross sections for derived processes
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
220public:
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
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
248
249 // Return the assigned fluctuation model
250 inline G4VEmFluctuationModel* FluctModel() const;
251
252 //------------------------------------------------------------------------
253 // Define and access particle type
254 //------------------------------------------------------------------------
255
256protected:
257 inline void SetParticle(const G4ParticleDefinition* p);
258 inline void SetSecondaryParticle(const G4ParticleDefinition* p);
259
260public:
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
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
281 const G4String& region,
282 G4bool flag = true);
283
284 void ActivateSecondaryBiasing(const G4String& region, G4double factor,
285 G4double energyLimit);
286
287 // obsolete method will be removed in the next major release
289
290 // interpolation of cross section
291 inline void SetSpline(G4bool val);
294
295 // Set/Get flag "isIonisation"
296 void SetIonisation(G4bool val);
297 inline G4bool IsIonisationProcess() const;
298
299 // Redefine parameteters of stepping control
301 void SetStepFunction(G4double v1, G4double v2);
303
304 inline G4int NumberOfSubCutoffRegions() const;
305
306 //------------------------------------------------------------------------
307 // Specific methods to path Physics Tables to the process
308 //------------------------------------------------------------------------
309
311 void SetCSDARangeTable(G4PhysicsTable* pRange);
315
316 // set properties of cross section shape
317 void SetTwoPeaksXS(std::vector<G4TwoPeaksXS*>*);
318 void SetEnergyOfCrossSectionMax(std::vector<G4double>*);
319
320 //------------------------------------------------------------------------
321 // Specific methods to define custom Physics Tables to the process
322 //------------------------------------------------------------------------
323
324 // Binning for dEdx, range, inverse range and lambda tables
325 void SetDEDXBinning(G4int nbins);
326
327 // Min kinetic energy for tables
329 inline G4double MinKinEnergy() const;
330
331 // Max kinetic energy for tables
333 inline G4double MaxKinEnergy() const;
334
335 // Biasing parameters
336 inline G4double CrossSectionBiasingFactor() const;
337
338 // Return values for given G4MaterialCutsCouple
339 inline G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple*);
340 inline G4double GetCSDADEDX(G4double kineticEnergy,
341 const G4MaterialCutsCouple*);
342 inline G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple*,
343 G4double logKineticEnergy);
344 inline G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple*);
345 inline G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple*,
346 G4double logKineticEnergy);
347 inline G4double GetCSDARange(G4double kineticEnergy,
348 const G4MaterialCutsCouple*);
349 inline G4double GetKineticEnergy(G4double range,
350 const G4MaterialCutsCouple*);
351 inline G4double GetLambda(G4double kineticEnergy,const G4MaterialCutsCouple*);
352 inline G4double GetLambda(G4double kineticEnergy,const G4MaterialCutsCouple*,
353 G4double logKineticEnergy);
354
355 inline G4bool TablesAreBuilt() const;
356
357 // Access to specific tables
358 inline G4PhysicsTable* DEDXTable() const;
360 inline G4PhysicsTable* IonisationTable() const;
361 inline G4PhysicsTable* CSDARangeTable() const;
362 inline G4PhysicsTable* RangeTableForLoss() const;
363 inline G4PhysicsTable* InverseRangeTable() const;
364 inline G4PhysicsTable* LambdaTable() const;
365 inline std::vector<G4TwoPeaksXS*>* TwoPeaksXS() const;
366 inline std::vector<G4double>* EnergyOfCrossSectionMax() const;
367
368 inline G4bool UseBaseMaterial() const;
369
370 //------------------------------------------------------------------------
371 // Run time method for simulation of ionisation
372 //------------------------------------------------------------------------
373
374 // access atom on which interaction happens
375 const G4Element* GetCurrentElement() const;
376
377 // Set scaling parameters for ions is needed to G4EmCalculator
378 void SetDynamicMassCharge(G4double massratio, G4double charge2ratio);
379
380private:
381
382 void FillSecondariesAlongStep(G4double weight);
383
384 void PrintWarning(const G4String&, G4double val) const;
385
386 // define material and indexes
387 inline void DefineMaterial(const G4MaterialCutsCouple* couple);
388
389 //------------------------------------------------------------------------
390 // Compute values using scaling relation, mass and charge of based particle
391 //------------------------------------------------------------------------
392 inline G4double GetDEDXForScaledEnergy(G4double scaledKinE);
393 inline G4double GetDEDXForScaledEnergy(G4double scaledKinE,
394 G4double logScaledKinE);
395 inline G4double GetIonisationForScaledEnergy(G4double scaledKinE);
396 inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinE);
397 inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinE,
398 G4double logScaledKinE);
399
400 inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinE);
401 inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinE,
402 G4double logScaledKinE);
403
404 inline G4double ScaledKinEnergyForLoss(G4double range);
405 inline G4double GetLambdaForScaledEnergy(G4double scaledKinE);
406 inline G4double GetLambdaForScaledEnergy(G4double scaledKinE,
407 G4double logScaledKinE);
408
409 inline G4double LogScaledEkin(const G4Track& aTrack);
410
411 void ComputeLambdaForScaledEnergy(G4double scaledKinE,
412 const G4Track& aTrack);
413
414 G4bool IsRegionForCubcutProcessor(const G4Track& aTrack);
415
416protected:
417
419 const G4Material* currentMaterial = nullptr;
421
422private:
423
424 G4LossTableManager* lManager;
425 G4EmModelManager* modelManager;
426 G4VEmModel* currentModel = nullptr;
427 G4EmBiasingManager* biasManager = nullptr;
428 G4SafetyHelper* safetyHelper;
429 G4EmParameters* theParameters;
430 G4VEmFluctuationModel* fluctModel = nullptr;
431 G4VAtomDeexcitation* atomDeexcitation = nullptr;
432 G4VSubCutProducer* subcutProducer = nullptr;
433
434 const G4ParticleDefinition* particle = nullptr;
435 const G4ParticleDefinition* baseParticle = nullptr;
436 const G4ParticleDefinition* secondaryParticle = nullptr;
437 G4EmDataHandler* theData = nullptr;
438
439 G4PhysicsTable* theDEDXTable = nullptr;
440 G4PhysicsTable* theDEDXunRestrictedTable = nullptr;
441 G4PhysicsTable* theIonisationTable = nullptr;
442 G4PhysicsTable* theRangeTableForLoss = nullptr;
443 G4PhysicsTable* theCSDARangeTable = nullptr;
444 G4PhysicsTable* theInverseRangeTable = nullptr;
445 G4PhysicsTable* theLambdaTable = nullptr;
446
447 std::vector<const G4Region*>* scoffRegions = nullptr;
448 std::vector<G4VEmModel*>* emModels = nullptr;
449 const std::vector<G4int>* theDensityIdx = nullptr;
450 const std::vector<G4double>* theDensityFactor = nullptr;
451 const G4DataVector* theCuts = nullptr;
452 const std::vector<G4bool>* theFluctuationFlags = nullptr;
453
454 std::vector<G4double>* theEnergyOfCrossSectionMax = nullptr;
455 std::vector<G4TwoPeaksXS*>* fXSpeaks = nullptr;
456
457 G4double lowestKinEnergy;
458 G4double minKinEnergy;
459 G4double maxKinEnergy;
460 G4double maxKinEnergyCSDA;
461
462 G4double linLossLimit = 0.01;
463 G4double dRoverRange = 0.2;
464 G4double finalRange;
465 G4double lambdaFactor = 0.8;
466 G4double invLambdaFactor;
467 G4double biasFactor = 1.0;
468
469 G4double massRatio = 1.0;
470 G4double logMassRatio = 0.0;
471 G4double fFactor = 1.0;
472 G4double reduceFactor = 1.0;
473 G4double chargeSqRatio = 1.0;
474 G4double fRange = 0.0;
475 G4double fRangeEnergy = 0.0;
476
477protected:
478
483
484 std::size_t currentCoupleIndex = 0;
485
486private:
487
488 G4int nBins;
489 G4int nBinsCSDA;
490 G4int numberOfModels = 0;
491 G4int nSCoffRegions = 0;
492 G4int secID = _DeltaElectron;
493 G4int tripletID = _TripletElectron;
494 G4int biasID = _DeltaEBelowCut;
495 G4int epixeID = _ePIXE;
496 G4int gpixeID = _GammaPIXE;
497 G4int mainSecondaries = 1;
498
499 std::size_t basedCoupleIndex = 0;
500 std::size_t coupleIdxRange = 0;
501 std::size_t idxDEDX = 0;
502 std::size_t idxDEDXunRestricted = 0;
503 std::size_t idxIonisation = 0;
504 std::size_t idxRange = 0;
505 std::size_t idxCSDA = 0;
506 std::size_t idxSecRange = 0;
507 std::size_t idxInverseRange = 0;
508 std::size_t idxLambda = 0;
509
510 G4GPILSelection aGPILSelection;
512
513 G4bool lossFluctuationFlag = true;
514 G4bool useCutAsFinalRange = false;
515 G4bool tablesAreBuilt = false;
516 G4bool spline = true;
517 G4bool isIon = false;
518 G4bool isIonisation = false;
519 G4bool useDeexcitation = false;
520 G4bool biasFlag = false;
521 G4bool weightFlag = false;
522 G4bool isMaster = false;
523 G4bool baseMat = false;
524
525 // flags allowing define table parameters individual for the process
526 G4bool actLinLossLimit = false;
527 G4bool actBinning = false;
528 G4bool actMinKinEnergy = false;
529 G4bool actMaxKinEnergy = false;
530
531 std::vector<G4DynamicParticle*> secParticles;
532 std::vector<G4Track*> scTracks;
533};
534
535// ======== Run time inline methods ================
536
538{
539 return currentCoupleIndex;
540}
541
542//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
543
545{
546 currentModel = modelManager->SelectModel(kinEnergy, currentCoupleIndex);
547 currentModel->SetCurrentCouple(currentCouple);
548}
549
550//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
551
553 G4double kinEnergy, std::size_t& idx) const
554{
555 return modelManager->SelectModel(kinEnergy, idx);
556}
557
558//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
559
560inline void
561G4VEnergyLossProcess::DefineMaterial(const G4MaterialCutsCouple* couple)
562{
563 if(couple != currentCouple) {
564 currentCouple = couple;
565 currentMaterial = couple->GetMaterial();
566 basedCoupleIndex = currentCoupleIndex = couple->GetIndex();
567 fFactor = chargeSqRatio*biasFactor;
569 idxLambda = 0;
570 if(baseMat) {
571 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
572 fFactor *= (*theDensityFactor)[currentCoupleIndex];
573 }
574 reduceFactor = 1.0/(fFactor*massRatio);
575 }
576}
577
578//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
579
580inline G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e)
581{
582 /*
583 G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
584 << basedCoupleIndex << " E(MeV)= " << e
585 << " Emin= " << minKinEnergy << " Factor= " << fFactor
586 << " " << theDEDXTable << G4endl; */
587 G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->Value(e, idxDEDX);
588 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
589 return x;
590}
591
592//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
593
594inline
595G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e, G4double loge)
596{
597 /*
598 G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
599 << basedCoupleIndex << " E(MeV)= " << e
600 << " Emin= " << minKinEnergy << " Factor= " << fFactor
601 << " " << theDEDXTable << G4endl; */
602 G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->LogVectorValue(e,loge);
603 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
604 return x;
605}
606
607//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
608
609inline G4double G4VEnergyLossProcess::GetIonisationForScaledEnergy(G4double e)
610{
611 G4double x =
612 fFactor*(*theIonisationTable)[basedCoupleIndex]->Value(e, idxIonisation);
613 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
614 return x;
615}
616
617//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
618
619inline G4double G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e)
620{
621 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
622 // << basedCoupleIndex << " E(MeV)= " << e
623 // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
624 if(currentCoupleIndex != coupleIdxRange || fRangeEnergy != e) {
625 coupleIdxRange = currentCoupleIndex;
626 fRangeEnergy = e;
627 fRange = reduceFactor*((*theRangeTableForLoss)[basedCoupleIndex])->Value(e, idxRange);
628 if (fRange < 0.0) { fRange = 0.0; }
629 else if (e < minKinEnergy) { fRange *= std::sqrt(e/minKinEnergy); }
630 }
631 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
632 // << basedCoupleIndex << " E(MeV)= " << e
633 // << " R= " << computedRange << " " << theRangeTableForLoss << G4endl;
634 return fRange;
635}
636
637inline G4double
638G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e, G4double loge)
639{
640 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
641 // << basedCoupleIndex << " E(MeV)= " << e
642 // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
643 if(currentCoupleIndex != coupleIdxRange || fRangeEnergy != e) {
644 coupleIdxRange = currentCoupleIndex;
645 fRangeEnergy = e;
646 fRange = reduceFactor*((*theRangeTableForLoss)[basedCoupleIndex])->LogVectorValue(e, loge);
647 if (fRange < 0.0) { fRange = 0.0; }
648 else if (e < minKinEnergy) { fRange *= std::sqrt(e/minKinEnergy); }
649 }
650 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
651 // << basedCoupleIndex << " E(MeV)= " << e
652 // << " R= " << fRange << " " << theRangeTableForLoss << G4endl;
653 return fRange;
654}
655
656//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
657
658inline G4double
659G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e)
660{
661 G4double x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e, idxCSDA);
662 if (x < 0.0) { x = 0.0; }
663 else if (e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
664 return x;
665}
666
667//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
668
669inline G4double
670G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e,
671 G4double loge)
672{
673 G4double x = ((*theCSDARangeTable)[basedCoupleIndex])->LogVectorValue(e, loge);
674 if (x < 0.0) { x = 0.0; }
675 else if (e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
676 return x;
677}
678
679//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
680
681inline G4double G4VEnergyLossProcess::ScaledKinEnergyForLoss(G4double r)
682{
683 //G4cout << "G4VEnergyLossProcess::GetEnergy: Idx= "
684 // << basedCoupleIndex << " R(mm)= " << r << " "
685 // << theInverseRangeTable << G4endl;
686 G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex];
687 G4double rmin = v->Energy(0);
688 G4double e = 0.0;
689 if(r >= rmin) { e = v->Value(r, idxInverseRange); }
690 else if(r > 0.0) {
691 G4double x = r/rmin;
692 e = minKinEnergy*x*x;
693 }
694 return e;
695}
696
697//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
698
699inline G4double G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e)
700{
701 return fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
702}
703
704//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
705
706inline G4double
707G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e, G4double loge)
708{
709 return fFactor*((*theLambdaTable)[basedCoupleIndex])->LogVectorValue(e, loge);
710}
711
712//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
713
714inline G4double G4VEnergyLossProcess::LogScaledEkin(const G4Track& track)
715{
716 return track.GetDynamicParticle()->GetLogKineticEnergy() + logMassRatio;
717}
718
719//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
720
721inline G4double
723 const G4MaterialCutsCouple* couple)
724{
725 DefineMaterial(couple);
726 return GetDEDXForScaledEnergy(kinEnergy*massRatio);
727}
728
729//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
730
731inline G4double
733 const G4MaterialCutsCouple* couple,
734 G4double logKinEnergy)
735{
736 DefineMaterial(couple);
737 return GetDEDXForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
738}
739
740//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
741
742inline G4double
744 const G4MaterialCutsCouple* couple)
745{
746 DefineMaterial(couple);
747 return GetScaledRangeForScaledEnergy(kinEnergy*massRatio);
748}
749
750//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
751
752inline G4double
754 const G4MaterialCutsCouple* couple,
755 G4double logKinEnergy)
756{
757 DefineMaterial(couple);
758 return GetScaledRangeForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
759}
760
761//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
762
763inline G4double
765 const G4MaterialCutsCouple* couple)
766{
767 DefineMaterial(couple);
768 return (nullptr == theCSDARangeTable) ? DBL_MAX :
769 GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
770}
771
772//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
773
774inline G4double
776 const G4MaterialCutsCouple* couple)
777{
778 DefineMaterial(couple);
779 return ScaledKinEnergyForLoss(range/reduceFactor)/massRatio;
780}
781
782//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
783
784inline G4double
786 const G4MaterialCutsCouple* couple)
787{
788 DefineMaterial(couple);
789 return (nullptr != theLambdaTable) ?
790 GetLambdaForScaledEnergy(kinEnergy*massRatio) : 0.0;
791}
792
793//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
794
795inline G4double
797 const G4MaterialCutsCouple* couple,
798 G4double logKinEnergy)
799{
800 DefineMaterial(couple);
801 return (nullptr != theLambdaTable) ?
802 GetLambdaForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio)
803 : 0.0;
804}
805
806// ======== Get/Set inline methods used at initialisation ================
807
809{
810 fluctModel = p;
811}
812
813//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
814
816{
817 return fluctModel;
818}
819
820//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
821
823{
824 particle = p;
825}
826
827//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
828
829inline void
831{
832 secondaryParticle = p;
833}
834
835//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
836
837inline void
839{
840 baseParticle = p;
841}
842
843//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
844
846{
847 return particle;
848}
849
850//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
851
853{
854 return baseParticle;
855}
856
857//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
858
859inline const G4ParticleDefinition*
861{
862 return secondaryParticle;
863}
864
865//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
866
868{
869 spline = val;
870}
871
872//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
873
875{
876 fXSType = val;
877}
878
879//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
880
882{
883 return fXSType;
884}
885
886//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
887
889{
890 return isIonisation;
891}
892
893//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
894
896{
897 return nSCoffRegions;
898}
899
900//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
901
903{
904 return minKinEnergy;
905}
906
907//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
908
910{
911 return maxKinEnergy;
912}
913
914//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
915
917{
918 return biasFactor;
919}
920
921//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
922
924{
925 return tablesAreBuilt;
926}
927
928//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
929
931{
932 return theDEDXTable;
933}
934
935//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
936
938{
939 return theDEDXunRestrictedTable;
940}
941
942//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
943
945{
946 return theIonisationTable;
947}
948
949//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
950
952{
953 return theCSDARangeTable;
954}
955
956//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
957
959{
960 return theRangeTableForLoss;
961}
962
963//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
964
966{
967 return theInverseRangeTable;
968}
969
970//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
971
973{
974 return theLambdaTable;
975}
976
977//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
978
980{
981 return baseMat;
982}
983
984//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
985
986inline std::vector<G4double>*
988{
989 return theEnergyOfCrossSectionMax;
990}
991
992//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
993
994inline std::vector<G4TwoPeaksXS*>* G4VEnergyLossProcess::TwoPeaksXS() const
995{
996 return fXSpeaks;
997}
998
999//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1000
1002{
1003 return numberOfModels;
1004}
1005
1006//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1007
1008inline G4VEmModel* G4VEnergyLossProcess::EmModel(std::size_t index) const
1009{
1010 return (index < emModels->size()) ? (*emModels)[index] : nullptr;
1011}
1012
1013//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1014
1015inline G4VEmModel*
1017{
1018 return modelManager->GetModel((G4int)idx, ver);
1019}
1020
1021//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1022
1023#endif
G4EmTableType
@ fRestricted
G4CrossSectionType
@ fEmOnePeak
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4GPILSelection
G4ProcessType
@ fElectromagnetic
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4double GetLogKineticEnergy() const
const G4Material * GetMaterial() const
G4double Energy(const std::size_t index) const
G4double Value(const G4double energy, std::size_t &lastidx) const
G4Region defines a region or a group of regions in the detector geometry setup, sharing properties as...
Definition G4Region.hh:90
G4SafetyHelper is a helper class for physics processes which require knowledge of the safety,...
const G4DynamicParticle * GetDynamicParticle() const
G4VContinuousDiscreteProcess(const G4String &, G4ProcessType aType=fNotDefined)
const G4ParticleDefinition * BaseParticle() const
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
G4PhysicsTable * RangeTableForLoss() const
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
G4ParticleChangeForLoss fParticleChange
void PreparePhysicsTable(const G4ParticleDefinition &) override
std::vector< G4double > * EnergyOfCrossSectionMax() const
G4double GetCSDADEDX(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4PhysicsTable * InverseRangeTable() const
G4double MeanFreePath(const G4Track &track)
void SetFluctModel(G4VEmFluctuationModel *)
G4VEnergyLossProcess(G4VEnergyLossProcess &)=delete
G4double GetKineticEnergy(G4double range, const G4MaterialCutsCouple *)
G4PhysicsTable * CSDARangeTable() const
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *)=0
void SelectModel(G4double kinEnergy)
void SetRangeTableForLoss(G4PhysicsTable *p)
G4int NumberOfSubCutoffRegions() const
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *, G4double cut)
void ProcessDescription(std::ostream &outFile) const override
std::size_t NumberOfModels() const
G4VEmModel * EmModel(std::size_t index=0) const
G4VEmModel * GetModelByIndex(std::size_t idx=0, G4bool ver=false) const
G4PhysicsTable * BuildDEDXTable(G4EmTableType tType=fRestricted)
void SetTwoPeaksXS(std::vector< G4TwoPeaksXS * > *)
const G4MaterialCutsCouple * currentCouple
std::size_t CurrentMaterialCutsCoupleIndex() const
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *, G4double cut)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
const G4Element * GetCurrentElement() const
void SetDEDXBinning(G4int nbins)
void SetStepFunction(G4double v1, G4double v2)
G4double GetLambda(G4double kineticEnergy, const G4MaterialCutsCouple *)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
const G4Material * currentMaterial
G4VEnergyLossProcess & operator=(const G4VEnergyLossProcess &right)=delete
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
std::vector< G4TwoPeaksXS * > * TwoPeaksXS() const
void SetParticle(const G4ParticleDefinition *p)
void SetCrossSectionType(G4CrossSectionType val)
const G4ParticleDefinition * Particle() const
void SetInverseRangeTable(G4PhysicsTable *p)
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
void ActivateForcedInteraction(G4double length, const G4String &region, G4bool flag=true)
G4CrossSectionType CrossSectionType() const
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
void SetBaseParticle(const G4ParticleDefinition *p)
G4double CrossSectionBiasingFactor() const
void SetEnergyOfCrossSectionMax(std::vector< G4double > *)
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple)
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
G4PhysicsTable * BuildLambdaTable(G4EmTableType tType=fRestricted)
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4double GetCSDARange(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
G4VEmFluctuationModel * FluctModel() const
void SetLinearLossLimit(G4double val)
void SetLambdaTable(G4PhysicsTable *p)
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4VEnergyLossProcess(const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
G4PhysicsTable * IonisationTable() const
void ActivateSubCutoff(const G4Region *region)
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, std::size_t &idxCouple) const
void SetSecondaryParticle(const G4ParticleDefinition *p)
G4PhysicsTable * LambdaTable() const
void SetCSDARangeTable(G4PhysicsTable *pRange)
G4PhysicsTable * DEDXunRestrictedTable() const
virtual void StreamProcessInfo(std::ostream &) const
G4PhysicsTable * DEDXTable() const
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
void StartTracking(G4Track *) override
const G4ParticleDefinition * SecondaryParticle() const
#define DBL_MAX
Definition templates.hh:62