source: trunk/source/processes/electromagnetic/utils/include/G4VEmProcess.hh @ 966

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

update processes

File size: 21.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//
[961]26// $Id: G4VEmProcess.hh,v 1.50 2009/02/19 09:57:36 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
[819]28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name:     G4VEmProcess
35//
36// Author:        Vladimir Ivanchenko
37//
38// Creation date: 01.10.2003
39//
40// Modifications:
41// 30-06-04 make destructor virtual (V.Ivanchenko)
42// 09-08-04 optimise integral option (V.Ivanchenko)
43// 11-08-04 add protected methods to access cuts (V.Ivanchenko)
44// 09-09-04 Bug fix for the integral mode with 2 peaks (V.Ivanchneko)
45// 16-09-04 Add flag for LambdaTable and method RecalculateLambda (VI)
46// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
47// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
48// 18-04-05 Use G4ParticleChangeForGamma (V.Ivantchenko)
49// 09-05-05 Fix problem in logic when path boundary between materials (VI)
50// 11-01-06 add A to parameters of ComputeCrossSectionPerAtom (VI)
51// 01-02-06 put default value A=0. to keep compatibility with v5.2 (mma)
52// 13-05-06 Add method to access model by index (V.Ivanchenko)
53// 12-09-06 add SetModel() (mma)
54// 25-09-07 More accurate handling zero xsect in
55//          PostStepGetPhysicalInteractionLength (V.Ivanchenko)
56// 27-10-07 Virtual functions moved to source (V.Ivanchenko)
[961]57// 15-07-08 Reorder class members for further multi-thread development (VI)
[819]58//
59// Class Description:
60//
61// It is the unified Discrete process
62
63// -------------------------------------------------------------------
64//
65
66#ifndef G4VEmProcess_h
67#define G4VEmProcess_h 1
68
69#include "G4VDiscreteProcess.hh"
70#include "globals.hh"
71#include "G4Material.hh"
72#include "G4MaterialCutsCouple.hh"
73#include "G4Track.hh"
74#include "G4EmModelManager.hh"
75#include "G4UnitsTable.hh"
76#include "G4ParticleDefinition.hh"
77#include "G4ParticleChangeForGamma.hh"
78
79class G4Step;
80class G4VEmModel;
81class G4DataVector;
82class G4VParticleChange;
83class G4PhysicsTable;
84class G4PhysicsVector;
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
87
88class G4VEmProcess : public G4VDiscreteProcess
89{
90public:
91
92  G4VEmProcess(const G4String& name,
[961]93               G4ProcessType type = fElectromagnetic);
[819]94
95  virtual ~G4VEmProcess();
96
97  //------------------------------------------------------------------------
98  // Virtual methods to be implemented in concrete processes
99  //------------------------------------------------------------------------
100
101  virtual G4bool IsApplicable(const G4ParticleDefinition& p) = 0;
102
103  virtual void PrintInfo() = 0;
104
105protected:
106
107  virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
108
109  //------------------------------------------------------------------------
[961]110  // Implementation of virtual methods common to all Discrete processes
[819]111  //------------------------------------------------------------------------
112
[961]113public:
[819]114
[961]115  // Initialise for build of tables
116  void PreparePhysicsTable(const G4ParticleDefinition&);
[819]117
[961]118  // Build physics table during initialisation
119  void BuildPhysicsTable(const G4ParticleDefinition&);
[819]120
121  void PrintInfoDefinition();
122
[961]123  // implementation of virtual method, specific for G4VEmProcess
124  G4double PostStepGetPhysicalInteractionLength(
125                             const G4Track& track,
126                             G4double   previousStepSize,
127                             G4ForceCondition* condition
128                            );
129
130  // implementation of virtual method, specific for G4VEmProcess
[819]131  G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
132
[961]133  // Store PhysicsTable in a file.
134  // Return false in case of failure at I/O
[819]135  G4bool StorePhysicsTable(const G4ParticleDefinition*,
136                           const G4String& directory,
137                           G4bool ascii = false);
138
[961]139  // Retrieve Physics from a file.
140  // (return true if the Physics Table can be build by using file)
141  // (return false if the process has no functionality or in case of failure)
142  // File name should is constructed as processName+particleName and the
143  // should be placed under the directory specifed by the argument.
[819]144  G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
145                              const G4String& directory,
146                              G4bool ascii);
147
[961]148  // deexcitation activated per G4Region
149  void ActivateDeexcitation(G4bool, const G4Region* r = 0);
150
[819]151  //------------------------------------------------------------------------
152  // Specific methods for Discrete EM post step simulation
153  //------------------------------------------------------------------------
154
[961]155  // It returns the cross section per volume for energy/ material
156  G4double CrossSectionPerVolume(G4double kineticEnergy,
157                                 const G4MaterialCutsCouple* couple);
[819]158
[961]159  // It returns the cross section of the process per atom
[819]160  inline G4double ComputeCrossSectionPerAtom(G4double kineticEnergy, 
161                                             G4double Z, G4double A=0., 
162                                             G4double cut=0.0);
163
164  inline G4double MeanFreePath(const G4Track& track);
165
[961]166  // It returns cross section per volume
[819]167  inline G4double GetLambda(G4double& kinEnergy, 
168                            const G4MaterialCutsCouple* couple);
169
170  //------------------------------------------------------------------------
171  // Specific methods to build and access Physics Tables
172  //------------------------------------------------------------------------
173
[961]174  // Binning for lambda table
[819]175  inline void SetLambdaBinning(G4int nbins);
176  inline G4int LambdaBinning() const;
177
[961]178  // Min kinetic energy for tables
[819]179  inline void SetMinKinEnergy(G4double e);
180  inline G4double MinKinEnergy() const;
181
[961]182  // Max kinetic energy for tables
[819]183  inline void SetMaxKinEnergy(G4double e);
184  inline G4double MaxKinEnergy() const;
185
[961]186  inline void SetPolarAngleLimit(G4double a);
187  inline G4double PolarAngleLimit() const;
188
[819]189  inline const G4PhysicsTable* LambdaTable() const;
190
191  //------------------------------------------------------------------------
192  // Define and access particle type
193  //------------------------------------------------------------------------
194
195  inline const G4ParticleDefinition* Particle() const;
196  inline const G4ParticleDefinition* SecondaryParticle() const;
197
198  //------------------------------------------------------------------------
[961]199  // Specific methods to set, access, modify models and basic parameters
[819]200  //------------------------------------------------------------------------
201
[961]202protected:
203  // Select model in run time
204  inline void SelectModel(G4double& kinEnergy);
205
206public:
207  // Select model by energy and region index
208  inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy, 
209                                            size_t& idxRegion) const;
210   
211  // Add model for region, smaller value of order defines which
212  // model will be selected for a given energy interval 
[819]213  inline void AddEmModel(G4int, G4VEmModel*, const G4Region* region = 0);
[961]214
[819]215  // Assign a model to a process
[961]216  inline void SetModel(G4VEmModel*, G4int index = 1);
[819]217 
218  // return the assigned model
[961]219  inline G4VEmModel* Model(G4int index = 1);
[819]220   
[961]221  // Define new energy range for the model identified by the name
[819]222  inline void UpdateEmModel(const G4String&, G4double, G4double);
223
224  // Access to models
[961]225  inline G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false);
[819]226
227  inline void SetLambdaFactor(G4double val);
228
229  inline void SetIntegral(G4bool val);
230  inline G4bool IsIntegral() const;
231
232  inline void SetApplyCuts(G4bool val);
[961]233
234  //------------------------------------------------------------------------
235  // Other generic methods
236  //------------------------------------------------------------------------
[819]237 
238protected:
239
240  G4double GetMeanFreePath(const G4Track& track,
241                           G4double previousStepSize,
242                           G4ForceCondition* condition);
243
244  G4PhysicsVector* LambdaPhysicsVector(const G4MaterialCutsCouple*);
245
[961]246  inline G4double RecalculateLambda(G4double kinEnergy,
247                                    const G4MaterialCutsCouple* couple);
248
249  inline G4ParticleChangeForGamma* GetParticleChange();
250
[819]251  inline void SetParticle(const G4ParticleDefinition* p);
252 
253  inline void SetSecondaryParticle(const G4ParticleDefinition* p);
254
255  inline size_t CurrentMaterialCutsCoupleIndex() const;
256
257  inline G4double GetGammaEnergyCut();
258
259  inline G4double GetElectronEnergyCut();
260
261  inline void SetBuildTableFlag(G4bool val);
262
263  inline void SetStartFromNullFlag(G4bool val);
264
265private:
266
267  void Clear();
268
269  void BuildLambdaTable();
270
271  void FindLambdaMax();
272
273  inline void InitialiseStep(const G4Track&);
274
275  inline void DefineMaterial(const G4MaterialCutsCouple* couple);
276
277  inline void ComputeIntegralLambda(G4double kinEnergy);
278
279  inline G4double GetLambdaFromTable(G4double kinEnergy);
280
281  inline G4double GetCurrentLambda(G4double kinEnergy);
282
283  inline G4double ComputeCurrentLambda(G4double kinEnergy);
284
[961]285  // copy constructor and hide assignment operator
[819]286  G4VEmProcess(G4VEmProcess &);
287  G4VEmProcess & operator=(const G4VEmProcess &right);
288
[961]289  // ======== Parameters of the class fixed at construction =========
[819]290
[961]291  G4EmModelManager*            modelManager;
292  const G4ParticleDefinition*  theGamma;
293  const G4ParticleDefinition*  theElectron;
294  const G4ParticleDefinition*  thePositron;
295  const G4ParticleDefinition*  secondaryParticle;
[819]296
[961]297  G4bool                       buildLambdaTable;
[819]298
[961]299  // ======== Parameters of the class fixed at initialisation =======
[819]300
[961]301  std::vector<G4VEmModel*>     emModels;
[819]302
303  // tables and vectors
304  G4PhysicsTable*              theLambdaTable;
305  G4double*                    theEnergyOfCrossSectionMax;
306  G4double*                    theCrossSectionMax;
307
308  const std::vector<G4double>* theCuts;
309  const std::vector<G4double>* theCutsGamma;
310  const std::vector<G4double>* theCutsElectron;
311  const std::vector<G4double>* theCutsPositron;
312
313  G4int                        nLambdaBins;
314
315  G4double                     minKinEnergy;
316  G4double                     maxKinEnergy;
317  G4double                     lambdaFactor;
[961]318  G4double                     polarAngleLimit;
[819]319
[961]320  G4bool                       integral;
321  G4bool                       applyCuts;
322  G4bool                       startFromNull;
323  G4bool                       useDeexcitation;
324
325  G4int                        nDERegions;
326  std::vector<const G4Region*> deRegions;
327  G4bool*                      idxDERegions;
328
329  // ======== Cashed values - may be state dependent ================
330
331protected:
332
333  G4ParticleChangeForGamma     fParticleChange;
334
335private:
336
337  std::vector<G4DynamicParticle*> secParticles;
338
339  G4VEmModel*                  currentModel; 
340
341  const G4ParticleDefinition*  particle;
342
[819]343  // cash
344  const G4Material*            currentMaterial;
345  const G4MaterialCutsCouple*  currentCouple;
346  size_t                       currentMaterialIndex;
347
348  G4double                     mfpKinEnergy;
349  G4double                     preStepKinEnergy;
350  G4double                     preStepLambda;
351
352};
353
354//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
355//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
356
[961]357inline G4double G4VEmProcess::ComputeCrossSectionPerAtom(
358                 G4double kineticEnergy, G4double Z, G4double A, G4double cut)
[819]359{
[961]360  SelectModel(kineticEnergy);
361  G4double x = 0.0;
362  if(currentModel) {
363   x = currentModel->ComputeCrossSectionPerAtom(particle,kineticEnergy,
364                                                 Z,A,cut);
[819]365  }
[961]366  return x;
[819]367}
368
369//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
370
[961]371inline G4double G4VEmProcess::MeanFreePath(const G4Track& track)
[819]372{
373  DefineMaterial(track.GetMaterialCutsCouple());
[961]374  preStepLambda = GetCurrentLambda(track.GetKineticEnergy());
375  G4double x = DBL_MAX;
376  if(DBL_MIN < preStepLambda) x = 1.0/preStepLambda;
377  return x;
[819]378}
379
380//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
381
382inline G4double G4VEmProcess::GetLambda(G4double& kineticEnergy,
[961]383                                        const G4MaterialCutsCouple* couple)
[819]384{
385  DefineMaterial(couple);
386  return GetCurrentLambda(kineticEnergy);
387}
388
389//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
390
[961]391inline void G4VEmProcess::SetLambdaBinning(G4int nbins)
[819]392{
[961]393  nLambdaBins = nbins;
[819]394}
395
396//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
397
[961]398inline G4int G4VEmProcess::LambdaBinning() const
[819]399{
[961]400  return nLambdaBins;
[819]401}
402
403//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
404
[961]405inline void G4VEmProcess::SetMinKinEnergy(G4double e)
[819]406{
[961]407  minKinEnergy = e;
[819]408}
409
410//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
411
[961]412inline G4double G4VEmProcess::MinKinEnergy() const
[819]413{
[961]414  return minKinEnergy;
[819]415}
416
417//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
418
[961]419inline void G4VEmProcess::SetMaxKinEnergy(G4double e)
[819]420{
[961]421  maxKinEnergy = e;
422}
[819]423
[961]424//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
425
426inline G4double G4VEmProcess::MaxKinEnergy() const
427{
428  return maxKinEnergy;
[819]429}
430
431//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
432
[961]433inline void G4VEmProcess::SetPolarAngleLimit(G4double val)
[819]434{
[961]435  if(val < 0.0)     polarAngleLimit = 0.0;
436  else if(val > pi) polarAngleLimit = pi;
437  else              polarAngleLimit = val;
[819]438}
439
440//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
441
[961]442inline G4double G4VEmProcess::PolarAngleLimit() const
[819]443{
[961]444  return polarAngleLimit;
[819]445}
446
447//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
448
[961]449inline const G4PhysicsTable* G4VEmProcess::LambdaTable() const
[819]450{
[961]451  return theLambdaTable;
[819]452}
453
454//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
455
456inline const G4ParticleDefinition* G4VEmProcess::Particle() const
457{
458  return particle;
459}
460
461//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
462
463inline const G4ParticleDefinition* G4VEmProcess::SecondaryParticle() const
464{
465  return secondaryParticle;
466}
467
468//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
469
[961]470inline void G4VEmProcess::SelectModel(G4double& kinEnergy)
[819]471{
[961]472  currentModel = modelManager->SelectModel(kinEnergy, currentMaterialIndex);
473  currentModel->SetCurrentCouple(currentCouple);
[819]474}
475
476//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
477
[961]478inline G4VEmModel* G4VEmProcess::SelectModelForMaterial(
479                                   G4double kinEnergy, size_t& idxRegion) const
[819]480{
[961]481  return modelManager->SelectModel(kinEnergy, idxRegion);
[819]482}
483
484//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
485
[961]486inline void G4VEmProcess::AddEmModel(G4int order, G4VEmModel* p, 
487                                     const G4Region* region)
[819]488{
[961]489  G4VEmFluctuationModel* fm = 0;
490  modelManager->AddEmModel(order, p, fm, region);
491  if(p) p->SetParticleChange(pParticleChange);
[819]492}
493
494//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
495
[961]496inline void G4VEmProcess::SetModel(G4VEmModel* p, G4int index)
[819]497{
[961]498  G4int n = emModels.size();
499  if(index >= n) for(G4int i=n; i<index+1; i++) {emModels.push_back(0);}
500  emModels[index] = p;
[819]501}
502
503//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
504
[961]505inline G4VEmModel* G4VEmProcess::Model(G4int index)
[819]506{
[961]507  G4VEmModel* p = 0;
508  if(index >= 0 && index <  G4int(emModels.size())) p = emModels[index];
509  return p;
[819]510}
511
512//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
513
[961]514inline void G4VEmProcess::UpdateEmModel(const G4String& nam, 
515                                        G4double emin, G4double emax)
[819]516{
[961]517  modelManager->UpdateEmModel(nam, emin, emax);
[819]518}
519
520//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
521
[961]522inline G4VEmModel* G4VEmProcess::GetModelByIndex(G4int idx, G4bool ver)
[819]523{
[961]524  return modelManager->GetModel(idx, ver);
[819]525}
526
527//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
528
[961]529inline void G4VEmProcess::SetLambdaFactor(G4double val)
[819]530{
[961]531  if(val > 0.0 && val <= 1.0) lambdaFactor = val;
[819]532}
533
534//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
535
[961]536inline void G4VEmProcess::SetIntegral(G4bool val)
[819]537{
[961]538  if(particle && particle != theGamma) integral = val;
539  if(integral) buildLambdaTable = true;
[819]540}
541
542//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
543
[961]544inline G4bool G4VEmProcess::IsIntegral() const
[819]545{
[961]546  return integral;
[819]547}
548
549//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
550
[961]551inline void G4VEmProcess::SetApplyCuts(G4bool val)
[819]552{
[961]553  applyCuts = val;
[819]554}
555
556//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
557
[961]558inline G4double G4VEmProcess::RecalculateLambda(G4double e, 
559                                                const G4MaterialCutsCouple* couple)
[819]560{
[961]561  DefineMaterial(couple);
562  return ComputeCurrentLambda(e);
[819]563}
564
565//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
566
[961]567inline G4ParticleChangeForGamma* G4VEmProcess::GetParticleChange()
[819]568{
[961]569  return &fParticleChange;
[819]570}
571
572//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
573
[961]574inline void G4VEmProcess::SetParticle(const G4ParticleDefinition* p)
[819]575{
[961]576  particle = p;
[819]577}
578
579//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
580
[961]581inline void G4VEmProcess::SetSecondaryParticle(const G4ParticleDefinition* p)
[819]582{
[961]583  secondaryParticle = p;
[819]584}
585
586//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
587
[961]588inline size_t G4VEmProcess::CurrentMaterialCutsCoupleIndex() const 
[819]589{
[961]590  return currentMaterialIndex;
[819]591}
592
593//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
594
[961]595inline G4double G4VEmProcess::GetGammaEnergyCut()
[819]596{
[961]597  return (*theCutsGamma)[currentMaterialIndex];
[819]598}
599
600//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
601
[961]602inline G4double G4VEmProcess::GetElectronEnergyCut()
603{
604  return (*theCutsElectron)[currentMaterialIndex];
605}
[819]606
607//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
608
[961]609inline void G4VEmProcess::SetBuildTableFlag(G4bool val)
[819]610{
[961]611  buildLambdaTable = val;
612  if(!val) integral = false;
[819]613}
614
615//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
616
[961]617inline void G4VEmProcess::SetStartFromNullFlag(G4bool val)
[819]618{
[961]619  startFromNull = val;
[819]620}
621
622//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
623
[961]624inline void G4VEmProcess::InitialiseStep(const G4Track& track)
[819]625{
[961]626  preStepKinEnergy = track.GetKineticEnergy();
627  DefineMaterial(track.GetMaterialCutsCouple());
628  if (theNumberOfInteractionLengthLeft < 0.0) mfpKinEnergy = DBL_MAX;
[819]629}
630
631//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
632
[961]633inline void G4VEmProcess::DefineMaterial(const G4MaterialCutsCouple* couple)
[819]634{
[961]635  if(couple != currentCouple) {
636    currentCouple   = couple;
637    currentMaterial = couple->GetMaterial();
638    currentMaterialIndex = couple->GetIndex();
639    mfpKinEnergy = DBL_MAX;
640  }
[819]641}
642
643//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
644
[961]645inline void G4VEmProcess::ComputeIntegralLambda(G4double e)
[819]646{
[961]647  mfpKinEnergy  = theEnergyOfCrossSectionMax[currentMaterialIndex];
648  if (e <= mfpKinEnergy) {
649    preStepLambda = GetLambdaFromTable(e);
650
651  } else {
652    G4double e1 = e*lambdaFactor;
653    if(e1 > mfpKinEnergy) {
654      preStepLambda  = GetLambdaFromTable(e);
655      G4double preStepLambda1 = GetLambdaFromTable(e1);
656      if(preStepLambda1 > preStepLambda) {
657        mfpKinEnergy = e1;
658        preStepLambda = preStepLambda1;
659      }
660    } else {
661      preStepLambda = theCrossSectionMax[currentMaterialIndex];
662    }
663  }
[819]664}
665
666//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
667
[961]668inline G4double G4VEmProcess::GetLambdaFromTable(G4double e)
[819]669{
[961]670  G4bool b;
671  return (((*theLambdaTable)[currentMaterialIndex])->GetValue(e, b));
[819]672}
673
674//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
675
[961]676inline G4double G4VEmProcess::GetCurrentLambda(G4double e)
[819]677{
[961]678  G4double x = 0.0;
679  if(theLambdaTable) x = GetLambdaFromTable(e);
680  else               x = ComputeCurrentLambda(e);
681  return x;
[819]682}
683
684//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
685
[961]686inline G4double G4VEmProcess::ComputeCurrentLambda(G4double e)
687{
688  SelectModel(e);
689  G4double x = 0.0;
690  if(currentModel) {
691    x = currentModel->CrossSectionPerVolume(currentMaterial,particle,
692                                            e,(*theCuts)[currentMaterialIndex]);
693  }
694  return x;
695}
696
697//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
698
[819]699#endif
Note: See TracBrowser for help on using the repository browser.