source: trunk/source/processes/electromagnetic/utils/include/G4VEnergyLossProcess.hh@ 1347

Last change on this file since 1347 was 1340, checked in by garnier, 15 years ago

update ti head

File size: 36.1 KB
RevLine 
[819]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//
[1340]26// $Id: G4VEnergyLossProcess.hh,v 1.93 2010/10/14 16:27:35 vnivanch Exp $
[819]27// GEANT4 tag $Name:
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4VEnergyLossProcess
35//
36// Author: Vladimir Ivanchenko on base of Laszlo Urban code
37//
38// Creation date: 03.01.2002
39//
40// Modifications:
41//
42// 26-12-02 Secondary production moved to derived classes (V.Ivanchenko)
43// 20-01-03 Migrade to cut per region (V.Ivanchenko)
44// 24-01-03 Make models region aware (V.Ivanchenko)
45// 05-02-03 Fix compilation warnings (V.Ivanchenko)
46// 13-02-03 SubCutoffProcessors defined for regions (V.Ivanchenko)
47// 17-02-03 Fix problem of store/restore tables (V.Ivanchenko)
48// 26-02-03 Region dependent step limit (V.Ivanchenko)
49// 26-03-03 Add GetDEDXDispersion (V.Ivanchenko)
50// 09-04-03 Fix problem of negative range limit for non integral (V.Ivanchenko)
51// 13-05-03 Add calculation of precise range (V.Ivanchenko)
52// 21-07-03 Add UpdateEmModel method (V.Ivanchenko)
53// 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
54// 14-01-04 Activate precise range calculation (V.Ivanchenko)
55// 10-03-04 Fix problem of step limit calculation (V.Ivanchenko)
56// 30-06-04 make destructor virtual (V.Ivanchenko)
57// 05-07-04 fix problem of GenericIons seen at small cuts (V.Ivanchenko)
58// 03-08-04 Add DEDX table to all processes for control on integral range(VI)
59// 06-08-04 Clear up names of member functions (V.Ivanchenko)
60// 27-08-04 Add NeedBuildTables method (V.Ivanchneko)
61// 09-09-04 Bug fix for the integral mode with 2 peaks (V.Ivanchneko)
62// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko)
63// 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
64// 11-04-05 Use MaxSecondaryEnergy from a model (V.Ivanchenko)
65// 10-01-05 Remove SetStepLimits (V.Ivanchenko)
66// 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko)
67// 13-01-06 Remove AddSubCutSecondaries and cleanup (V.Ivantchenko)
68// 20-01-06 Introduce G4EmTableType and reducing number of methods (VI)
69// 26-01-06 Add public method GetCSDARange (V.Ivanchenko)
70// 22-03-06 Add SetDynamicMassCharge (V.Ivanchenko)
71// 23-03-06 Use isIonisation flag (V.Ivanchenko)
72// 13-05-06 Add method to access model by index (V.Ivanchenko)
73// 14-01-07 add SetEmModel(index) and SetFluctModel() (mma)
74// 15-01-07 Add separate ionisation tables and reorganise get/set methods for
75// dedx tables (V.Ivanchenko)
76// 13-03-07 use SafetyHelper instead of navigator (V.Ivanchenko)
77// 27-07-07 use stl vector for emModels instead of C-array (V.Ivanchenko)
78// 25-09-07 More accurate handling zero xsect in
79// PostStepGetPhysicalInteractionLength (V.Ivanchenko)
80// 27-10-07 Virtual functions moved to source (V.Ivanchenko)
[961]81// 15-07-08 Reorder class members for further multi-thread development (VI)
[819]82//
83// Class Description:
84//
85// It is the unified energy loss process it calculates the continuous
86// energy loss for charged particles using a set of Energy Loss
87// models valid for different energy regions. There are a possibility
88// to create and access to dE/dx and range tables, or to calculate
89// that information on fly.
90
91// -------------------------------------------------------------------
92//
93
94#ifndef G4VEnergyLossProcess_h
95#define G4VEnergyLossProcess_h 1
96
97#include "G4VContinuousDiscreteProcess.hh"
98#include "globals.hh"
99#include "G4Material.hh"
100#include "G4MaterialCutsCouple.hh"
101#include "G4Track.hh"
102#include "G4EmModelManager.hh"
103#include "G4UnitsTable.hh"
104#include "G4ParticleChangeForLoss.hh"
105#include "G4EmTableType.hh"
106#include "G4PhysicsTable.hh"
107#include "G4PhysicsVector.hh"
108
109class G4Step;
110class G4ParticleDefinition;
111class G4VEmModel;
112class G4VEmFluctuationModel;
113class G4DataVector;
114class G4Region;
115class G4SafetyHelper;
[1340]116class G4VAtomDeexcitation;
[819]117
118//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
119
120class G4VEnergyLossProcess : public G4VContinuousDiscreteProcess
121{
122public:
123
124 G4VEnergyLossProcess(const G4String& name = "EnergyLoss",
125 G4ProcessType type = fElectromagnetic);
126
127 virtual ~G4VEnergyLossProcess();
128
[1055]129private:
130 // clean vectors and arrays
131 void Clean();
132
[819]133 //------------------------------------------------------------------------
134 // Virtual methods to be implemented in concrete processes
135 //------------------------------------------------------------------------
136
[1055]137public:
[819]138 virtual G4bool IsApplicable(const G4ParticleDefinition& p) = 0;
139
140 virtual void PrintInfo() = 0;
141
142protected:
143
144 virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition*,
145 const G4ParticleDefinition*) = 0;
146
147 //------------------------------------------------------------------------
148 // Methods with standard implementation; may be overwritten if needed
149 //------------------------------------------------------------------------
150
151 virtual G4double MinPrimaryEnergy(const G4ParticleDefinition*,
152 const G4Material*, G4double cut);
153
154 //------------------------------------------------------------------------
[1055]155 // Virtual methods implementation common to all EM ContinuousDiscrete
156 // processes. Further inheritance is not assumed
[819]157 //------------------------------------------------------------------------
[961]158
[819]159public:
160
[1055]161 // prepare all tables
[819]162 void PreparePhysicsTable(const G4ParticleDefinition&);
163
[1055]164 // build all tables
[819]165 void BuildPhysicsTable(const G4ParticleDefinition&);
166
[1055]167 // build a table
168 G4PhysicsTable* BuildDEDXTable(G4EmTableType tType = fRestricted);
169
170 // build a table
171 G4PhysicsTable* BuildLambdaTable(G4EmTableType tType = fRestricted);
172
173 // summary printout after initialisation
174 void PrintInfoDefinition();
175
176 // Add subcutoff option for the region
177 void ActivateSubCutoff(G4bool val, const G4Region* region = 0);
178
179 // Activate deexcitation code for region
180 void ActivateDeexcitation(G4bool, const G4Region* region = 0);
181
182 // Step limit from AlongStep
[961]183 G4double AlongStepGetPhysicalInteractionLength(const G4Track&,
184 G4double previousStepSize,
185 G4double currentMinimumStep,
186 G4double& currentSafety,
187 G4GPILSelection* selection);
188
[1055]189 // Step limit from cross section
[961]190 G4double PostStepGetPhysicalInteractionLength(const G4Track& track,
191 G4double previousStepSize,
192 G4ForceCondition* condition);
193
[1055]194 // AlongStep computations
[819]195 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&);
196
[1055]197 // Sampling of secondaries in vicinity of geometrical boundary
198 void SampleSubCutSecondaries(std::vector<G4Track*>&, const G4Step&,
[1315]199 G4VEmModel* model, G4int matIdx);
[1055]200
201 // PostStep sampling of secondaries
[819]202 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
203
[1055]204 // Store all PhysicsTable in files.
205 // Return false in case of any fatal failure at I/O
[819]206 G4bool StorePhysicsTable(const G4ParticleDefinition*,
207 const G4String& directory,
208 G4bool ascii = false);
209
[1055]210 // Retrieve all Physics from a files.
211 // Return true if all the Physics Table are built.
212 // Return false if any fatal failure.
[819]213 G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
214 const G4String& directory,
215 G4bool ascii);
216
[1055]217private:
218 // store a table
219 G4bool StoreTable(const G4ParticleDefinition* p,
220 G4PhysicsTable*, G4bool ascii,
221 const G4String& directory,
222 const G4String& tname);
[819]223
[1055]224 // retrieve a table
225 G4bool RetrieveTable(const G4ParticleDefinition* p,
226 G4PhysicsTable*, G4bool ascii,
227 const G4String& directory,
228 const G4String& tname,
229 G4bool mandatory);
[819]230
231 //------------------------------------------------------------------------
[1055]232 // Public interface to cross section, mfp and sampling of fluctuations
233 // These methods are not used in run time
[819]234 //------------------------------------------------------------------------
235
236public:
[1055]237 // access to dispersion of restricted energy loss
[819]238 G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple,
239 const G4DynamicParticle* dp,
240 G4double length);
241
[1055]242 // Access to cross section table
243 G4double CrossSectionPerVolume(G4double kineticEnergy,
244 const G4MaterialCutsCouple* couple);
[819]245
[1055]246 // access to cross section
247 G4double MeanFreePath(const G4Track& track);
[819]248
[1055]249 // access to step limit
250 G4double ContinuousStepLimit(const G4Track& track,
251 G4double previousStepSize,
252 G4double currentMinimumStep,
253 G4double& currentSafety);
[819]254
[1055]255protected:
[819]256
[1055]257 // implementation of the pure virtual method
258 G4double GetMeanFreePath(const G4Track& track,
259 G4double previousStepSize,
260 G4ForceCondition* condition);
[819]261
[1055]262 // implementation of the pure virtual method
263 G4double GetContinuousStepLimit(const G4Track& track,
264 G4double previousStepSize,
265 G4double currentMinimumStep,
266 G4double& currentSafety);
[819]267
[1055]268 //------------------------------------------------------------------------
269 // Run time method which may be also used by derived processes
270 //------------------------------------------------------------------------
[819]271
[1055]272 // creeation of an empty vector for cross section
273 G4PhysicsVector* LambdaPhysicsVector(const G4MaterialCutsCouple*,
274 G4double cut);
[819]275
[1055]276 inline size_t CurrentMaterialCutsCoupleIndex() const;
[819]277
[1055]278 inline G4double GetCurrentRange() const;
[819]279
280 //------------------------------------------------------------------------
[1055]281 // Specific methods to set, access, modify models
[819]282 //------------------------------------------------------------------------
283
[1055]284 // Select model in run time
285 inline void SelectModel(G4double kinEnergy);
[819]286
[1055]287public:
288 // Select model by energy and region index
289 inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
290 size_t& idx) const;
[819]291
[1055]292 // Add EM model coupled with fluctuation model for region, smaller value
293 // of order defines which pair of models will be selected for a given
294 // energy interval
295 void AddEmModel(G4int, G4VEmModel*,
296 G4VEmFluctuationModel* fluc = 0,
297 const G4Region* region = 0);
[819]298
[1055]299 // Define new energy range for the model identified by the name
300 void UpdateEmModel(const G4String&, G4double, G4double);
301
[819]302 // Assign a model to a process
[1055]303 void SetEmModel(G4VEmModel*, G4int index=1);
[819]304
305 // return the assigned model
[1055]306 G4VEmModel* EmModel(G4int index=1);
[819]307
[1055]308 // Access to models
309 G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false);
310
311 G4int NumberOfModels();
312
[819]313 // Assign a fluctuation model to a process
[1055]314 void SetFluctModel(G4VEmFluctuationModel*);
[819]315
316 // return the assigned fluctuation model
317 inline G4VEmFluctuationModel* FluctModel();
318
[1055]319 //------------------------------------------------------------------------
320 // Define and access particle type
321 //------------------------------------------------------------------------
[819]322
[1055]323protected:
324 inline void SetParticle(const G4ParticleDefinition* p);
325 inline void SetSecondaryParticle(const G4ParticleDefinition* p);
[819]326
[1055]327public:
328 inline void SetBaseParticle(const G4ParticleDefinition* p);
329 inline const G4ParticleDefinition* Particle() const;
330 inline const G4ParticleDefinition* BaseParticle() const;
331 inline const G4ParticleDefinition* SecondaryParticle() const;
[819]332
333 //------------------------------------------------------------------------
[1055]334 // Get/set parameters to configure the process at initialisation time
[819]335 //------------------------------------------------------------------------
336
[1055]337 // Add subcutoff process (bremsstrahlung) to sample secondary
338 // particle production in vicinity of the geometry boundary
339 void AddCollaborativeProcess(G4VEnergyLossProcess*);
340
[819]341 inline void SetLossFluctuations(G4bool val);
342 inline void SetRandomStep(G4bool val);
[1055]343
[819]344 inline void SetIntegral(G4bool val);
345 inline G4bool IsIntegral() const;
346
347 // Set/Get flag "isIonisation"
348 inline void SetIonisation(G4bool val);
349 inline G4bool IsIonisationProcess() const;
350
351 // Redefine parameteters for stepping control
352 //
353 inline void SetLinearLossLimit(G4double val);
354 inline void SetMinSubRange(G4double val);
[1055]355 inline void SetLambdaFactor(G4double val);
[1007]356 inline void SetStepFunction(G4double v1, G4double v2);
[1196]357 inline void SetLowestEnergyLimit(G4double);
[819]358
359 inline G4int NumberOfSubCutoffRegions() const;
[1055]360 inline G4int NumberOfDERegions() const;
[819]361
362 //------------------------------------------------------------------------
[1055]363 // Specific methods to path Physics Tables to the process
[819]364 //------------------------------------------------------------------------
365
[1055]366 void SetDEDXTable(G4PhysicsTable* p, G4EmTableType tType);
367 void SetCSDARangeTable(G4PhysicsTable* pRange);
368 void SetRangeTableForLoss(G4PhysicsTable* p);
369 void SetSecondaryRangeTable(G4PhysicsTable* p);
370 void SetInverseRangeTable(G4PhysicsTable* p);
[819]371
[1055]372 void SetLambdaTable(G4PhysicsTable* p);
373 void SetSubLambdaTable(G4PhysicsTable* p);
[819]374
[1055]375 // Binning for dEdx, range, inverse range and labda tables
376 inline void SetDEDXBinning(G4int nbins);
377 inline void SetLambdaBinning(G4int nbins);
[819]378
[1055]379 // Binning for dEdx, range, and inverse range tables
380 inline void SetDEDXBinningForCSDARange(G4int nbins);
[819]381
[1055]382 // Min kinetic energy for tables
383 inline void SetMinKinEnergy(G4double e);
384 inline G4double MinKinEnergy() const;
[819]385
[1055]386 // Max kinetic energy for tables
387 inline void SetMaxKinEnergy(G4double e);
388 inline G4double MaxKinEnergy() const;
[819]389
[1055]390 // Max kinetic energy for tables
391 inline void SetMaxKinEnergyForCSDARange(G4double e);
[819]392
[1055]393 // Return values for given G4MaterialCutsCouple
394 inline G4double GetDEDX(G4double& kineticEnergy, const G4MaterialCutsCouple*);
395 inline G4double GetDEDXForSubsec(G4double& kineticEnergy,
396 const G4MaterialCutsCouple*);
397 inline G4double GetRange(G4double& kineticEnergy, const G4MaterialCutsCouple*);
398 inline G4double GetCSDARange(G4double& kineticEnergy, const G4MaterialCutsCouple*);
399 inline G4double GetRangeForLoss(G4double& kineticEnergy, const G4MaterialCutsCouple*);
400 inline G4double GetKineticEnergy(G4double& range, const G4MaterialCutsCouple*);
401 inline G4double GetLambda(G4double& kineticEnergy, const G4MaterialCutsCouple*);
[819]402
[1055]403 inline G4bool TablesAreBuilt() const;
[819]404
[1055]405 // Access to specific tables
406 inline G4PhysicsTable* DEDXTable() const;
407 inline G4PhysicsTable* DEDXTableForSubsec() const;
408 inline G4PhysicsTable* DEDXunRestrictedTable() const;
409 inline G4PhysicsTable* IonisationTable() const;
410 inline G4PhysicsTable* IonisationTableForSubsec() const;
411 inline G4PhysicsTable* CSDARangeTable() const;
412 inline G4PhysicsTable* RangeTableForLoss() const;
413 inline G4PhysicsTable* InverseRangeTable() const;
414 inline G4PhysicsTable* LambdaTable();
415 inline G4PhysicsTable* SubLambdaTable();
[819]416
[961]417 //------------------------------------------------------------------------
[1055]418 // Run time method for simulation of ionisation
[961]419 //------------------------------------------------------------------------
[819]420
[1315]421 // access atom on which interaction happens
422 const G4Element* GetCurrentElement() const;
423
[1055]424 // sample range at the end of a step
425 inline G4double SampleRange();
[819]426
[1055]427 // Set scaling parameters for ions is needed to G4EmCalculator
428 inline void SetDynamicMassCharge(G4double massratio, G4double charge2ratio);
[819]429
[1055]430private:
[819]431
[961]432 // define material and indexes
433 inline void DefineMaterial(const G4MaterialCutsCouple* couple);
[819]434
[1055]435 //------------------------------------------------------------------------
436 // Compute values using scaling relation, mass and charge of based particle
437 //------------------------------------------------------------------------
438
[819]439 inline G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy);
440 inline G4double GetSubDEDXForScaledEnergy(G4double scaledKinEnergy);
441 inline G4double GetIonisationForScaledEnergy(G4double scaledKinEnergy);
442 inline G4double GetSubIonisationForScaledEnergy(G4double scaledKinEnergy);
443 inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy);
444 inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinEnergy);
[1055]445 inline G4double ScaledKinEnergyForLoss(G4double range);
[1007]446 inline G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy);
[819]447 inline void ComputeLambdaForScaledEnergy(G4double scaledKinEnergy);
448
449 // hide assignment operator
450 G4VEnergyLossProcess(G4VEnergyLossProcess &);
451 G4VEnergyLossProcess & operator=(const G4VEnergyLossProcess &right);
452
[961]453 // ======== Parameters of the class fixed at construction =========
[819]454
[961]455 G4EmModelManager* modelManager;
456 G4SafetyHelper* safetyHelper;
[819]457
[961]458 const G4ParticleDefinition* secondaryParticle;
459 const G4ParticleDefinition* theElectron;
460 const G4ParticleDefinition* thePositron;
461 const G4ParticleDefinition* theGenericIon;
[819]462
[961]463 G4PhysicsVector* vstrag;
[819]464
[961]465 // ======== Parameters of the class fixed at initialisation =======
466
[819]467 std::vector<G4VEmModel*> emModels;
468 G4VEmFluctuationModel* fluctModel;
[1340]469 G4VAtomDeexcitation* atomDeexcitation;
[819]470 std::vector<const G4Region*> scoffRegions;
[1055]471 std::vector<const G4Region*> deRegions;
[819]472 G4int nSCoffRegions;
[1055]473 G4int nDERegions;
474 G4bool* idxSCoffRegions;
475 G4bool* idxDERegions;
[961]476
[819]477 std::vector<G4VEnergyLossProcess*> scProcesses;
478 G4int nProcesses;
479
480 // tables and vectors
481 G4PhysicsTable* theDEDXTable;
482 G4PhysicsTable* theDEDXSubTable;
483 G4PhysicsTable* theDEDXunRestrictedTable;
484 G4PhysicsTable* theIonisationTable;
485 G4PhysicsTable* theIonisationSubTable;
486 G4PhysicsTable* theRangeTableForLoss;
487 G4PhysicsTable* theCSDARangeTable;
488 G4PhysicsTable* theSecondaryRangeTable;
489 G4PhysicsTable* theInverseRangeTable;
490 G4PhysicsTable* theLambdaTable;
491 G4PhysicsTable* theSubLambdaTable;
492 G4double* theDEDXAtMaxEnergy;
493 G4double* theRangeAtMaxEnergy;
494 G4double* theEnergyOfCrossSectionMax;
495 G4double* theCrossSectionMax;
496
497 const G4DataVector* theCuts;
498 const G4DataVector* theSubCuts;
499
500 const G4ParticleDefinition* baseParticle;
501
502 G4int nBins;
503 G4int nBinsCSDA;
504
505 G4double lowestKinEnergy;
506 G4double minKinEnergy;
507 G4double maxKinEnergy;
508 G4double maxKinEnergyCSDA;
509
510 G4double linLossLimit;
511 G4double minSubRange;
512 G4double dRoverRange;
513 G4double finalRange;
514 G4double lambdaFactor;
515
516 G4bool lossFluctuationFlag;
517 G4bool rndmStepFlag;
518 G4bool tablesAreBuilt;
519 G4bool integral;
[961]520 G4bool isIon;
[819]521 G4bool isIonisation;
522 G4bool useSubCutoff;
[1055]523 G4bool useDeexcitation;
[819]524
[961]525protected:
[819]526
[961]527 G4ParticleChangeForLoss fParticleChange;
[819]528
[961]529 // ======== Cashed values - may be state dependent ================
[819]530
[961]531private:
[819]532
[961]533 std::vector<G4DynamicParticle*> secParticles;
534 std::vector<G4Track*> scTracks;
[819]535
[961]536 const G4ParticleDefinition* particle;
[819]537
[961]538 G4VEmModel* currentModel;
539 const G4Material* currentMaterial;
540 const G4MaterialCutsCouple* currentCouple;
541 size_t currentMaterialIndex;
542
543 G4int nWarnings;
544
545 G4double massRatio;
546 G4double reduceFactor;
547 G4double chargeSqRatio;
548
549 G4double preStepLambda;
550 G4double fRange;
551 G4double preStepKinEnergy;
552 G4double preStepScaledEnergy;
553 G4double mfpKinEnergy;
554
555 G4GPILSelection aGPILSelection;
556
557};
558
[1315]559// ======== Run time inline methods ================
[819]560
[1055]561inline size_t G4VEnergyLossProcess::CurrentMaterialCutsCoupleIndex() const
[819]562{
[1055]563 return currentMaterialIndex;
[819]564}
565
566//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[1055]567
568inline G4double G4VEnergyLossProcess::GetCurrentRange() const
[819]569{
[1055]570 return fRange;
[819]571}
572
573//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[1007]574
[1055]575inline void G4VEnergyLossProcess::SelectModel(G4double kinEnergy)
[819]576{
[1055]577 currentModel = modelManager->SelectModel(kinEnergy, currentMaterialIndex);
578 currentModel->SetCurrentCouple(currentCouple);
[819]579}
580
581//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
582
[1055]583inline G4VEmModel* G4VEnergyLossProcess::SelectModelForMaterial(
584 G4double kinEnergy, size_t& idx) const
[819]585{
[1055]586 return modelManager->SelectModel(kinEnergy, idx);
[819]587}
588
589//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
590
[1315]591inline void
592G4VEnergyLossProcess::DefineMaterial(const G4MaterialCutsCouple* couple)
[819]593{
[1315]594 if(couple != currentCouple) {
595 currentCouple = couple;
596 currentMaterial = couple->GetMaterial();
597 currentMaterialIndex = couple->GetIndex();
598 mfpKinEnergy = DBL_MAX;
599 }
[819]600}
601
602//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
603
[1315]604inline void G4VEnergyLossProcess::SetDynamicMassCharge(G4double massratio,
605 G4double charge2ratio)
[819]606{
[1315]607 massRatio = massratio;
608 chargeSqRatio = charge2ratio;
609 reduceFactor = 1.0/(chargeSqRatio*massRatio);
[819]610}
611
612//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
613
[1315]614inline G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e)
[819]615{
[1315]616 G4double x = ((*theDEDXTable)[currentMaterialIndex]->Value(e))*chargeSqRatio;
617 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
618 return x;
[819]619}
620
621//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
622
[1315]623inline G4double G4VEnergyLossProcess::GetSubDEDXForScaledEnergy(G4double e)
[819]624{
[1315]625 G4double x = ((*theDEDXSubTable)[currentMaterialIndex]->Value(e))*chargeSqRatio;
626 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
627 return x;
[819]628}
629
630//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
631
[1315]632inline G4double G4VEnergyLossProcess::GetIonisationForScaledEnergy(G4double e)
[819]633{
[1315]634 //G4double x = 0.0;
635 // if(theIonisationTable) {
636 G4double x = ((*theIonisationTable)[currentMaterialIndex]->Value(e))*chargeSqRatio;
637 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
638 //}
639 return x;
[819]640}
641
642//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
643
[1315]644inline
645G4double G4VEnergyLossProcess::GetSubIonisationForScaledEnergy(G4double e)
[819]646{
[1315]647 // G4double x = 0.0;
648 //if(theIonisationSubTable) {
649 G4double x = ((*theIonisationSubTable)[currentMaterialIndex]->Value(e))*chargeSqRatio;
650 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
651 //}
652 return x;
[819]653}
654
655//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
656
[1315]657inline G4double G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e)
[819]658{
[1315]659 G4double x = ((*theRangeTableForLoss)[currentMaterialIndex])->Value(e);
660 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
661 return x;
[819]662}
663
664//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
665
[1315]666inline G4double
667G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e)
[819]668{
[1315]669 G4double x;
670
671 if (e < maxKinEnergyCSDA) {
672 x = ((*theCSDARangeTable)[currentMaterialIndex])->Value(e);
673 if(e < minKinEnergy) x *= std::sqrt(e/minKinEnergy);
674 } else {
675 x = theRangeAtMaxEnergy[currentMaterialIndex] +
676 (e - maxKinEnergyCSDA)/theDEDXAtMaxEnergy[currentMaterialIndex];
677 }
678 return x;
[819]679}
680
681//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
682
[1315]683inline G4double G4VEnergyLossProcess::ScaledKinEnergyForLoss(G4double r)
[819]684{
[1315]685 G4PhysicsVector* v = (*theInverseRangeTable)[currentMaterialIndex];
686 G4double rmin = v->Energy(0);
687 G4double e = 0.0;
688 if(r >= rmin) { e = v->Value(r); }
689 else if(r > 0.0) {
690 G4double x = r/rmin;
691 e = minKinEnergy*x*x;
692 }
693 return e;
[819]694}
695
696//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
697
[1315]698inline G4double G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e)
[819]699{
[1315]700 return chargeSqRatio*(((*theLambdaTable)[currentMaterialIndex])->Value(e));
[819]701}
702
703//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
704
[1315]705inline G4double
706G4VEnergyLossProcess::GetDEDX(G4double& kineticEnergy,
707 const G4MaterialCutsCouple* couple)
[819]708{
[1315]709 DefineMaterial(couple);
710 return GetDEDXForScaledEnergy(kineticEnergy*massRatio);
[819]711}
712
713//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[1315]714
715inline G4double
716G4VEnergyLossProcess::GetDEDXForSubsec(G4double& kineticEnergy,
717 const G4MaterialCutsCouple* couple)
[819]718{
[1315]719 DefineMaterial(couple);
720 return GetSubDEDXForScaledEnergy(kineticEnergy*massRatio);
[819]721}
722
723//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
724
[1315]725inline G4double
726G4VEnergyLossProcess::GetRange(G4double& kineticEnergy,
727 const G4MaterialCutsCouple* couple)
[819]728{
[1315]729 G4double x = fRange;
730 if(kineticEnergy != preStepKinEnergy || couple != currentCouple) {
731 DefineMaterial(couple);
732 if(theCSDARangeTable)
733 x = GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)
734 * reduceFactor;
735 else if(theRangeTableForLoss)
736 x = GetScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
737 }
738 return x;
[819]739}
740
741//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
742
[1315]743inline G4double
744G4VEnergyLossProcess::GetCSDARange(G4double& kineticEnergy,
745 const G4MaterialCutsCouple* couple)
[819]746{
[1315]747 DefineMaterial(couple);
748 G4double x = DBL_MAX;
749 if(theCSDARangeTable)
750 x = GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)
751 * reduceFactor;
752 return x;
[819]753}
754
755//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
756
[1315]757inline G4double
758G4VEnergyLossProcess::GetRangeForLoss(G4double& kineticEnergy,
759 const G4MaterialCutsCouple* couple)
[819]760{
[1315]761 DefineMaterial(couple);
762 G4double x = DBL_MAX;
763 if(theRangeTableForLoss)
764 x = GetScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
765 // G4cout << "Range from " << GetProcessName()
766 // << " e= " << kineticEnergy << " r= " << x << G4endl;
767 return x;
[819]768}
769
770//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
771
[1315]772inline G4double
773G4VEnergyLossProcess::GetKineticEnergy(G4double& range,
774 const G4MaterialCutsCouple* couple)
[819]775{
[1315]776 DefineMaterial(couple);
777 G4double r = range/reduceFactor;
778 G4double e = ScaledKinEnergyForLoss(r)/massRatio;
779 return e;
[819]780}
781
782//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
783
[1315]784inline G4double
785G4VEnergyLossProcess::GetLambda(G4double& kineticEnergy,
786 const G4MaterialCutsCouple* couple)
[819]787{
[1315]788 DefineMaterial(couple);
789 G4double x = 0.0;
790 if(theLambdaTable) { x = GetLambdaForScaledEnergy(kineticEnergy*massRatio); }
791 return x;
[819]792}
793
794//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
795
[1315]796inline void G4VEnergyLossProcess::ComputeLambdaForScaledEnergy(G4double e)
[819]797{
[1315]798 mfpKinEnergy = theEnergyOfCrossSectionMax[currentMaterialIndex];
799 if (e <= mfpKinEnergy) {
800 preStepLambda = GetLambdaForScaledEnergy(e);
801
802 } else {
803 G4double e1 = e*lambdaFactor;
804 if(e1 > mfpKinEnergy) {
805 preStepLambda = GetLambdaForScaledEnergy(e);
806 G4double preStepLambda1 = GetLambdaForScaledEnergy(e1);
807 if(preStepLambda1 > preStepLambda) {
808 mfpKinEnergy = e1;
809 preStepLambda = preStepLambda1;
810 }
811 } else {
812 preStepLambda = chargeSqRatio*theCrossSectionMax[currentMaterialIndex];
813 }
814 }
[819]815}
816
817//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[1007]818
[1315]819inline G4double G4VEnergyLossProcess::SampleRange()
[1196]820{
[1315]821 G4double e = amu_c2*preStepKinEnergy/particle->GetPDGMass();
822 G4double s = fRange*std::pow(10.,vstrag->Value(e));
823 G4double x = fRange + G4RandGauss::shoot(0.0,s);
824 if(x > 0.0) { fRange = x; }
825 return fRange;
[1196]826}
827
[1315]828// ======== Get/Set inline methods used at initialisation ================
[1196]829
[1315]830inline void G4VEnergyLossProcess::SetFluctModel(G4VEmFluctuationModel* p)
[819]831{
[1315]832 fluctModel = p;
[819]833}
834
835//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
836
[1315]837inline G4VEmFluctuationModel* G4VEnergyLossProcess::FluctModel()
[819]838{
[1315]839 return fluctModel;
[819]840}
841
842//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
843
[1315]844inline void G4VEnergyLossProcess::SetParticle(const G4ParticleDefinition* p)
[819]845{
[1315]846 particle = p;
[819]847}
848
849//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
850
[1315]851inline void G4VEnergyLossProcess::SetSecondaryParticle(const G4ParticleDefinition* p)
[961]852{
[1315]853 secondaryParticle = p;
[961]854}
[819]855
856//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
857
[1315]858inline void G4VEnergyLossProcess::SetBaseParticle(const G4ParticleDefinition* p)
[819]859{
[1315]860 baseParticle = p;
[819]861}
862
863//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
864
[1315]865inline const G4ParticleDefinition* G4VEnergyLossProcess::Particle() const
[819]866{
[1315]867 return particle;
[819]868}
869
870//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
871
[1315]872inline const G4ParticleDefinition* G4VEnergyLossProcess::BaseParticle() const
[819]873{
[1315]874 return baseParticle;
[819]875}
876
877//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
878
[1315]879inline const G4ParticleDefinition* G4VEnergyLossProcess::SecondaryParticle() const
[819]880{
[1315]881 return secondaryParticle;
[819]882}
883
884//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
885
[1315]886inline void G4VEnergyLossProcess::SetLossFluctuations(G4bool val)
[819]887{
[1315]888 lossFluctuationFlag = val;
[819]889}
890
891//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
892
[1315]893inline void G4VEnergyLossProcess::SetRandomStep(G4bool val)
[819]894{
[1315]895 rndmStepFlag = val;
[819]896}
897
898//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
899
[1315]900inline void G4VEnergyLossProcess::SetIntegral(G4bool val)
[819]901{
[1315]902 integral = val;
[819]903}
904
905//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[1315]906
907inline G4bool G4VEnergyLossProcess::IsIntegral() const
[819]908{
[1315]909 return integral;
[819]910}
911
912//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
913
[1315]914inline void G4VEnergyLossProcess::SetIonisation(G4bool val)
[819]915{
[1315]916 isIonisation = val;
917 if(val) { aGPILSelection = CandidateForSelection; }
918 else { aGPILSelection = NotCandidateForSelection; }
[819]919}
920
921//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
922
[1315]923inline G4bool G4VEnergyLossProcess::IsIonisationProcess() const
[819]924{
[1315]925 return isIonisation;
[819]926}
927
928//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[961]929
[1315]930inline void G4VEnergyLossProcess::SetLinearLossLimit(G4double val)
[819]931{
[1315]932 linLossLimit = val;
[819]933}
934
935//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[1005]936
[1315]937inline void G4VEnergyLossProcess::SetMinSubRange(G4double val)
[991]938{
[1315]939 minSubRange = val;
[991]940}
[819]941
[991]942//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
943
[1315]944inline void G4VEnergyLossProcess::SetLambdaFactor(G4double val)
[819]945{
[1315]946 if(val > 0.0 && val <= 1.0) { lambdaFactor = val; }
[819]947}
948
949//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
950
[1315]951void G4VEnergyLossProcess::SetStepFunction(G4double v1, G4double v2)
[819]952{
[1315]953 dRoverRange = v1;
954 finalRange = v2;
955 if (dRoverRange > 0.999) { dRoverRange = 1.0; }
956 currentCouple = 0;
957 mfpKinEnergy = DBL_MAX;
[819]958}
959
960//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[1005]961
[1315]962inline void G4VEnergyLossProcess::SetLowestEnergyLimit(G4double val)
[991]963{
[1315]964 lowestKinEnergy = val;
[991]965}
[819]966
[991]967//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
968
[1315]969inline G4int G4VEnergyLossProcess::NumberOfSubCutoffRegions() const
[819]970{
[1315]971 return nSCoffRegions;
[819]972}
973
974//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
975
[1315]976inline G4int G4VEnergyLossProcess::NumberOfDERegions() const
[819]977{
[1315]978 return nDERegions;
[819]979}
980
981//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
982
[1315]983inline void G4VEnergyLossProcess::SetDEDXBinning(G4int nbins)
[819]984{
[1315]985 nBins = nbins;
[819]986}
987
988//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
989
[1315]990inline void G4VEnergyLossProcess::SetLambdaBinning(G4int nbins)
[819]991{
[1315]992 nBins = nbins;
[819]993}
994
995//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
996
[1315]997inline void G4VEnergyLossProcess::SetDEDXBinningForCSDARange(G4int nbins)
[819]998{
[1315]999 nBinsCSDA = nbins;
[819]1000}
1001
1002//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1003
[1315]1004inline void G4VEnergyLossProcess::SetMinKinEnergy(G4double e)
[819]1005{
[1315]1006 minKinEnergy = e;
[819]1007}
1008
1009//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1010
[1315]1011inline G4double G4VEnergyLossProcess::MinKinEnergy() const
[819]1012{
[1315]1013 return minKinEnergy;
[819]1014}
1015
1016//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1017
[1315]1018inline void G4VEnergyLossProcess::SetMaxKinEnergy(G4double e)
[819]1019{
[1315]1020 maxKinEnergy = e;
1021 if(e < maxKinEnergyCSDA) { maxKinEnergyCSDA = e; }
[819]1022}
1023
1024//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1025
[1315]1026inline G4double G4VEnergyLossProcess::MaxKinEnergy() const
[819]1027{
[1315]1028 return maxKinEnergy;
[819]1029}
1030
1031//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1032
[1315]1033inline void G4VEnergyLossProcess::SetMaxKinEnergyForCSDARange(G4double e)
[819]1034{
[1315]1035 maxKinEnergyCSDA = e;
[819]1036}
1037
[1315]1038
[819]1039//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1040
[1315]1041inline G4bool G4VEnergyLossProcess::TablesAreBuilt() const
[819]1042{
[1315]1043 return tablesAreBuilt;
[819]1044}
1045
1046//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1047
[1315]1048inline G4PhysicsTable* G4VEnergyLossProcess::DEDXTable() const
[819]1049{
[1315]1050 return theDEDXTable;
[819]1051}
1052
1053//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1054
[1315]1055inline G4PhysicsTable* G4VEnergyLossProcess::DEDXTableForSubsec() const
[819]1056{
[1315]1057 return theDEDXSubTable;
[819]1058}
1059
1060//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1061
[1315]1062inline G4PhysicsTable* G4VEnergyLossProcess::DEDXunRestrictedTable() const
[819]1063{
[1315]1064 return theDEDXunRestrictedTable;
[819]1065}
1066
1067//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1068
[1315]1069inline G4PhysicsTable* G4VEnergyLossProcess::IonisationTable() const
[819]1070{
[1315]1071 G4PhysicsTable* t = theDEDXTable;
1072 if(theIonisationTable) { t = theIonisationTable; }
1073 return t;
[819]1074}
1075
1076//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1077
[1315]1078inline G4PhysicsTable* G4VEnergyLossProcess::IonisationTableForSubsec() const
[819]1079{
[1315]1080 G4PhysicsTable* t = theDEDXSubTable;
1081 if(theIonisationSubTable) { t = theIonisationSubTable; }
1082 return t;
[819]1083}
1084
1085//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1086
[1315]1087inline G4PhysicsTable* G4VEnergyLossProcess::CSDARangeTable() const
[819]1088{
[1315]1089 return theCSDARangeTable;
[819]1090}
1091
1092//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1093
[1315]1094inline G4PhysicsTable* G4VEnergyLossProcess::RangeTableForLoss() const
[819]1095{
[1315]1096 return theRangeTableForLoss;
[819]1097}
1098
1099//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1100
[1315]1101inline G4PhysicsTable* G4VEnergyLossProcess::InverseRangeTable() const
[819]1102{
[1315]1103 return theInverseRangeTable;
[819]1104}
1105
1106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1107
[1315]1108inline G4PhysicsTable* G4VEnergyLossProcess::LambdaTable()
[819]1109{
[1315]1110 return theLambdaTable;
[819]1111}
1112
1113//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1114
[1315]1115inline G4PhysicsTable* G4VEnergyLossProcess::SubLambdaTable()
[819]1116{
[1315]1117 return theSubLambdaTable;
[819]1118}
1119
1120//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1121
1122#endif
Note: See TracBrowser for help on using the repository browser.