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

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

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

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