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

Last change on this file since 1047 was 1007, checked in by garnier, 15 years ago

update to geant4.9.2

File size: 37.4 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//
[1007]26// $Id: G4VEnergyLossProcess.hh,v 1.83 2008/09/12 16:19:01 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;
116
117//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
118
119class G4VEnergyLossProcess : public G4VContinuousDiscreteProcess
120{
121public:
122
123  G4VEnergyLossProcess(const G4String& name = "EnergyLoss",
124                       G4ProcessType type = fElectromagnetic);
125
126  virtual ~G4VEnergyLossProcess();
127
128  //------------------------------------------------------------------------
129  // Virtual methods to be implemented in concrete processes
130  //------------------------------------------------------------------------
131
132  virtual G4bool IsApplicable(const G4ParticleDefinition& p) = 0;
133 
134  virtual void PrintInfo() = 0;
135
136protected:
137
138  virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition*,
139                                           const G4ParticleDefinition*) = 0;
140
141  //------------------------------------------------------------------------
142  // Methods with standard implementation; may be overwritten if needed
143  //------------------------------------------------------------------------
144
[1007]145protected:
146
[819]147  virtual G4double MinPrimaryEnergy(const G4ParticleDefinition*,
148                                    const G4Material*, G4double cut);
149
150  //------------------------------------------------------------------------
[1007]151  // Virtual methods common to all EM ContinuousDiscrete processes
152  // Further inheritance is not assumed
[819]153  //------------------------------------------------------------------------
[961]154
[819]155public:
156
[1007]157  void PrintInfoDefinition();
158
[819]159  void PreparePhysicsTable(const G4ParticleDefinition&);
160
161  void BuildPhysicsTable(const G4ParticleDefinition&);
162
[961]163  G4double AlongStepGetPhysicalInteractionLength(const G4Track&,
164                                                 G4double  previousStepSize,
165                                                 G4double  currentMinimumStep,
166                                                 G4double& currentSafety,
167                                                 G4GPILSelection* selection);
168
169  G4double PostStepGetPhysicalInteractionLength(const G4Track& track,
170                                                G4double   previousStepSize,
171                                                G4ForceCondition* condition);
172
[819]173  G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&);
174
175  G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
176
[1007]177  // Store PhysicsTable in a file.
178  // Return false in case of failure at I/O
[819]179  G4bool StorePhysicsTable(const G4ParticleDefinition*,
180                           const G4String& directory,
181                           G4bool ascii = false);
182
[1007]183  // Retrieve Physics from a file.
184  // (return true if the Physics Table can be build by using file)
185  // (return false if the process has no functionality or in case of failure)
186  // File name should is constructed as processName+particleName and the
187  // should be placed under the directory specifed by the argument.
[819]188  G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
189                              const G4String& directory,
190                              G4bool ascii);
191
[1007]192protected:
[819]193
[1007]194  G4double GetMeanFreePath(const G4Track& track,
195                           G4double previousStepSize,
196                           G4ForceCondition* condition);
[819]197
[1007]198  G4double GetContinuousStepLimit(const G4Track& track,
199                                  G4double previousStepSize,
200                                  G4double currentMinimumStep,
201                                  G4double& currentSafety);
202
[819]203  //------------------------------------------------------------------------
[1007]204  // Specific methods for along/post step EM processes
[819]205  //------------------------------------------------------------------------
206
207public:
[1007]208
209  void AddCollaborativeProcess(G4VEnergyLossProcess*);
210
211  void SampleSubCutSecondaries(std::vector<G4Track*>&, const G4Step&, 
212                               G4VEmModel* model, G4int matIdx,
213                               G4double& extraEdep); 
214
[819]215  G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple,
216                             const G4DynamicParticle* dp,
217                             G4double length);
218
[1007]219  //------------------------------------------------------------------------
220  // Specific methods to build and access Physics Tables
221  //------------------------------------------------------------------------
[819]222
[1007]223  G4PhysicsTable* BuildDEDXTable(G4EmTableType tType = fRestricted);
[819]224
[1007]225  G4PhysicsTable* BuildLambdaTable(G4EmTableType tType = fRestricted);
[819]226
[1007]227  void SetDEDXTable(G4PhysicsTable* p, G4EmTableType tType);
228  void SetCSDARangeTable(G4PhysicsTable* pRange);
229  void SetRangeTableForLoss(G4PhysicsTable* p);
230  void SetInverseRangeTable(G4PhysicsTable* p);
231  void SetSecondaryRangeTable(G4PhysicsTable* p);
[819]232
[1007]233  void SetLambdaTable(G4PhysicsTable* p);
234  void SetSubLambdaTable(G4PhysicsTable* p);
[819]235
[1007]236  // Binning for dEdx, range, inverse range and labda tables
237  inline void SetDEDXBinning(G4int nbins);
238  inline void SetLambdaBinning(G4int nbins);
[819]239
[1007]240  // Binning for dEdx, range, and inverse range tables
241  inline void SetDEDXBinningForCSDARange(G4int nbins);
[819]242
[1007]243  // Min kinetic energy for tables
244  inline void SetMinKinEnergy(G4double e);
245  inline G4double MinKinEnergy() const;
[819]246
[1007]247  // Max kinetic energy for tables
248  inline void SetMaxKinEnergy(G4double e);
249  inline G4double MaxKinEnergy() const;
[819]250
[1007]251  // Max kinetic energy for tables
252  inline void SetMaxKinEnergyForCSDARange(G4double e);
[819]253
[1007]254  // Access to specific tables
255  inline G4PhysicsTable* DEDXTable() const;
256  inline G4PhysicsTable* DEDXTableForSubsec() const;
257  inline G4PhysicsTable* DEDXunRestrictedTable() const;
258  inline G4PhysicsTable* IonisationTable() const;
259  inline G4PhysicsTable* IonisationTableForSubsec() const;
260  inline G4PhysicsTable* CSDARangeTable() const;
261  inline G4PhysicsTable* RangeTableForLoss() const;
262  inline G4PhysicsTable* InverseRangeTable() const;
263  inline G4PhysicsTable* LambdaTable();
264  inline G4PhysicsTable* SubLambdaTable();
[819]265
[1007]266  // Return values for given G4MaterialCutsCouple
267  inline G4double GetDEDX(G4double& kineticEnergy, const G4MaterialCutsCouple*);
268  inline G4double GetDEDXForSubsec(G4double& kineticEnergy, 
269                                   const G4MaterialCutsCouple*);
270  inline G4double GetRange(G4double& kineticEnergy, const G4MaterialCutsCouple*);
271  inline G4double GetCSDARange(G4double& kineticEnergy, const G4MaterialCutsCouple*);
272  inline G4double GetRangeForLoss(G4double& kineticEnergy, const G4MaterialCutsCouple*);
273  inline G4double GetKineticEnergy(G4double& range, const G4MaterialCutsCouple*);
274  inline G4double GetLambda(G4double& kineticEnergy, const G4MaterialCutsCouple*);
275
276  inline G4bool TablesAreBuilt() const;
277
[819]278  //------------------------------------------------------------------------
[1007]279  // Define and access particle type
[819]280  //------------------------------------------------------------------------
281
[1007]282  inline void SetBaseParticle(const G4ParticleDefinition* p);
283  inline const G4ParticleDefinition* Particle() const;
284  inline const G4ParticleDefinition* BaseParticle() const;
285  inline const G4ParticleDefinition* SecondaryParticle() const;
[819]286
[1007]287  //------------------------------------------------------------------------
288  // Specific methods to set, access, modify models
289  //------------------------------------------------------------------------
[819]290
[1007]291  // Add EM model coupled with fluctuation model for the region
[961]292  inline void AddEmModel(G4int, G4VEmModel*, 
293                         G4VEmFluctuationModel* fluc = 0,
294                         const G4Region* region = 0);
[819]295
296  // Assign a model to a process
297  inline void SetEmModel(G4VEmModel*, G4int index=1);
298 
299  // return the assigned model
300  inline G4VEmModel* EmModel(G4int index=1);
301 
302  // Assign a fluctuation model to a process
303  inline void SetFluctModel(G4VEmFluctuationModel*);
304 
305  // return the assigned fluctuation model
306  inline G4VEmFluctuationModel* FluctModel();
307   
[1007]308  // Define new energy range for the model identified by the name
309  inline void UpdateEmModel(const G4String&, G4double, G4double);
[819]310
[1007]311  // Access to models
312  inline G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false);
[819]313
[1007]314  inline G4int NumberOfModels();
[819]315
316  //------------------------------------------------------------------------
[1007]317  // Get/set parameters used for simulation of energy loss
[819]318  //------------------------------------------------------------------------
319
320  inline void SetLossFluctuations(G4bool val);
321  inline void SetRandomStep(G4bool val);
322  inline void SetIntegral(G4bool val);
323  inline G4bool IsIntegral() const;
324
325  // Set/Get flag "isIonisation"
326  inline void SetIonisation(G4bool val);
327  inline G4bool IsIonisationProcess() const;
328
329  // Redefine parameteters for stepping control
330  //
331  inline void SetLinearLossLimit(G4double val);
332  inline void SetMinSubRange(G4double val);
[1007]333  inline void SetStepFunction(G4double v1, G4double v2);
[1005]334  inline void SetLambdaFactor(G4double val);
[819]335
[1007]336
337  // Add subcutoff option for the region
338  void ActivateSubCutoff(G4bool val, const G4Region* region = 0);
339
[819]340  inline G4int NumberOfSubCutoffRegions() const;
341
[1007]342  // Activate deexcitation code
343  virtual void ActivateDeexcitation(G4bool, const G4Region* region = 0);
344
[819]345  //------------------------------------------------------------------------
[1007]346  // Public interface to helper functions
[819]347  //------------------------------------------------------------------------
348
[1007]349  inline 
350  G4VEmModel* SelectModelForMaterial(G4double kinEnergy, size_t& idx) const;
[819]351
[1007]352  inline G4double MeanFreePath(const G4Track& track);
[819]353
[1007]354  inline G4double ContinuousStepLimit(const G4Track& track,
355                                      G4double previousStepSize,
356                                      G4double currentMinimumStep,
357                                      G4double& currentSafety);
[819]358
[1007]359  //------------------------------------------------------------------------
360  // Run time method for simulation of ionisation
361  //------------------------------------------------------------------------
[819]362
[1007]363  // sample range at the end of a step
364  inline G4double SampleRange();
[819]365
[1007]366  // Set scaling parameters for ions is needed to G4EmCalculator
367  inline void SetDynamicMassCharge(G4double massratio, G4double charge2ratio);
[819]368
[1007]369  // Access to cross section table
370  G4double CrossSectionPerVolume(G4double kineticEnergy,
371                                 const G4MaterialCutsCouple* couple);
[819]372
[1007]373protected:
[819]374
[1007]375  G4PhysicsVector* LambdaPhysicsVector(const G4MaterialCutsCouple*, 
376                                       G4double cut);
[819]377
[1007]378  inline G4ParticleChangeForLoss* GetParticleChange();
[819]379
[1007]380  inline void SetParticle(const G4ParticleDefinition* p);
381
382  inline void SetSecondaryParticle(const G4ParticleDefinition* p);
383
384  inline void SelectModel(G4double kinEnergy);
385
386  inline size_t CurrentMaterialCutsCoupleIndex() const;
387
388  inline G4double GetCurrentRange() const;
389
390private:
391
[961]392  //------------------------------------------------------------------------
[1007]393  // Management of tables
[961]394  //------------------------------------------------------------------------
[819]395
[1007]396  void Clear();
[819]397
[1007]398  G4bool StoreTable(const G4ParticleDefinition* p, 
399                    G4PhysicsTable*, G4bool ascii,
400                    const G4String& directory, 
401                    const G4String& tname);
[819]402
[1007]403  G4bool RetrieveTable(const G4ParticleDefinition* p, 
404                       G4PhysicsTable*, G4bool ascii,
405                       const G4String& directory, 
406                       const G4String& tname, 
407                       G4bool mandatory);
[819]408
[961]409  // define material and indexes
410  inline void DefineMaterial(const G4MaterialCutsCouple* couple);
[819]411
[1007]412  // Returnd values for scaled energy using mass of the base particle
413  //
[819]414  inline G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy);
415  inline G4double GetSubDEDXForScaledEnergy(G4double scaledKinEnergy);
416  inline G4double GetIonisationForScaledEnergy(G4double scaledKinEnergy);
417  inline G4double GetSubIonisationForScaledEnergy(G4double scaledKinEnergy);
418  inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy);
419  inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinEnergy);
[1007]420  inline G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy);
[1005]421  inline G4double ScaledKinEnergyForLoss(G4double range);
[819]422  inline void ComputeLambdaForScaledEnergy(G4double scaledKinEnergy);
423
424  // hide  assignment operator
[1007]425
[819]426  G4VEnergyLossProcess(G4VEnergyLossProcess &);
427  G4VEnergyLossProcess & operator=(const G4VEnergyLossProcess &right);
428
[961]429  // ======== Parameters of the class fixed at construction =========
[819]430
[961]431  G4EmModelManager*           modelManager;
432  G4SafetyHelper*             safetyHelper;
[819]433
[961]434  const G4ParticleDefinition* secondaryParticle;
435  const G4ParticleDefinition* theElectron;
436  const G4ParticleDefinition* thePositron;
437  const G4ParticleDefinition* theGenericIon;
[819]438
[961]439  G4PhysicsVector*            vstrag;
[819]440
[961]441  // ======== Parameters of the class fixed at initialisation =======
442
[819]443  std::vector<G4VEmModel*>              emModels;
444  G4VEmFluctuationModel*                fluctModel;
445  std::vector<const G4Region*>          scoffRegions;
446  G4int                                 nSCoffRegions;
[1007]447  G4int*                                idxSCoffRegions;
[961]448
[819]449  std::vector<G4VEnergyLossProcess*>    scProcesses;
450  G4int                                 nProcesses;
451
452  // tables and vectors
453  G4PhysicsTable*             theDEDXTable;
454  G4PhysicsTable*             theDEDXSubTable;
455  G4PhysicsTable*             theDEDXunRestrictedTable;
456  G4PhysicsTable*             theIonisationTable;
457  G4PhysicsTable*             theIonisationSubTable;
458  G4PhysicsTable*             theRangeTableForLoss;
459  G4PhysicsTable*             theCSDARangeTable;
460  G4PhysicsTable*             theSecondaryRangeTable;
461  G4PhysicsTable*             theInverseRangeTable;
462  G4PhysicsTable*             theLambdaTable;
463  G4PhysicsTable*             theSubLambdaTable;
464  G4double*                   theDEDXAtMaxEnergy;
465  G4double*                   theRangeAtMaxEnergy;
466  G4double*                   theEnergyOfCrossSectionMax;
467  G4double*                   theCrossSectionMax;
468
469  const G4DataVector*         theCuts;
470  const G4DataVector*         theSubCuts;
471
472  const G4ParticleDefinition* baseParticle;
473
474  G4int    nBins;
475  G4int    nBinsCSDA;
476
477  G4double lowestKinEnergy;
478  G4double minKinEnergy;
479  G4double maxKinEnergy;
480  G4double maxKinEnergyCSDA;
481
482  G4double linLossLimit;
483  G4double minSubRange;
484  G4double dRoverRange;
485  G4double finalRange;
486  G4double lambdaFactor;
487
488  G4bool   lossFluctuationFlag;
489  G4bool   rndmStepFlag;
490  G4bool   tablesAreBuilt;
491  G4bool   integral;
[961]492  G4bool   isIon;
[819]493  G4bool   isIonisation;
494  G4bool   useSubCutoff;
495
[961]496protected:
[819]497
[961]498  G4ParticleChangeForLoss          fParticleChange;
[819]499
[961]500  // ======== Cashed values - may be state dependent ================
[819]501
[961]502private:
[819]503
[961]504  std::vector<G4DynamicParticle*>  secParticles;
505  std::vector<G4Track*>            scTracks;
[819]506
[961]507  const G4ParticleDefinition* particle;
[819]508
[961]509  G4VEmModel*                 currentModel;
510  const G4Material*           currentMaterial;
511  const G4MaterialCutsCouple* currentCouple;
512  size_t                      currentMaterialIndex;
513
514  G4int    nWarnings;
515
516  G4double massRatio;
517  G4double reduceFactor;
518  G4double chargeSqRatio;
519
520  G4double preStepLambda;
521  G4double fRange;
522  G4double preStepKinEnergy;
523  G4double preStepScaledEnergy;
524  G4double mfpKinEnergy;
525
526  G4GPILSelection  aGPILSelection;
527
528};
529
[819]530//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[961]531//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[819]532
[1007]533inline void G4VEnergyLossProcess::DefineMaterial(
534            const G4MaterialCutsCouple* couple)
[819]535{
[1007]536  if(couple != currentCouple) {
537    currentCouple   = couple;
538    currentMaterial = couple->GetMaterial();
539    currentMaterialIndex = couple->GetIndex();
540    mfpKinEnergy = DBL_MAX;
541  }
[819]542}
543
544//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
545
[1007]546inline G4double G4VEnergyLossProcess::GetDEDX(G4double& kineticEnergy,
547                                        const G4MaterialCutsCouple* couple)
[819]548{
[1007]549  DefineMaterial(couple);
550  return GetDEDXForScaledEnergy(kineticEnergy*massRatio);
[819]551}
552
553//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[1007]554
555inline G4double G4VEnergyLossProcess::GetDEDXForSubsec(G4double& kineticEnergy,
556                                        const G4MaterialCutsCouple* couple)
[819]557{
[1007]558  DefineMaterial(couple);
559  return GetSubDEDXForScaledEnergy(kineticEnergy*massRatio);
[819]560}
561
562//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
563
[1007]564inline G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e)
[819]565{
[1007]566  G4bool b;
567  G4double x = 
568    ((*theDEDXTable)[currentMaterialIndex]->GetValue(e, b))*chargeSqRatio;
569  if(e < minKinEnergy) x *= std::sqrt(e/minKinEnergy);
570  return x;
[819]571}
572
573//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
574
[1007]575inline G4double G4VEnergyLossProcess::GetSubDEDXForScaledEnergy(G4double e)
[819]576{
[1007]577  G4bool b;
578  G4double x = 
579    ((*theDEDXSubTable)[currentMaterialIndex]->GetValue(e, b))*chargeSqRatio;
580  if(e < minKinEnergy) x *= std::sqrt(e/minKinEnergy);
581  return x;
[819]582}
583
584//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
585
[1007]586inline G4double G4VEnergyLossProcess::GetIonisationForScaledEnergy(G4double e)
[819]587{
[1007]588  G4bool b;
589  G4double x = 0.0;
590  //  if(theIonisationTable) {
591  x = ((*theIonisationTable)[currentMaterialIndex]->GetValue(e, b))
592    *chargeSqRatio;
593  if(e < minKinEnergy) x *= std::sqrt(e/minKinEnergy);
594  //}
595  return x;
[819]596}
597
598//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
599
[1007]600inline 
601G4double G4VEnergyLossProcess::GetSubIonisationForScaledEnergy(G4double e)
[819]602{
[1007]603  G4bool b;
604  G4double x = 0.0;
605  //if(theIonisationSubTable) {
606  x = ((*theIonisationSubTable)[currentMaterialIndex]->GetValue(e, b))
607    *chargeSqRatio;
608  if(e < minKinEnergy) x *= std::sqrt(e/minKinEnergy);
609  //}
610  return x;
[819]611}
612
613//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
614
[1007]615inline G4double G4VEnergyLossProcess::GetRange(G4double& kineticEnergy,
616                                         const G4MaterialCutsCouple* couple)
[819]617{
[1007]618  G4double x = fRange;
619  if(kineticEnergy != preStepKinEnergy || couple != currentCouple) { 
620    DefineMaterial(couple);
621    if(theCSDARangeTable)
622      x = GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)
623        * reduceFactor;
624    else if(theRangeTableForLoss)
625      x = GetScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
626  }
627  return x;
[819]628}
629
630//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
631
[1007]632inline G4double G4VEnergyLossProcess::GetCSDARange(
633       G4double& kineticEnergy, const G4MaterialCutsCouple* couple)
[819]634{
[1007]635  DefineMaterial(couple);
636  G4double x = DBL_MAX;
637  if(theCSDARangeTable)
638    x = GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)
639      * reduceFactor;
640  return x;
[819]641}
642
643//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
644
[1007]645inline G4double G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(
646                G4double e)
[819]647{
[1007]648  G4bool b;
649  G4double x;
650
651  if (e < maxKinEnergyCSDA) {
652    x = ((*theCSDARangeTable)[currentMaterialIndex])->GetValue(e, b);
653    if(e < minKinEnergy) x *= std::sqrt(e/minKinEnergy);
654  } else {
655    x = theRangeAtMaxEnergy[currentMaterialIndex] +
656         (e - maxKinEnergyCSDA)/theDEDXAtMaxEnergy[currentMaterialIndex];
657  }
658  return x;
[819]659}
660
661//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
662
[1007]663inline G4double G4VEnergyLossProcess::GetRangeForLoss(
664                G4double& kineticEnergy,
665                const G4MaterialCutsCouple* couple)
[819]666{
[1007]667  DefineMaterial(couple);
668  G4double x = DBL_MAX;
669  if(theRangeTableForLoss) 
670    x = GetScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
671  //  G4cout << "Range from " << GetProcessName()
672  //         << "  e= " << kineticEnergy << " r= " << x << G4endl;
673  return x;
[819]674}
675
676//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
677
[1007]678inline G4double G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e)
[819]679{
[1007]680  G4bool b;
681  G4double x = ((*theRangeTableForLoss)[currentMaterialIndex])->GetValue(e, b);
682  if(e < minKinEnergy) x *= std::sqrt(e/minKinEnergy);
683  return x;
[819]684}
685
686//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
687
[1007]688inline G4double G4VEnergyLossProcess::GetKineticEnergy(
689                G4double& range,
690                const G4MaterialCutsCouple* couple)
[819]691{
[1007]692  DefineMaterial(couple);
693  G4double r = range/reduceFactor;
694  G4double e = ScaledKinEnergyForLoss(r)/massRatio;
695  return e;
[819]696}
697
698//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
699
[1007]700inline G4double G4VEnergyLossProcess::ScaledKinEnergyForLoss(G4double r)
[819]701{
[1007]702  G4PhysicsVector* v = (*theInverseRangeTable)[currentMaterialIndex];
703  G4double rmin = v->GetLowEdgeEnergy(0);
704  G4double e = 0.0; 
705  if(r >= rmin) {
706    G4bool b;
707    e = v->GetValue(r, b);
708  } else if(r > 0.0) {
709    G4double x = r/rmin;
710    e = minKinEnergy*x*x;
711  }
712  return e;
[819]713}
714
715//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
716
[1007]717inline G4double G4VEnergyLossProcess::GetLambda(G4double& kineticEnergy,
718                                          const G4MaterialCutsCouple* couple)
[819]719{
[1007]720  DefineMaterial(couple);
721  G4double x = 0.0;
722  if(theLambdaTable) x = GetLambdaForScaledEnergy(kineticEnergy*massRatio);
723  return x;
[819]724}
725
726//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
727
[1007]728inline G4double G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e)
[819]729{
[1007]730  G4bool b;
731  return 
732    chargeSqRatio*(((*theLambdaTable)[currentMaterialIndex])->GetValue(e, b));
[819]733}
734
735//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
736
[1007]737inline void G4VEnergyLossProcess::ComputeLambdaForScaledEnergy(G4double e)
[819]738{
[1007]739  mfpKinEnergy  = theEnergyOfCrossSectionMax[currentMaterialIndex];
740  if (e <= mfpKinEnergy) {
741    preStepLambda = GetLambdaForScaledEnergy(e);
742
743  } else {
744    G4double e1 = e*lambdaFactor;
745    if(e1 > mfpKinEnergy) {
746      preStepLambda  = GetLambdaForScaledEnergy(e);
747      G4double preStepLambda1 = GetLambdaForScaledEnergy(e1);
748      if(preStepLambda1 > preStepLambda) {
749        mfpKinEnergy = e1;
750        preStepLambda = preStepLambda1;
751      }
752    } else {
753      preStepLambda = chargeSqRatio*theCrossSectionMax[currentMaterialIndex];
754    }
755  }
[819]756}
757
758//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
759
[1007]760inline G4double G4VEnergyLossProcess::ContinuousStepLimit(
761         const G4Track& track, G4double x, G4double y, G4double& z)
[819]762{
[1007]763  G4GPILSelection sel;
764  return AlongStepGetPhysicalInteractionLength(track, x, y, z, &sel);
[819]765}
766
767//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
768
[1007]769inline G4double G4VEnergyLossProcess::SampleRange()
[819]770{
[1007]771  G4double e = amu_c2*preStepKinEnergy/particle->GetPDGMass();
772  G4bool b;
773  G4double s = fRange*std::pow(10.,vstrag->GetValue(e,b));
774  G4double x = fRange + G4RandGauss::shoot(0.0,s);
775  if(x > 0.0) fRange = x;
776  return fRange;
[819]777}
778
779//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
780
[1007]781inline G4double G4VEnergyLossProcess::MeanFreePath(const G4Track& track)
[819]782{
[1007]783  DefineMaterial(track.GetMaterialCutsCouple());
784  preStepLambda = GetLambdaForScaledEnergy(track.GetKineticEnergy()*massRatio);
785  G4double x = DBL_MAX;
786  if(DBL_MIN < preStepLambda) x = 1.0/preStepLambda;
787  return x;
[819]788}
789
790//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
791
[1007]792inline G4double G4VEnergyLossProcess::MinPrimaryEnergy(
793                const G4ParticleDefinition*, const G4Material*, G4double cut)
[819]794{
[1007]795  return cut;
[819]796}
797
798//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
799
[1007]800inline void G4VEnergyLossProcess::SelectModel(G4double kinEnergy)
[819]801{
[1007]802  currentModel = modelManager->SelectModel(kinEnergy, currentMaterialIndex);
[819]803}
804
805//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[1007]806
807inline G4VEmModel* G4VEnergyLossProcess::SelectModelForMaterial(
808                   G4double kinEnergy, size_t& idx) const
[819]809{
[1007]810  return modelManager->SelectModel(kinEnergy, idx);
[819]811}
812
813//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
814
[1007]815inline G4ParticleChangeForLoss* G4VEnergyLossProcess::GetParticleChange()
[819]816{
[1007]817  return &fParticleChange;
[819]818}
819
820//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
821
[1007]822inline const G4ParticleDefinition* G4VEnergyLossProcess::Particle() const
[819]823{
[1007]824  return particle;
[819]825}
826
827//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
828
[1007]829inline const G4ParticleDefinition* G4VEnergyLossProcess::BaseParticle() const
[961]830{
[1007]831  return baseParticle;
[961]832}
[819]833
834//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
835
[1007]836inline const G4ParticleDefinition* G4VEnergyLossProcess::SecondaryParticle() const
[819]837{
[1007]838  return secondaryParticle;
[819]839}
840
841//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
842
[1007]843inline G4PhysicsTable* G4VEnergyLossProcess::DEDXTable() const
[819]844{
[1007]845  return theDEDXTable;
[819]846}
847
848//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
849
[1007]850inline G4PhysicsTable* G4VEnergyLossProcess::DEDXTableForSubsec() const
[819]851{
[1007]852  return theDEDXSubTable;
[819]853}
854
855//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
856
[1007]857inline G4PhysicsTable* G4VEnergyLossProcess::DEDXunRestrictedTable() const
[819]858{
[1007]859  return theDEDXunRestrictedTable;
[819]860}
861
862//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
863
[1007]864inline G4PhysicsTable* G4VEnergyLossProcess::IonisationTable() const
[819]865{
[1007]866  G4PhysicsTable* t = theDEDXTable;
867  if(theIonisationTable) t = theIonisationTable; 
868  return t;
[819]869}
870
871//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
872
[1007]873inline G4PhysicsTable* G4VEnergyLossProcess::IonisationTableForSubsec() const
[819]874{
[1007]875  G4PhysicsTable* t = theDEDXSubTable;
876  if(theIonisationSubTable) t = theIonisationSubTable; 
877  return t;
[819]878}
879
880//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
881
[1007]882inline G4PhysicsTable* G4VEnergyLossProcess::CSDARangeTable() const
[819]883{
[1007]884  return theCSDARangeTable;
[819]885}
886
887//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
888
[1007]889inline G4PhysicsTable* G4VEnergyLossProcess::RangeTableForLoss() const
[819]890{
[1007]891  return theRangeTableForLoss;
[819]892}
893
894//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
895
[1007]896inline G4PhysicsTable* G4VEnergyLossProcess::InverseRangeTable() const
[819]897{
[1007]898  return theInverseRangeTable;
[819]899}
900
901//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
902
[1007]903inline G4PhysicsTable* G4VEnergyLossProcess::LambdaTable()
[819]904{
[1007]905  return theLambdaTable;
[819]906}
907
908//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[961]909
[1007]910inline G4PhysicsTable* G4VEnergyLossProcess::SubLambdaTable()
[819]911{
[1007]912  return theSubLambdaTable;
[819]913}
914
915//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[1007]916 
917inline G4bool G4VEnergyLossProcess::IsIntegral() const 
918{
919  return integral;
920}
[1005]921
[1007]922//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
923
924inline size_t G4VEnergyLossProcess::CurrentMaterialCutsCoupleIndex() const 
[991]925{
[1007]926  return currentMaterialIndex;
[991]927}
[819]928
[991]929//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
930
[1007]931inline void G4VEnergyLossProcess::SetDynamicMassCharge(G4double massratio,
932                                                       G4double charge2ratio)
[819]933{
[1007]934  massRatio     = massratio;
935  chargeSqRatio = charge2ratio;
936  reduceFactor  = 1.0/(chargeSqRatio*massRatio);
[819]937}
938
939//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[1007]940 
941inline G4double G4VEnergyLossProcess::GetCurrentRange() const
942{
943  return fRange;
944}
[819]945
[1007]946//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
947
948inline
949void G4VEnergyLossProcess::AddEmModel(G4int order, G4VEmModel* p, 
950                                      G4VEmFluctuationModel* fluc,
951                                      const G4Region* region)
[819]952{
[1007]953  modelManager->AddEmModel(order, p, fluc, region);
954  if(p) p->SetParticleChange(pParticleChange, fluc);
[819]955}
956
957//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[1005]958
[1007]959inline 
960G4VEmModel* G4VEnergyLossProcess::GetModelByIndex(G4int idx, G4bool ver)
[991]961{
[1007]962  return modelManager->GetModel(idx, ver);
[991]963}
[819]964
[991]965//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
966
[1007]967inline G4int G4VEnergyLossProcess::NumberOfModels()
[819]968{
[1007]969  return modelManager->NumberOfModels();
[819]970}
971
972//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
973
[1007]974inline void G4VEnergyLossProcess::SetEmModel(G4VEmModel* p, G4int index)
[819]975{
[1007]976  G4int n = emModels.size();
977  if(index >= n) for(G4int i=n; i<index+1; i++) {emModels.push_back(0);}
978  emModels[index] = p;
[819]979}
980
981//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
982
[1007]983inline G4VEmModel* G4VEnergyLossProcess::EmModel(G4int index)
[819]984{
[1007]985  G4VEmModel* p = 0;
986  if(index >= 0 && index <  G4int(emModels.size())) p = emModels[index];
987  return p;
[819]988}
989
990//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
991
[1007]992inline void G4VEnergyLossProcess::SetFluctModel(G4VEmFluctuationModel* p)
[819]993{
[1007]994  fluctModel = p;
[819]995}
996
997//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
998
[1007]999inline G4VEmFluctuationModel* G4VEnergyLossProcess::FluctModel()
[819]1000{
[1007]1001  return fluctModel;
[819]1002}
1003
1004//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1005
[1007]1006inline void G4VEnergyLossProcess::UpdateEmModel(const G4String& nam, 
1007                                                G4double emin, G4double emax)
[819]1008{
[1007]1009  modelManager->UpdateEmModel(nam, emin, emax);
[819]1010}
1011
1012//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1013
[1007]1014inline void G4VEnergyLossProcess::SetIntegral(G4bool val)
[819]1015{
[1007]1016  integral = val;
[819]1017}
1018
1019//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1020
[1007]1021inline void G4VEnergyLossProcess::SetParticle(const G4ParticleDefinition* p)
[819]1022{
[1007]1023  particle = p;
[819]1024}
1025
1026//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1027
[1007]1028inline void G4VEnergyLossProcess::SetBaseParticle(const G4ParticleDefinition* p)
[819]1029{
[1007]1030  baseParticle = p;
[819]1031}
1032
1033//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1034
[1007]1035inline void G4VEnergyLossProcess::SetSecondaryParticle(const G4ParticleDefinition* p)
[819]1036{
[1007]1037  secondaryParticle = p;
[819]1038}
1039
1040//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1041
[1007]1042inline void G4VEnergyLossProcess::SetLinearLossLimit(G4double val)
[819]1043{
[1007]1044  linLossLimit = val;
[819]1045}
1046
1047//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1048
[1007]1049inline void G4VEnergyLossProcess::SetLossFluctuations(G4bool val)
[819]1050{
[1007]1051  lossFluctuationFlag = val;
[819]1052}
1053
1054//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1055
[1007]1056inline void G4VEnergyLossProcess::SetRandomStep(G4bool val)
[819]1057{
[1007]1058  rndmStepFlag = val;
[819]1059}
1060
1061//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1062
[1007]1063inline void G4VEnergyLossProcess::SetMinSubRange(G4double val)
[819]1064{
[1007]1065  minSubRange = val;
[819]1066}
1067
1068//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1069
[1007]1070inline G4bool G4VEnergyLossProcess::TablesAreBuilt() const
[819]1071{
[1007]1072  return  tablesAreBuilt;
[819]1073}
1074
1075//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1076
[1007]1077inline G4int G4VEnergyLossProcess::NumberOfSubCutoffRegions() const
[819]1078{
[1007]1079  return nSCoffRegions;
[819]1080}
1081
1082//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1083
[1007]1084inline void G4VEnergyLossProcess::SetDEDXBinning(G4int nbins)
[819]1085{
[1007]1086  nBins = nbins;
[819]1087}
1088
1089//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1090
[1007]1091inline void G4VEnergyLossProcess::SetLambdaBinning(G4int nbins)
[819]1092{
[1007]1093  nBins = nbins;
[819]1094}
1095
1096//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1097
[1007]1098inline void G4VEnergyLossProcess::SetDEDXBinningForCSDARange(G4int nbins)
[819]1099{
[1007]1100  nBinsCSDA = nbins;
[819]1101}
1102
1103//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1104
[1007]1105inline G4double G4VEnergyLossProcess::MinKinEnergy() const
[819]1106{
[1007]1107  return minKinEnergy;
[819]1108}
1109
1110//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1111
[1007]1112inline void G4VEnergyLossProcess::SetMinKinEnergy(G4double e)
[819]1113{
[1007]1114  minKinEnergy = e;
[819]1115}
1116
1117//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1118
[1007]1119inline void G4VEnergyLossProcess::SetMaxKinEnergy(G4double e)
[819]1120{
[1007]1121  maxKinEnergy = e;
1122  if(e < maxKinEnergyCSDA) maxKinEnergyCSDA = e;
[819]1123}
1124
1125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1126
[1007]1127inline void G4VEnergyLossProcess::SetMaxKinEnergyForCSDARange(G4double e)
[819]1128{
[1007]1129  maxKinEnergyCSDA = e;
[819]1130}
1131
1132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1133
[1007]1134inline G4double G4VEnergyLossProcess::MaxKinEnergy() const
[819]1135{
[1007]1136  return maxKinEnergy;
[819]1137}
1138
1139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1140
[1007]1141inline void G4VEnergyLossProcess::SetLambdaFactor(G4double val)
[819]1142{
[1007]1143  if(val > 0.0 && val <= 1.0) lambdaFactor = val;
[819]1144}
1145
1146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1147
[1007]1148inline void G4VEnergyLossProcess::SetIonisation(G4bool val)
[819]1149{
[1007]1150  isIonisation = val;
1151  if(val) aGPILSelection = CandidateForSelection;
1152  else    aGPILSelection = NotCandidateForSelection;
[819]1153}
1154
1155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1156
[1007]1157inline G4bool G4VEnergyLossProcess::IsIonisationProcess() const
[819]1158{
[1007]1159  return isIonisation;
[819]1160}
1161
1162//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1163
[1007]1164void G4VEnergyLossProcess::SetStepFunction(G4double v1, G4double v2)
[819]1165{
[1007]1166  dRoverRange = v1;
1167  finalRange = v2;
1168  if (dRoverRange > 0.999) dRoverRange = 1.0;
1169  currentCouple = 0;
1170  mfpKinEnergy  = DBL_MAX;
[819]1171}
1172
1173//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1174
1175#endif
Note: See TracBrowser for help on using the repository browser.