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

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

update processes

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