Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 inline void SetLossFluctuations(G4bool val);
288
289 inline void SetSpline(G4bool val);
292
293 // Set/Get flag "isIonisation"
294 void SetIonisation(G4bool val);
295 inline G4bool IsIonisationProcess() const;
296
297 // Redefine parameteters for stepping control
299 void SetStepFunction(G4double v1, G4double v2);
301
302 inline G4int NumberOfSubCutoffRegions() const;
303
304 //------------------------------------------------------------------------
305 // Specific methods to path Physics Tables to the process
306 //------------------------------------------------------------------------
307
309 void SetCSDARangeTable(G4PhysicsTable* pRange);
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
326 inline G4double MinKinEnergy() const;
327
328 // Max kinetic energy for tables
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;
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
377private:
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 void ComputeLambdaForScaledEnergy(G4double scaledKinE,
407 G4double logScaledKinE);
408
409 G4bool IsRegionForCubcutProcessor(const G4Track& aTrack);
410
411protected:
412
414 const G4Material* currentMaterial = nullptr;
416
417private:
418
419 G4LossTableManager* lManager;
420 G4EmModelManager* modelManager;
421 G4VEmModel* currentModel = nullptr;
422 G4EmBiasingManager* biasManager = nullptr;
423 G4SafetyHelper* safetyHelper;
424 G4EmParameters* theParameters;
425 G4VEmFluctuationModel* fluctModel = nullptr;
426 G4VAtomDeexcitation* atomDeexcitation = nullptr;
427 G4VSubCutProducer* subcutProducer = nullptr;
428
429 const G4ParticleDefinition* particle = nullptr;
430 const G4ParticleDefinition* baseParticle = nullptr;
431 const G4ParticleDefinition* secondaryParticle = nullptr;
432 G4EmDataHandler* theData = nullptr;
433
434 G4PhysicsTable* theDEDXTable = nullptr;
435 G4PhysicsTable* theDEDXunRestrictedTable = nullptr;
436 G4PhysicsTable* theIonisationTable = nullptr;
437 G4PhysicsTable* theRangeTableForLoss = nullptr;
438 G4PhysicsTable* theCSDARangeTable = nullptr;
439 G4PhysicsTable* theInverseRangeTable = nullptr;
440 G4PhysicsTable* theLambdaTable = nullptr;
441
442 std::vector<const G4Region*>* scoffRegions = nullptr;
443 std::vector<G4VEmModel*>* emModels = nullptr;
444 const std::vector<G4int>* theDensityIdx = nullptr;
445 const std::vector<G4double>* theDensityFactor = nullptr;
446 const G4DataVector* theCuts = nullptr;
447
448 std::vector<G4double>* theEnergyOfCrossSectionMax = nullptr;
449 std::vector<G4TwoPeaksXS*>* fXSpeaks = nullptr;
450
451 G4double lowestKinEnergy;
452 G4double minKinEnergy;
453 G4double maxKinEnergy;
454 G4double maxKinEnergyCSDA;
455
456 G4double linLossLimit = 0.01;
457 G4double dRoverRange = 0.2;
458 G4double finalRange;
459 G4double lambdaFactor = 0.8;
460 G4double invLambdaFactor;
461 G4double biasFactor = 1.0;
462
463 G4double massRatio = 1.0;
464 G4double logMassRatio = 0.0;
465 G4double fFactor = 1.0;
466 G4double reduceFactor = 1.0;
467 G4double chargeSqRatio = 1.0;
468 G4double fRange = 0.0;
469 G4double fRangeEnergy = 0.0;
470
471protected:
472
479
480 std::size_t currentCoupleIndex = 0;
481
482private:
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 mainSecondaries = 1;
492
493 std::size_t basedCoupleIndex = 0;
494 std::size_t coupleIdxRange = 0;
495 std::size_t idxDEDX = 0;
496 std::size_t idxDEDXunRestricted = 0;
497 std::size_t idxIonisation = 0;
498 std::size_t idxRange = 0;
499 std::size_t idxCSDA = 0;
500 std::size_t idxSecRange = 0;
501 std::size_t idxInverseRange = 0;
502 std::size_t idxLambda = 0;
503
504 G4GPILSelection aGPILSelection;
506
507 G4bool lossFluctuationFlag = true;
508 G4bool rndmStepFlag = false;
509 G4bool tablesAreBuilt = false;
510 G4bool spline = true;
511 G4bool isIon = false;
512 G4bool isIonisation = true;
513 G4bool useDeexcitation = false;
514 G4bool biasFlag = false;
515 G4bool weightFlag = false;
516 G4bool isMaster = true;
517 G4bool baseMat = false;
518 G4bool actLinLossLimit = false;
519 G4bool actLossFluc = false;
520 G4bool actBinning = false;
521 G4bool actMinKinEnergy = false;
522 G4bool actMaxKinEnergy = false;
523
524 std::vector<G4DynamicParticle*> secParticles;
525 std::vector<G4Track*> scTracks;
526};
527
528// ======== Run time inline methods ================
529
531{
532 return currentCoupleIndex;
533}
534
535//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
536
538{
539 currentModel = modelManager->SelectModel(kinEnergy, currentCoupleIndex);
540 currentModel->SetCurrentCouple(currentCouple);
541}
542
543//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
544
546 G4double kinEnergy, std::size_t& idx) const
547{
548 return modelManager->SelectModel(kinEnergy, idx);
549}
550
551//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
552
553inline void
554G4VEnergyLossProcess::DefineMaterial(const G4MaterialCutsCouple* couple)
555{
556 if(couple != currentCouple) {
557 currentCouple = couple;
558 currentMaterial = couple->GetMaterial();
559 basedCoupleIndex = currentCoupleIndex = couple->GetIndex();
560 fFactor = chargeSqRatio*biasFactor;
562 idxLambda = 0;
563 if(baseMat) {
564 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
565 fFactor *= (*theDensityFactor)[currentCoupleIndex];
566 }
567 reduceFactor = 1.0/(fFactor*massRatio);
568 }
569}
570
571//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
572
573inline G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e)
574{
575 /*
576 G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
577 << basedCoupleIndex << " E(MeV)= " << e
578 << " Emin= " << minKinEnergy << " Factor= " << fFactor
579 << " " << theDEDXTable << G4endl; */
580 G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->Value(e, idxDEDX);
581 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
582 return x;
583}
584
585inline
586G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e, G4double loge)
587{
588 /*
589 G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
590 << basedCoupleIndex << " E(MeV)= " << e
591 << " Emin= " << minKinEnergy << " Factor= " << fFactor
592 << " " << theDEDXTable << G4endl; */
593 G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->LogVectorValue(e,loge);
594 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
595 return x;
596}
597
598//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
599
600inline G4double G4VEnergyLossProcess::GetIonisationForScaledEnergy(G4double e)
601{
602 G4double x =
603 fFactor*(*theIonisationTable)[basedCoupleIndex]->Value(e, idxIonisation);
604 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
605 return x;
606}
607
608//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
609
610inline G4double G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e)
611{
612 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
613 // << basedCoupleIndex << " E(MeV)= " << e
614 // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
615 if(currentCoupleIndex != coupleIdxRange || fRangeEnergy != e) {
616 coupleIdxRange = currentCoupleIndex;
617 fRangeEnergy = e;
618 fRange = reduceFactor*((*theRangeTableForLoss)[basedCoupleIndex])->Value(e, idxRange);
619 if(e < minKinEnergy) { fRange *= std::sqrt(e/minKinEnergy); }
620 }
621 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
622 // << basedCoupleIndex << " E(MeV)= " << e
623 // << " R= " << computedRange << " " << theRangeTableForLoss << G4endl;
624 return fRange;
625}
626
627inline G4double
628G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e, G4double loge)
629{
630 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
631 // << basedCoupleIndex << " E(MeV)= " << e
632 // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
633 if(currentCoupleIndex != coupleIdxRange || fRangeEnergy != e) {
634 coupleIdxRange = currentCoupleIndex;
635 fRangeEnergy = e;
636 fRange = reduceFactor*((*theRangeTableForLoss)[basedCoupleIndex])->LogVectorValue(e, loge);
637 if(e < minKinEnergy) { fRange *= std::sqrt(e/minKinEnergy); }
638 }
639 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
640 // << basedCoupleIndex << " E(MeV)= " << e
641 // << " R= " << fRange << " " << theRangeTableForLoss << G4endl;
642 return fRange;
643}
644
645//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
646
647inline G4double
648G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e)
649{
650 G4double x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e, idxCSDA);
651 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
652 return x;
653}
654
655//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
656
657inline G4double
658G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e,
659 G4double loge)
660{
661 G4double x = ((*theCSDARangeTable)[basedCoupleIndex])->LogVectorValue(e, loge);
662 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
663 return x;
664}
665
666//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
667
668inline G4double G4VEnergyLossProcess::ScaledKinEnergyForLoss(G4double r)
669{
670 //G4cout << "G4VEnergyLossProcess::GetEnergy: Idx= "
671 // << basedCoupleIndex << " R(mm)= " << r << " "
672 // << theInverseRangeTable << G4endl;
673 G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex];
674 G4double rmin = v->Energy(0);
675 G4double e = 0.0;
676 if(r >= rmin) { e = v->Value(r, idxInverseRange); }
677 else if(r > 0.0) {
678 G4double x = r/rmin;
679 e = minKinEnergy*x*x;
680 }
681 return e;
682}
683
684//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
685
686inline G4double G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e)
687{
688 return fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
689}
690
691//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
692
693inline G4double
694G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e, G4double loge)
695{
696 return fFactor*((*theLambdaTable)[basedCoupleIndex])->LogVectorValue(e, loge);
697}
698
699//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
700
701inline G4double
703 const G4MaterialCutsCouple* couple)
704{
705 DefineMaterial(couple);
706 return GetDEDXForScaledEnergy(kinEnergy*massRatio);
707}
708
709//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
710
711inline G4double
713 const G4MaterialCutsCouple* couple,
714 G4double logKinEnergy)
715{
716 DefineMaterial(couple);
717 return GetDEDXForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
718}
719
720//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
721
722inline G4double
724 const G4MaterialCutsCouple* couple)
725{
726 DefineMaterial(couple);
727 return GetScaledRangeForScaledEnergy(kinEnergy*massRatio);
728}
729
730//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
731
732inline G4double
734 const G4MaterialCutsCouple* couple,
735 G4double logKinEnergy)
736{
737 DefineMaterial(couple);
738 return GetScaledRangeForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
739}
740
741//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
742
743inline G4double
745 const G4MaterialCutsCouple* couple)
746{
747 DefineMaterial(couple);
748 return (nullptr == theCSDARangeTable) ? DBL_MAX :
749 GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
750}
751
752//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
753
754inline G4double
756 const G4MaterialCutsCouple* couple)
757{
758 DefineMaterial(couple);
759 return ScaledKinEnergyForLoss(range/reduceFactor)/massRatio;
760}
761
762//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
763
764inline G4double
766 const G4MaterialCutsCouple* couple)
767{
768 DefineMaterial(couple);
769 return (nullptr != theLambdaTable) ?
770 GetLambdaForScaledEnergy(kinEnergy*massRatio) : 0.0;
771}
772
773//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
774
775inline G4double
777 const G4MaterialCutsCouple* couple,
778 G4double logKinEnergy)
779{
780 DefineMaterial(couple);
781 return (nullptr != theLambdaTable) ?
782 GetLambdaForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio)
783 : 0.0;
784}
785
786// ======== Get/Set inline methods used at initialisation ================
787
789{
790 fluctModel = p;
791}
792
793//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
794
796{
797 return fluctModel;
798}
799
800//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
801
803{
804 particle = p;
805}
806
807//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
808
809inline void
811{
812 secondaryParticle = p;
813}
814
815//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
816
817inline void
819{
820 baseParticle = p;
821}
822
823//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
824
826{
827 return particle;
828}
829
830//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
831
833{
834 return baseParticle;
835}
836
837//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
838
839inline const G4ParticleDefinition*
841{
842 return secondaryParticle;
843}
844
845//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
846
848{
849 lossFluctuationFlag = val;
850 actLossFluc = true;
851}
852
853//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
854
856{
857 spline = val;
858}
859
860//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
861
863{
864 fXSType = val;
865}
866
867//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
868
870{
871 return fXSType;
872}
873
874//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
875
877{
878 return isIonisation;
879}
880
881//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
882
884{
885 return nSCoffRegions;
886}
887
888//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
889
891{
892 return minKinEnergy;
893}
894
895//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
896
898{
899 return maxKinEnergy;
900}
901
902//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
903
905{
906 return biasFactor;
907}
908
909//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
910
912{
913 return tablesAreBuilt;
914}
915
916//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
917
919{
920 return theDEDXTable;
921}
922
923//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
924
926{
927 return theDEDXunRestrictedTable;
928}
929
930//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
931
933{
934 return theIonisationTable;
935}
936
937//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
938
940{
941 return theCSDARangeTable;
942}
943
944//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
945
947{
948 return theRangeTableForLoss;
949}
950
951//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
952
954{
955 return theInverseRangeTable;
956}
957
958//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
959
961{
962 return theLambdaTable;
963}
964
965//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
966
968{
969 return baseMat;
970}
971
972//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
973
974inline std::vector<G4double>*
976{
977 return theEnergyOfCrossSectionMax;
978}
979
980//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
981
982inline std::vector<G4TwoPeaksXS*>* G4VEnergyLossProcess::TwoPeaksXS() const
983{
984 return fXSpeaks;
985}
986
987//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
988
989inline std::size_t G4VEnergyLossProcess::NumberOfModels() const
990{
991 return numberOfModels;
992}
993
994//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
995
996inline G4VEmModel* G4VEnergyLossProcess::EmModel(std::size_t index) const
997{
998 return (index < emModels->size()) ? (*emModels)[index] : nullptr;
999}
1000
1001//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1002
1003inline G4VEmModel*
1005{
1006 return modelManager->GetModel((G4int)idx, ver);
1007}
1008
1009//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1010
1011#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
G4VEmModel * SelectModel(G4double energy, std::size_t index)
G4VEmModel * GetModel(G4int idx, G4bool ver=false) const
const G4Material * GetMaterial() const
G4double Energy(const std::size_t index) const
G4double Value(const G4double energy, std::size_t &lastidx) const
Definition: G4Step.hh:62
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:468
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
void SetMaxKinEnergy(G4double e)
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
G4double MaxKinEnergy() const
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 SetLossFluctuations(G4bool val)
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)
G4double MinKinEnergy() const
void SetBaseParticle(const G4ParticleDefinition *p)
G4double CrossSectionBiasingFactor() const
G4bool IsIonisationProcess() 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 *)
void SetIonisation(G4bool val)
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 SetLowestEnergyLimit(G4double)
void SetLambdaTable(G4PhysicsTable *p)
void BuildPhysicsTable(const G4ParticleDefinition &) override
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
void SetMinKinEnergy(G4double e)
const G4ParticleDefinition * SecondaryParticle() const
#define LOG_EKIN_MIN
Definition: templates.hh:98
#define DBL_MAX
Definition: templates.hh:62