Geant4 10.7.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"
60#include "G4UnitsTable.hh"
62#include "G4EmTableType.hh"
63#include "G4PhysicsTable.hh"
64#include "G4PhysicsVector.hh"
65#include "G4EmParameters.hh"
66
67class G4Step;
69class G4VEmModel;
71class G4DataVector;
72class G4Region;
73class G4SafetyHelper;
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80
82{
83public:
84
85 G4VEnergyLossProcess(const G4String& name = "EnergyLoss",
87
88 virtual ~G4VEnergyLossProcess();
89
90private:
91 // clean vectors and arrays
92 void Clean();
93
94 //------------------------------------------------------------------------
95 // Virtual methods to be implemented in concrete processes
96 //------------------------------------------------------------------------
97
98public:
99 virtual G4bool IsApplicable(const G4ParticleDefinition& p) override = 0;
100
101 // obsolete to be removed
102 virtual void PrintInfo() {};
103
104 virtual void ProcessDescription(std::ostream& outFile) const override;
105
106protected:
107
108 virtual void StreamProcessInfo(std::ostream&) const {};
109
111 const G4ParticleDefinition*) = 0;
112
113 //------------------------------------------------------------------------
114 // Methods with standard implementation; may be overwritten if needed
115 //------------------------------------------------------------------------
116
118 const G4Material*, G4double cut);
119
120 //------------------------------------------------------------------------
121 // Virtual methods implementation common to all EM ContinuousDiscrete
122 // processes. Further inheritance is not assumed
123 //------------------------------------------------------------------------
124
125public:
126
127 // prepare all tables
128 virtual void PreparePhysicsTable(const G4ParticleDefinition&) override;
129
130 // build all tables
131 virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
132
133 // build a table
135
136 // build a table
138
139 // Called before tracking of each new G4Track
140 virtual void StartTracking(G4Track*) override;
141
142 // Step limit from AlongStep
144 const G4Track&,
145 G4double previousStepSize,
146 G4double currentMinimumStep,
147 G4double& currentSafety,
148 G4GPILSelection* selection) override;
149
150 // Step limit from cross section
152 const G4Track& track,
153 G4double previousStepSize,
154 G4ForceCondition* condition) override;
155
156 // AlongStep computations
157 virtual G4VParticleChange* AlongStepDoIt(const G4Track&,
158 const G4Step&) override;
159
160 // Sampling of secondaries in vicinity of geometrical boundary
161 // Return sum of secodaries energy
162 G4double SampleSubCutSecondaries(std::vector<G4Track*>&, const G4Step&,
163 G4VEmModel* model, G4int matIdx);
164
165 // PostStep sampling of secondaries
166 virtual G4VParticleChange* PostStepDoIt(const G4Track&,
167 const G4Step&) override;
168
169 // Store all PhysicsTable in files.
170 // Return false in case of any fatal failure at I/O
172 const G4String& directory,
173 G4bool ascii = false) override;
174
175 // Retrieve all Physics from a files.
176 // Return true if all the Physics Table are built.
177 // Return false if any fatal failure.
179 const G4String& directory,
180 G4bool ascii) override;
181
182private:
183
184 // summary printout after initialisation
185 void StreamInfo(std::ostream& out, const G4ParticleDefinition& part,
186 G4bool rst=false) const;
187
188 // store a table
189 G4bool StoreTable(const G4ParticleDefinition* p,
190 G4PhysicsTable*, G4bool ascii,
191 const G4String& directory,
192 const G4String& tname);
193
194 // retrieve a table
195 G4bool RetrieveTable(const G4ParticleDefinition* p,
196 G4PhysicsTable*, G4bool ascii,
197 const G4String& directory,
198 const G4String& tname,
199 G4bool mandatory);
200
201 //------------------------------------------------------------------------
202 // Public interface to cross section, mfp and sampling of fluctuations
203 // These methods are not used in run time
204 //------------------------------------------------------------------------
205
206public:
207
208 // access to dispersion of restricted energy loss
210 const G4DynamicParticle* dp,
211 G4double length);
212
213 // Access to cross section table
215 const G4MaterialCutsCouple* couple);
217 const G4MaterialCutsCouple* couple,
218 G4double logKineticEnergy);
219
220 // access to cross section
221 G4double MeanFreePath(const G4Track& track);
222
223 // access to step limit
225 G4double previousStepSize,
226 G4double currentMinimumStep,
227 G4double& currentSafety);
228
229protected:
230
231 // implementation of the pure virtual method
232 virtual G4double GetMeanFreePath(const G4Track& track,
233 G4double previousStepSize,
234 G4ForceCondition* condition) override;
235
236 // implementation of the pure virtual method
237 virtual G4double GetContinuousStepLimit(const G4Track& track,
238 G4double previousStepSize,
239 G4double currentMinimumStep,
240 G4double& currentSafety) override;
241
242 //------------------------------------------------------------------------
243 // Run time method which may be also used by derived processes
244 //------------------------------------------------------------------------
245
246 // creeation of an empty vector for cross section
248 G4double cut);
249
250 inline size_t CurrentMaterialCutsCoupleIndex() const;
251
252 //------------------------------------------------------------------------
253 // Specific methods to set, access, modify models
254 //------------------------------------------------------------------------
255
256 // Select model in run time
257 inline void SelectModel(G4double kinEnergy);
258
259public:
260 // Select model by energy and region index
261 inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
262 size_t& idx) const;
263
264 // Add EM model coupled with fluctuation model for region, smaller value
265 // of order defines which pair of models will be selected for a given
266 // energy interval
268 G4VEmFluctuationModel* fluc = 0,
269 const G4Region* region = nullptr);
270
271 // Define new energy range for the model identified by the name
273
274 // Assign a model to a process local list, to enable the list in run time
275 // the derived process should execute AddEmModel(..) for all such models
276 void SetEmModel(G4VEmModel*, G4int index=0);
277
278 // return a model from the local list
279 G4VEmModel* EmModel(size_t index=0) const;
280
281 // Access to models
282 G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
283
284 G4int NumberOfModels() const;
285
286 // Assign a fluctuation model to a process
288
289 // return the assigned fluctuation model
291
292 //------------------------------------------------------------------------
293 // Define and access particle type
294 //------------------------------------------------------------------------
295
296protected:
297 inline void SetParticle(const G4ParticleDefinition* p);
298 inline void SetSecondaryParticle(const G4ParticleDefinition* p);
299
300public:
301 inline void SetBaseParticle(const G4ParticleDefinition* p);
302 inline const G4ParticleDefinition* Particle() const;
303 inline const G4ParticleDefinition* BaseParticle() const;
304 inline const G4ParticleDefinition* SecondaryParticle() const;
305
306 //------------------------------------------------------------------------
307 // Get/set parameters to configure the process at initialisation time
308 //------------------------------------------------------------------------
309
310 // Add subcutoff option for the region
311 void ActivateSubCutoff(G4bool val, const G4Region* region = nullptr);
312
313 // Activate biasing
314 void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
315
317 const G4String& region,
318 G4bool flag = true);
319
320 void ActivateSecondaryBiasing(const G4String& region, G4double factor,
321 G4double energyLimit);
322
323 // Add subcutoff process (bremsstrahlung) to sample secondary
324 // particle production in vicinity of the geometry boundary
326
327 inline void SetLossFluctuations(G4bool val);
328
329 inline void SetIntegral(G4bool val);
330 inline G4bool IsIntegral() const;
331
332 // Set/Get flag "isIonisation"
333 void SetIonisation(G4bool val);
334 inline G4bool IsIonisationProcess() const;
335
336 // Redefine parameteters for stepping control
338 void SetStepFunction(G4double v1, G4double v2);
340
341 inline G4int NumberOfSubCutoffRegions() const;
342
343 //------------------------------------------------------------------------
344 // Specific methods to path Physics Tables to the process
345 //------------------------------------------------------------------------
346
348 void SetCSDARangeTable(G4PhysicsTable* pRange);
352
355
356 // Binning for dEdx, range, inverse range and labda tables
357 void SetDEDXBinning(G4int nbins);
358
359 // Min kinetic energy for tables
361 inline G4double MinKinEnergy() const;
362
363 // Max kinetic energy for tables
365 inline G4double MaxKinEnergy() const;
366
367 // Biasing parameters
368 inline G4double CrossSectionBiasingFactor() const;
369
370 // Return values for given G4MaterialCutsCouple
371 inline G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple*);
372 inline G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple*,
373 G4double logKineticEnergy);
374 inline G4double GetDEDXForSubsec(G4double kineticEnergy,
375 const G4MaterialCutsCouple*);
376 inline G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple*);
377 inline G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple*,
378 G4double logKineticEnergy);
379 inline G4double GetCSDARange(G4double kineticEnergy,
380 const G4MaterialCutsCouple*);
381 inline G4double GetRangeForLoss(G4double kineticEnergy,
382 const G4MaterialCutsCouple*);
383 inline G4double GetRangeForLoss(G4double kineticEnergy,
385 G4double logKineticEnergy);
386 inline G4double GetKineticEnergy(G4double range,
387 const G4MaterialCutsCouple*);
388 inline G4double GetLambda(G4double kineticEnergy,const G4MaterialCutsCouple*);
389 inline G4double GetLambda(G4double kineticEnergy,const G4MaterialCutsCouple*,
390 G4double logKineticEnergy);
391
392 inline G4bool TablesAreBuilt() const;
393
394 // Access to specific tables
395 inline G4PhysicsTable* DEDXTable() const;
396 inline G4PhysicsTable* DEDXTableForSubsec() const;
398 inline G4PhysicsTable* IonisationTable() const;
400 inline G4PhysicsTable* CSDARangeTable() const;
401 inline G4PhysicsTable* SecondaryRangeTable() const;
402 inline G4PhysicsTable* RangeTableForLoss() const;
403 inline G4PhysicsTable* InverseRangeTable() const;
404 inline G4PhysicsTable* LambdaTable() const;
405 inline G4PhysicsTable* SubLambdaTable() const;
406
407 //------------------------------------------------------------------------
408 // Run time method for simulation of ionisation
409 //------------------------------------------------------------------------
410
411 // access atom on which interaction happens
412 const G4Element* GetCurrentElement() const;
413
414 // Set scaling parameters for ions is needed to G4EmCalculator
415 inline void SetDynamicMassCharge(G4double massratio, G4double charge2ratio);
416
417private:
418
419 void FillSecondariesAlongStep(G4double& eloss, G4double& weight);
420
421 void PrintWarning(G4String, G4double val);
422
423 // define material and indexes
424 inline void DefineMaterial(const G4MaterialCutsCouple* couple);
425
426 //------------------------------------------------------------------------
427 // Compute values using scaling relation, mass and charge of based particle
428 //------------------------------------------------------------------------
429 inline G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy);
430 inline G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy,
431 G4double logScaledKinEnergy);
432 inline G4double GetSubDEDXForScaledEnergy(G4double scaledKinEnergy);
433 inline G4double GetIonisationForScaledEnergy(G4double scaledKinEnergy);
434 inline G4double GetSubIonisationForScaledEnergy(G4double scaledKinEnergy);
435 inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy);
436 inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy,
437 G4double logScaledKinEnergy);
438 inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinEnergy);
439 inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinEnergy,
440 G4double logScaledKinEnergy);
441 inline G4double ScaledKinEnergyForLoss(G4double range);
442 inline G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy);
443 inline G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy,
444 G4double logScaledKinEnergy);
445 void ComputeLambdaForScaledEnergy(G4double scaledKinEnergy,
446 G4double logScaledKinEnergy);
447
448 // hide assignment operator
450 G4VEnergyLossProcess & operator=(const G4VEnergyLossProcess &right) = delete;
451
452 // ======== Parameters of the class fixed at construction =========
453
454 G4LossTableManager* lManager;
455 G4EmModelManager* modelManager;
456 G4EmBiasingManager* biasManager;
457 G4SafetyHelper* safetyHelper;
458 G4EmParameters* theParameters;
459
460 const G4ParticleDefinition* secondaryParticle;
461 const G4ParticleDefinition* theElectron;
462 const G4ParticleDefinition* thePositron;
463 const G4ParticleDefinition* theGamma;
464 const G4ParticleDefinition* theGenericIon;
465
466 // ======== Parameters of the class fixed at initialisation =======
467
468 std::vector<G4VEmModel*> emModels;
469 G4VEmFluctuationModel* fluctModel;
470 G4VAtomDeexcitation* atomDeexcitation;
471 G4VSubCutProducer* subcutProducer;
472 std::vector<const G4Region*> scoffRegions;
473 G4int nSCoffRegions;
474 G4bool* idxSCoffRegions;
475
476 std::vector<G4VEnergyLossProcess*> scProcesses;
477 G4int nProcesses;
478
479 // tables and vectors
480 G4PhysicsTable* theDEDXTable;
481 G4PhysicsTable* theDEDXSubTable;
482 G4PhysicsTable* theDEDXunRestrictedTable;
483 G4PhysicsTable* theIonisationTable;
484 G4PhysicsTable* theIonisationSubTable;
485 G4PhysicsTable* theRangeTableForLoss;
486 G4PhysicsTable* theCSDARangeTable;
487 G4PhysicsTable* theSecondaryRangeTable;
488 G4PhysicsTable* theInverseRangeTable;
489 G4PhysicsTable* theLambdaTable;
490 G4PhysicsTable* theSubLambdaTable;
491
492 size_t idxDEDX = 0;
493 size_t idxDEDXSub = 0;
494 size_t idxDEDXunRestricted = 0;
495 size_t idxIonisation = 0;
496 size_t idxIonisationSub = 0;
497 size_t idxRange = 0;
498 size_t idxCSDA = 0;
499 size_t idxSecRange = 0;
500 size_t idxInverseRange = 0;
501 size_t idxLambda = 0;
502 size_t idxSubLambda = 0;
503
504 std::vector<G4double> theDEDXAtMaxEnergy;
505 std::vector<G4double> theRangeAtMaxEnergy;
506 std::vector<G4double> theEnergyOfCrossSectionMax;
507 std::vector<G4double> theCrossSectionMax;
508
509 const std::vector<G4double>* theDensityFactor;
510 const std::vector<G4int>* theDensityIdx;
511
512 const G4DataVector* theCuts;
513 const G4DataVector* theSubCuts;
514
515 const G4ParticleDefinition* baseParticle;
516
517 G4int nBins;
518 G4int nBinsCSDA;
519
520 G4double lowestKinEnergy;
521 G4double minKinEnergy;
522 G4double maxKinEnergy;
523 G4double maxKinEnergyCSDA;
524
525 G4double linLossLimit;
526 G4double dRoverRange;
527 G4double finalRange;
528 G4double lambdaFactor;
529 G4double logLambdafactor;
530 G4double biasFactor;
531
532 G4bool lossFluctuationFlag;
533 G4bool rndmStepFlag;
534 G4bool tablesAreBuilt;
535 G4bool integral;
536 G4bool isIon;
537 G4bool isIonisation;
538 G4bool useSubCutoff;
539 G4bool useDeexcitation;
540 G4bool biasFlag;
541 G4bool weightFlag;
542 G4bool isMaster;
543 G4bool actIntegral;
544 G4bool actLinLossLimit;
545 G4bool actLossFluc;
546 G4bool actBinning;
547 G4bool actMinKinEnergy;
548 G4bool actMaxKinEnergy;
549
550protected:
551
556
566
567 // ======== Cached values - may be state dependent ================
568
569private:
570
571 std::vector<G4DynamicParticle*> secParticles;
572 std::vector<G4Track*> scTracks;
573
574 const G4ParticleDefinition* particle;
575
576 G4VEmModel* currentModel;
577 size_t basedCoupleIndex;
578 size_t lastIdx;
579
580 G4double massRatio;
581 G4double logMassRatio;
582 G4double fFactor;
583 G4double reduceFactor;
584 G4double chargeSqRatio;
585
586 G4GPILSelection aGPILSelection;
587
588 G4int secID;
589 G4int subsecID;
590 G4int biasID;
591};
592
593// ======== Run time inline methods ================
594
596{
597 return currentCoupleIndex;
598}
599
600//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
601
603{
604 currentModel = modelManager->SelectModel(kinEnergy, currentCoupleIndex);
605 currentModel->SetCurrentCouple(currentCouple);
606}
607
608//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
609
611 G4double kinEnergy, size_t& idx) const
612{
613 return modelManager->SelectModel(kinEnergy, idx);
614}
615
616//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
617
618inline void
619G4VEnergyLossProcess::DefineMaterial(const G4MaterialCutsCouple* couple)
620{
621 if(couple != currentCouple) {
622 currentCouple = couple;
623 currentMaterial = couple->GetMaterial();
624 currentCoupleIndex = couple->GetIndex();
625 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
626 fFactor = chargeSqRatio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
627 reduceFactor = 1.0/(fFactor*massRatio);
629 idxLambda = idxSubLambda = 0;
630 }
631}
632
633//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
634
636 G4double charge2ratio)
637{
638 massRatio = massratio;
639 logMassRatio = G4Log(massRatio);
640 fFactor = charge2ratio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
641 chargeSqRatio = charge2ratio;
642 reduceFactor = 1.0/(fFactor*massRatio);
643}
644
645//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
646
647inline G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e)
648{
649 /*
650 G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
651 << basedCoupleIndex << " E(MeV)= " << e
652 << " Emin= " << minKinEnergy << " Factor= " << fFactor
653 << " " << theDEDXTable << G4endl; */
654 G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->Value(e, idxDEDX);
655 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
656 return x;
657}
658
659inline
660G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e, G4double loge)
661{
662 /*
663 G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= "
664 << basedCoupleIndex << " E(MeV)= " << e
665 << " Emin= " << minKinEnergy << " Factor= " << fFactor
666 << " " << theDEDXTable << G4endl; */
667 G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->LogVectorValue(e,loge);
668 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
669 return x;
670}
671
672//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
673
674inline G4double G4VEnergyLossProcess::GetSubDEDXForScaledEnergy(G4double e)
675{
676 G4double x =
677 fFactor*(*theDEDXSubTable)[basedCoupleIndex]->Value(e, idxDEDXSub);
678 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
679 return x;
680}
681
682//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
683
684inline G4double G4VEnergyLossProcess::GetIonisationForScaledEnergy(G4double e)
685{
686 G4double x =
687 fFactor*(*theIonisationTable)[basedCoupleIndex]->Value(e, idxIonisation);
688 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
689 return x;
690}
691
692//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
693
694inline
695G4double G4VEnergyLossProcess::GetSubIonisationForScaledEnergy(G4double e)
696{
697 G4double x = fFactor*
698 (*theIonisationSubTable)[basedCoupleIndex]->Value(e, idxIonisationSub);
699 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
700 return x;
701}
702
703//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
704
705inline G4double G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e)
706{
707 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
708 // << basedCoupleIndex << " E(MeV)= " << e
709 // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
710 if(basedCoupleIndex != lastIdx || preStepRangeEnergy != e) {
711 lastIdx = basedCoupleIndex;
714 ((*theRangeTableForLoss)[basedCoupleIndex])->Value(e, idxRange);
715 if(e < minKinEnergy) { computedRange *= std::sqrt(e/minKinEnergy); }
716 }
717 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
718 // << basedCoupleIndex << " E(MeV)= " << e
719 // << " R= " << fRange << " " << theRangeTableForLoss << G4endl;
720
721 return computedRange;
722}
723
724inline G4double
725G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e, G4double loge)
726{
727 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
728 // << basedCoupleIndex << " E(MeV)= " << e
729 // << " lastIdx= " << lastIdx << " " << theRangeTableForLoss << G4endl;
730 if(basedCoupleIndex != lastIdx || preStepRangeEnergy != e) {
731 lastIdx = basedCoupleIndex;
734 ((*theRangeTableForLoss)[basedCoupleIndex])->LogVectorValue(e, loge);
735 if(e < minKinEnergy) { computedRange *= std::sqrt(e/minKinEnergy); }
736 }
737 //G4cout << "G4VEnergyLossProcess::GetScaledRange: Idx= "
738 // << basedCoupleIndex << " E(MeV)= " << e
739 // << " R= " << fRange << " " << theRangeTableForLoss << G4endl;
740
741 return computedRange;
742}
743
744//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
745
746inline G4double
747G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e)
748{
749 G4double x;
750 if (e < maxKinEnergyCSDA) {
751 x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e, idxCSDA);
752 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
753 } else {
754 x = theRangeAtMaxEnergy[basedCoupleIndex] +
755 (e - maxKinEnergyCSDA)/theDEDXAtMaxEnergy[basedCoupleIndex];
756 }
757 return x;
758}
759
760inline G4double
761G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e,
762 G4double loge)
763{
764 G4double x;
765 if (e < maxKinEnergyCSDA) {
766 x = ((*theCSDARangeTable)[basedCoupleIndex])->LogVectorValue(e, loge);
767 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
768 } else {
769 x = theRangeAtMaxEnergy[basedCoupleIndex] +
770 (e - maxKinEnergyCSDA)/theDEDXAtMaxEnergy[basedCoupleIndex];
771 }
772 return x;
773}
774
775//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
776
777inline G4double G4VEnergyLossProcess::ScaledKinEnergyForLoss(G4double r)
778{
779 //G4cout << "G4VEnergyLossProcess::GetEnergy: Idx= "
780 // << basedCoupleIndex << " R(mm)= " << r << " "
781 // << theInverseRangeTable << G4endl;
782 G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex];
783 G4double rmin = v->Energy(0);
784 G4double e = 0.0;
785 if(r >= rmin) { e = v->Value(r, idxInverseRange); }
786 else if(r > 0.0) {
787 G4double x = r/rmin;
788 e = minKinEnergy*x*x;
789 }
790 return e;
791}
792
793//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
794
795inline G4double G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e)
796{
797 return fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
798}
799
800inline G4double
801G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e, G4double loge)
802{
803 return fFactor*((*theLambdaTable)[basedCoupleIndex])->LogVectorValue(e,loge);
804}
805
806//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
807
808inline G4double
810 const G4MaterialCutsCouple* couple)
811{
812 DefineMaterial(couple);
813 return GetDEDXForScaledEnergy(kinEnergy*massRatio);
814}
815
816inline G4double
818 const G4MaterialCutsCouple* couple,
819 G4double logKinEnergy)
820{
821 DefineMaterial(couple);
822 return GetDEDXForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
823}
824
825//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
826
827inline G4double
829 const G4MaterialCutsCouple* couple)
830{
831 DefineMaterial(couple);
832 return GetSubDEDXForScaledEnergy(kineticEnergy*massRatio);
833}
834
835//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
836
837inline G4double
839 const G4MaterialCutsCouple* couple)
840{
841 G4double x = fRange;
842 DefineMaterial(couple);
843 if(theCSDARangeTable) {
844 x = reduceFactor * GetLimitScaledRangeForScaledEnergy(kinEnergy*massRatio);
845 } else if(theRangeTableForLoss) {
846 x = reduceFactor * GetScaledRangeForScaledEnergy(kinEnergy*massRatio);
847 }
848 return x;
849}
850
851inline G4double
853 const G4MaterialCutsCouple* couple,
854 G4double logKinEnergy)
855{
856 G4double x = fRange;
857 DefineMaterial(couple);
858 if(theCSDARangeTable) {
859 x = reduceFactor * GetLimitScaledRangeForScaledEnergy(kinEnergy*massRatio,
860 logKinEnergy+logMassRatio);
861 } else if(theRangeTableForLoss) {
862 x = reduceFactor * GetScaledRangeForScaledEnergy(kinEnergy*massRatio,
863 logKinEnergy+logMassRatio);
864 }
865 return x;
866}
867
868//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
869
870inline G4double
872 const G4MaterialCutsCouple* couple)
873{
874 DefineMaterial(couple);
875 return (theCSDARangeTable) ?
876 GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor
877 : DBL_MAX;
878}
879
880//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
881
882inline G4double
884 const G4MaterialCutsCouple* couple)
885{
886 // G4cout << "GetRangeForLoss: Range from " << GetProcessName() << G4endl;
887 DefineMaterial(couple);
888 return reduceFactor * GetScaledRangeForScaledEnergy(kinEnergy*massRatio);
889}
890
891
892inline G4double
894 const G4MaterialCutsCouple* couple,
895 G4double logKinEnergy)
896{
897 // G4cout << "GetRangeForLoss: Range from " << GetProcessName() << G4endl;
898 DefineMaterial(couple);
899 return reduceFactor * GetScaledRangeForScaledEnergy(kinEnergy*massRatio,
900 logKinEnergy+logMassRatio);
901}
902
903//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
904
905inline G4double
907 const G4MaterialCutsCouple* couple)
908{
909 DefineMaterial(couple);
910 return ScaledKinEnergyForLoss(range/reduceFactor)/massRatio;
911}
912
913//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
914
915inline G4double
917 const G4MaterialCutsCouple* couple)
918{
919 DefineMaterial(couple);
920 return theLambdaTable ? GetLambdaForScaledEnergy(kinEnergy*massRatio) : 0.0;
921}
922
923inline G4double
925 const G4MaterialCutsCouple* couple,
926 G4double logKinEnergy)
927{
928 DefineMaterial(couple);
929 return theLambdaTable
930 ? GetLambdaForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio)
931 : 0.0;
932}
933
934// ======== Get/Set inline methods used at initialisation ================
935
937{
938 fluctModel = p;
939}
940
941//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
942
944{
945 return fluctModel;
946}
947
948//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
949
951{
952 particle = p;
953}
954
955//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
956
957inline void
959{
960 secondaryParticle = p;
961}
962
963//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
964
965inline void
967{
968 baseParticle = p;
969}
970
971//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
972
974{
975 return particle;
976}
977
978//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
979
981{
982 return baseParticle;
983}
984
985//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
986
987inline const G4ParticleDefinition*
989{
990 return secondaryParticle;
991}
992
993//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
994
996{
997 lossFluctuationFlag = val;
998 actLossFluc = true;
999}
1000
1001//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1002
1004{
1005 integral = val;
1006 actIntegral = true;
1007}
1008
1009//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1010
1012{
1013 return integral;
1014}
1015
1016//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1017
1019{
1020 return isIonisation;
1021}
1022
1023//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1024
1026{
1027 return nSCoffRegions;
1028}
1029
1030//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1031
1033{
1034 return minKinEnergy;
1035}
1036
1037//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1038
1040{
1041 return maxKinEnergy;
1042}
1043
1044//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1045
1047{
1048 return biasFactor;
1049}
1050
1051//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1052
1054{
1055 return tablesAreBuilt;
1056}
1057
1058//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1059
1061{
1062 return theDEDXTable;
1063}
1064
1065//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1066
1068{
1069 return theDEDXSubTable;
1070}
1071
1072//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1073
1075{
1076 return theDEDXunRestrictedTable;
1077}
1078
1079//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1080
1082{
1083 return theIonisationTable;
1084}
1085
1086//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1087
1089{
1090 return theIonisationSubTable;
1091}
1092
1093//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1094
1096{
1097 return theCSDARangeTable;
1098}
1099
1100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1101
1103{
1104 return theSecondaryRangeTable;
1105}
1106
1107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1108
1110{
1111 return theRangeTableForLoss;
1112}
1113
1114//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1115
1117{
1118 return theInverseRangeTable;
1119}
1120
1121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1122
1124{
1125 return theLambdaTable;
1126}
1127
1128//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1129
1131{
1132 return theSubLambdaTable;
1133}
1134
1135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1136
1137#endif
G4EmTableType
@ fRestricted
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4GPILSelection
G4double G4Log(G4double x)
Definition: G4Log.hh:226
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4VEmModel * SelectModel(G4double &energy, size_t &index)
const G4Material * GetMaterial() const
G4double Energy(std::size_t index) const
G4double Value(G4double theEnergy, std::size_t &lastidx) const
Definition: G4Step.hh:62
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:465
const G4ParticleDefinition * BaseParticle() const
G4PhysicsTable * RangeTableForLoss() const
G4double GetRangeForLoss(G4double kineticEnergy, const G4MaterialCutsCouple *)
virtual G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
void SetMaxKinEnergy(G4double e)
G4ParticleChangeForLoss fParticleChange
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
G4PhysicsTable * InverseRangeTable() const
void ActivateSubCutoff(G4bool val, const G4Region *region=nullptr)
G4double MeanFreePath(const G4Track &track)
void SetFluctModel(G4VEmFluctuationModel *)
G4double GetKineticEnergy(G4double range, const G4MaterialCutsCouple *)
G4PhysicsTable * CSDARangeTable() const
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) 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)
virtual void ProcessDescription(std::ostream &outFile) const override
void UpdateEmModel(const G4String &, G4double, G4double)
G4PhysicsTable * BuildDEDXTable(G4EmTableType tType=fRestricted)
const G4MaterialCutsCouple * currentCouple
void AddCollaborativeProcess(G4VEnergyLossProcess *)
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idx) const
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *, G4double cut)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
const G4Element * GetCurrentElement() const
void SetDEDXBinning(G4int nbins)
void SetStepFunction(G4double v1, G4double v2)
G4PhysicsTable * DEDXTableForSubsec() const
G4double GetLambda(G4double kineticEnergy, const G4MaterialCutsCouple *)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
G4PhysicsTable * IonisationTableForSubsec() const
const G4Material * currentMaterial
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
void SetParticle(const G4ParticleDefinition *p)
G4VEmModel * EmModel(size_t index=0) const
G4PhysicsTable * SecondaryRangeTable() const
void SetLossFluctuations(G4bool val)
const G4ParticleDefinition * Particle() const
void SetInverseRangeTable(G4PhysicsTable *p)
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
G4double SampleSubCutSecondaries(std::vector< G4Track * > &, const G4Step &, G4VEmModel *model, G4int matIdx)
void ActivateForcedInteraction(G4double length, const G4String &region, G4bool flag=true)
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
G4PhysicsTable * SubLambdaTable() const
void SetBaseParticle(const G4ParticleDefinition *p)
G4double GetDEDXForSubsec(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double CrossSectionBiasingFactor() const
G4bool IsIonisationProcess() const
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple)
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
virtual G4bool IsApplicable(const G4ParticleDefinition &p) override=0
void SetSecondaryRangeTable(G4PhysicsTable *p)
G4PhysicsTable * BuildLambdaTable(G4EmTableType tType=fRestricted)
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4double GetCSDARange(G4double kineticEnergy, const G4MaterialCutsCouple *)
void SetIonisation(G4bool val)
void SetSubLambdaTable(G4PhysicsTable *p)
G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple *)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
void SetLinearLossLimit(G4double val)
void SetLowestEnergyLimit(G4double)
void SetLambdaTable(G4PhysicsTable *p)
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
size_t CurrentMaterialCutsCoupleIndex() const
G4PhysicsTable * IonisationTable() 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
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
G4VEmFluctuationModel * FluctModel()
virtual void StartTracking(G4Track *) override
void SetMinKinEnergy(G4double e)
const G4ParticleDefinition * SecondaryParticle() const
#define DBL_MAX
Definition: templates.hh:62