Geant4 11.2.2
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 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 inline G4double LogScaledEkin(const G4Track& aTrack);
407
408 void ComputeLambdaForScaledEnergy(G4double scaledKinE,
409 const G4Track& aTrack);
410
411 G4bool IsRegionForCubcutProcessor(const G4Track& aTrack);
412
413protected:
414
416 const G4Material* currentMaterial = nullptr;
418
419private:
420
421 G4LossTableManager* lManager;
422 G4EmModelManager* modelManager;
423 G4VEmModel* currentModel = nullptr;
424 G4EmBiasingManager* biasManager = nullptr;
425 G4SafetyHelper* safetyHelper;
426 G4EmParameters* theParameters;
427 G4VEmFluctuationModel* fluctModel = nullptr;
428 G4VAtomDeexcitation* atomDeexcitation = nullptr;
429 G4VSubCutProducer* subcutProducer = nullptr;
430
431 const G4ParticleDefinition* particle = nullptr;
432 const G4ParticleDefinition* baseParticle = nullptr;
433 const G4ParticleDefinition* secondaryParticle = nullptr;
434 G4EmDataHandler* theData = nullptr;
435
436 G4PhysicsTable* theDEDXTable = nullptr;
437 G4PhysicsTable* theDEDXunRestrictedTable = nullptr;
438 G4PhysicsTable* theIonisationTable = nullptr;
439 G4PhysicsTable* theRangeTableForLoss = nullptr;
440 G4PhysicsTable* theCSDARangeTable = nullptr;
441 G4PhysicsTable* theInverseRangeTable = nullptr;
442 G4PhysicsTable* theLambdaTable = nullptr;
443
444 std::vector<const G4Region*>* scoffRegions = nullptr;
445 std::vector<G4VEmModel*>* emModels = nullptr;
446 const std::vector<G4int>* theDensityIdx = nullptr;
447 const std::vector<G4double>* theDensityFactor = nullptr;
448 const G4DataVector* theCuts = nullptr;
449
450 std::vector<G4double>* theEnergyOfCrossSectionMax = nullptr;
451 std::vector<G4TwoPeaksXS*>* fXSpeaks = nullptr;
452
453 G4double lowestKinEnergy;
454 G4double minKinEnergy;
455 G4double maxKinEnergy;
456 G4double maxKinEnergyCSDA;
457
458 G4double linLossLimit = 0.01;
459 G4double dRoverRange = 0.2;
460 G4double finalRange;
461 G4double lambdaFactor = 0.8;
462 G4double invLambdaFactor;
463 G4double biasFactor = 1.0;
464
465 G4double massRatio = 1.0;
466 G4double logMassRatio = 0.0;
467 G4double fFactor = 1.0;
468 G4double reduceFactor = 1.0;
469 G4double chargeSqRatio = 1.0;
470 G4double fRange = 0.0;
471 G4double fRangeEnergy = 0.0;
472
473protected:
474
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 useCutAsFinalRange = 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
585//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
586
587inline
588G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e, G4double loge)
589{
590 /*
591 G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
592 << basedCoupleIndex << " E(MeV)= " << e
593 << " Emin= " << minKinEnergy << " Factor= " << fFactor
594 << " " << theDEDXTable << G4endl; */
595 G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->LogVectorValue(e,loge);
596 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
597 return x;
598}
599
600//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
601
602inline G4double G4VEnergyLossProcess::GetIonisationForScaledEnergy(G4double e)
603{
604 G4double x =
605 fFactor*(*theIonisationTable)[basedCoupleIndex]->Value(e, idxIonisation);
606 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
607 return x;
608}
609
610//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
611
612inline G4double G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e)
613{
614 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
615 // << basedCoupleIndex << " E(MeV)= " << e
616 // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
617 if(currentCoupleIndex != coupleIdxRange || fRangeEnergy != e) {
618 coupleIdxRange = currentCoupleIndex;
619 fRangeEnergy = e;
620 fRange = reduceFactor*((*theRangeTableForLoss)[basedCoupleIndex])->Value(e, idxRange);
621 if (fRange < 0.0) { fRange = 0.0; }
622 else if (e < minKinEnergy) { fRange *= std::sqrt(e/minKinEnergy); }
623 }
624 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
625 // << basedCoupleIndex << " E(MeV)= " << e
626 // << " R= " << computedRange << " " << theRangeTableForLoss << G4endl;
627 return fRange;
628}
629
630inline G4double
631G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e, G4double loge)
632{
633 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
634 // << basedCoupleIndex << " E(MeV)= " << e
635 // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
636 if(currentCoupleIndex != coupleIdxRange || fRangeEnergy != e) {
637 coupleIdxRange = currentCoupleIndex;
638 fRangeEnergy = e;
639 fRange = reduceFactor*((*theRangeTableForLoss)[basedCoupleIndex])->LogVectorValue(e, loge);
640 if (fRange < 0.0) { fRange = 0.0; }
641 else if (e < minKinEnergy) { fRange *= std::sqrt(e/minKinEnergy); }
642 }
643 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
644 // << basedCoupleIndex << " E(MeV)= " << e
645 // << " R= " << fRange << " " << theRangeTableForLoss << G4endl;
646 return fRange;
647}
648
649//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
650
651inline G4double
652G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e)
653{
654 G4double x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e, idxCSDA);
655 if (x < 0.0) { x = 0.0; }
656 else if (e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
657 return x;
658}
659
660//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
661
662inline G4double
663G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e,
664 G4double loge)
665{
666 G4double x = ((*theCSDARangeTable)[basedCoupleIndex])->LogVectorValue(e, loge);
667 if (x < 0.0) { x = 0.0; }
668 else if (e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
669 return x;
670}
671
672//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
673
674inline G4double G4VEnergyLossProcess::ScaledKinEnergyForLoss(G4double r)
675{
676 //G4cout << "G4VEnergyLossProcess::GetEnergy: Idx= "
677 // << basedCoupleIndex << " R(mm)= " << r << " "
678 // << theInverseRangeTable << G4endl;
679 G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex];
680 G4double rmin = v->Energy(0);
681 G4double e = 0.0;
682 if(r >= rmin) { e = v->Value(r, idxInverseRange); }
683 else if(r > 0.0) {
684 G4double x = r/rmin;
685 e = minKinEnergy*x*x;
686 }
687 return e;
688}
689
690//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
691
692inline G4double G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e)
693{
694 return fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
695}
696
697//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
698
699inline G4double
700G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e, G4double loge)
701{
702 return fFactor*((*theLambdaTable)[basedCoupleIndex])->LogVectorValue(e, loge);
703}
704
705//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
706
707inline G4double G4VEnergyLossProcess::LogScaledEkin(const G4Track& track)
708{
709 return track.GetDynamicParticle()->GetLogKineticEnergy() + logMassRatio;
710}
711
712//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
713
714inline G4double
716 const G4MaterialCutsCouple* couple)
717{
718 DefineMaterial(couple);
719 return GetDEDXForScaledEnergy(kinEnergy*massRatio);
720}
721
722//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
723
724inline G4double
726 const G4MaterialCutsCouple* couple,
727 G4double logKinEnergy)
728{
729 DefineMaterial(couple);
730 return GetDEDXForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
731}
732
733//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
734
735inline G4double
737 const G4MaterialCutsCouple* couple)
738{
739 DefineMaterial(couple);
740 return GetScaledRangeForScaledEnergy(kinEnergy*massRatio);
741}
742
743//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
744
745inline G4double
747 const G4MaterialCutsCouple* couple,
748 G4double logKinEnergy)
749{
750 DefineMaterial(couple);
751 return GetScaledRangeForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
752}
753
754//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
755
756inline G4double
758 const G4MaterialCutsCouple* couple)
759{
760 DefineMaterial(couple);
761 return (nullptr == theCSDARangeTable) ? DBL_MAX :
762 GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
763}
764
765//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
766
767inline G4double
769 const G4MaterialCutsCouple* couple)
770{
771 DefineMaterial(couple);
772 return ScaledKinEnergyForLoss(range/reduceFactor)/massRatio;
773}
774
775//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
776
777inline G4double
779 const G4MaterialCutsCouple* couple)
780{
781 DefineMaterial(couple);
782 return (nullptr != theLambdaTable) ?
783 GetLambdaForScaledEnergy(kinEnergy*massRatio) : 0.0;
784}
785
786//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
787
788inline G4double
790 const G4MaterialCutsCouple* couple,
791 G4double logKinEnergy)
792{
793 DefineMaterial(couple);
794 return (nullptr != theLambdaTable) ?
795 GetLambdaForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio)
796 : 0.0;
797}
798
799// ======== Get/Set inline methods used at initialisation ================
800
802{
803 fluctModel = p;
804}
805
806//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
807
809{
810 return fluctModel;
811}
812
813//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
814
816{
817 particle = p;
818}
819
820//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
821
822inline void
824{
825 secondaryParticle = p;
826}
827
828//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
829
830inline void
832{
833 baseParticle = p;
834}
835
836//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
837
839{
840 return particle;
841}
842
843//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
844
846{
847 return baseParticle;
848}
849
850//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
851
852inline const G4ParticleDefinition*
854{
855 return secondaryParticle;
856}
857
858//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
859
861{
862 lossFluctuationFlag = val;
863 actLossFluc = true;
864}
865
866//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
867
869{
870 spline = val;
871}
872
873//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
874
876{
877 fXSType = val;
878}
879
880//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
881
883{
884 return fXSType;
885}
886
887//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
888
890{
891 return isIonisation;
892}
893
894//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
895
897{
898 return nSCoffRegions;
899}
900
901//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
902
904{
905 return minKinEnergy;
906}
907
908//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
909
911{
912 return maxKinEnergy;
913}
914
915//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
916
918{
919 return biasFactor;
920}
921
922//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
923
925{
926 return tablesAreBuilt;
927}
928
929//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
930
932{
933 return theDEDXTable;
934}
935
936//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
937
939{
940 return theDEDXunRestrictedTable;
941}
942
943//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
944
946{
947 return theIonisationTable;
948}
949
950//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
951
953{
954 return theCSDARangeTable;
955}
956
957//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
958
960{
961 return theRangeTableForLoss;
962}
963
964//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
965
967{
968 return theInverseRangeTable;
969}
970
971//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
972
974{
975 return theLambdaTable;
976}
977
978//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
979
981{
982 return baseMat;
983}
984
985//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
986
987inline std::vector<G4double>*
989{
990 return theEnergyOfCrossSectionMax;
991}
992
993//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
994
995inline std::vector<G4TwoPeaksXS*>* G4VEnergyLossProcess::TwoPeaksXS() const
996{
997 return fXSpeaks;
998}
999
1000//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1001
1003{
1004 return numberOfModels;
1005}
1006
1007//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1008
1009inline G4VEmModel* G4VEnergyLossProcess::EmModel(std::size_t index) const
1010{
1011 return (index < emModels->size()) ? (*emModels)[index] : nullptr;
1012}
1013
1014//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1015
1016inline G4VEmModel*
1018{
1019 return modelManager->GetModel((G4int)idx, ver);
1020}
1021
1022//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1023
1024#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
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
const G4DynamicParticle * GetDynamicParticle() const
void SetCurrentCouple(const G4MaterialCutsCouple *)
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 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)
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