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