source: trunk/source/processes/electromagnetic/utils/include/G4VMultipleScattering.hh@ 1199

Last change on this file since 1199 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 18.3 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: G4VMultipleScattering.hh,v 1.62 2009/10/29 17:56:04 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4VMultipleScattering
35//
36// Author: Vladimir Ivanchenko on base of Laszlo Urban code
37//
38// Creation date: 12.03.2002
39//
40// Modifications:
41//
42// 16-07-03 Update GetRange interface (V.Ivanchenko)
43//
44//
45// Class Description:
46//
47// It is the generic process of multiple scattering it includes common
48// part of calculations for all charged particles
49//
50// 26-11-03 bugfix in AlongStepDoIt (L.Urban)
51// 25-05-04 add protection against case when range is less than steplimit (VI)
52// 30-06-04 make destructor virtual (V.Ivanchenko)
53// 27-08-04 Add InitialiseForRun method (V.Ivanchneko)
54// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko)
55// 15-04-05 optimize internal interfaces (V.Ivanchenko)
56// 15-04-05 remove boundary flag (V.Ivanchenko)
57// 07-10-05 error in a protection in GetContinuousStepLimit corrected (L.Urban)
58// 27-10-05 introduce virtual function MscStepLimitation() (V.Ivanchenko)
59// 26-01-06 Rename GetRange -> GetRangeFromRestricteDEDX (V.Ivanchenko)
60// 17-02-06 Save table of transport cross sections not mfp (V.Ivanchenko)
61// 07-03-06 Move step limit calculation to model (V.Ivanchenko)
62// 13-05-06 Add method to access model by index (V.Ivanchenko)
63// 12-02-07 Add get/set skin (V.Ivanchenko)
64// 27-10-07 Virtual functions moved to source (V.Ivanchenko)
65// 15-07-08 Reorder class members for further multi-thread development (VI)
66// 07-04-09 Moved msc methods from G4VEmModel to G4VMscModel (VI)
67//
68
69// -------------------------------------------------------------------
70//
71
72#ifndef G4VMultipleScattering_h
73#define G4VMultipleScattering_h 1
74
75#include "G4VContinuousDiscreteProcess.hh"
76#include "globals.hh"
77#include "G4Material.hh"
78#include "G4MaterialCutsCouple.hh"
79#include "G4ParticleChangeForMSC.hh"
80#include "G4Track.hh"
81#include "G4Step.hh"
82#include "G4EmModelManager.hh"
83#include "G4VMscModel.hh"
84#include "G4MscStepLimitType.hh"
85
86class G4ParticleDefinition;
87class G4DataVector;
88class G4PhysicsTable;
89class G4PhysicsVector;
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92
93class G4VMultipleScattering : public G4VContinuousDiscreteProcess
94{
95public:
96
97 G4VMultipleScattering(const G4String& name = "msc",
98 G4ProcessType type = fElectromagnetic);
99
100 virtual ~G4VMultipleScattering();
101
102 //------------------------------------------------------------------------
103 // Virtual methods to be implemented for the concrete model
104 //------------------------------------------------------------------------
105
106 virtual G4bool IsApplicable(const G4ParticleDefinition& p) = 0;
107
108 virtual void PrintInfo() = 0;
109
110protected:
111
112 virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
113
114public:
115
116 //------------------------------------------------------------------------
117 // Generic methods common to all ContinuousDiscrete processes
118 //------------------------------------------------------------------------
119
120 // Initialise for build of tables
121 void PreparePhysicsTable(const G4ParticleDefinition&);
122
123 // Build physics table during initialisation
124 void BuildPhysicsTable(const G4ParticleDefinition&);
125
126 // Print out of generic class parameters
127 void PrintInfoDefinition();
128
129 // Store PhysicsTable in a file.
130 // Return false in case of failure at I/O
131 G4bool StorePhysicsTable(const G4ParticleDefinition*,
132 const G4String& directory,
133 G4bool ascii = false);
134
135 // Retrieve Physics from a file.
136 // (return true if the Physics Table can be build by using file)
137 // (return false if the process has no functionality or in case of failure)
138 // File name should is constructed as processName+particleName and the
139 // should be placed under the directory specifed by the argument.
140 G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
141 const G4String& directory,
142 G4bool ascii);
143
144 // The function overloads the corresponding function of the base
145 // class.It limits the step near to boundaries only
146 // and invokes the method GetMscContinuousStepLimit at every step.
147 G4double AlongStepGetPhysicalInteractionLength(
148 const G4Track&,
149 G4double previousStepSize,
150 G4double currentMinimalStep,
151 G4double& currentSafety,
152 G4GPILSelection* selection);
153
154 // The function overloads the corresponding function of the base
155 // class.
156 inline G4double PostStepGetPhysicalInteractionLength(
157 const G4Track&,
158 G4double previousStepSize,
159 G4ForceCondition* condition);
160
161 // Along step actions
162 inline G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&);
163
164 // Post step actions
165 inline G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
166
167 // This method does not used for tracking, it is intended only for tests
168 inline G4double ContinuousStepLimit(const G4Track& track,
169 G4double previousStepSize,
170 G4double currentMinimalStep,
171 G4double& currentSafety);
172
173 //------------------------------------------------------------------------
174 // Specific methods to build and access Physics Tables
175 //------------------------------------------------------------------------
176
177 // Build empty Physics Vector
178 G4PhysicsVector* PhysicsVector(const G4MaterialCutsCouple*);
179
180 inline void SetBinning(G4int nbins);
181 inline G4int Binning() const;
182
183 inline void SetMinKinEnergy(G4double e);
184 inline G4double MinKinEnergy() const;
185
186 inline void SetMaxKinEnergy(G4double e);
187 inline G4double MaxKinEnergy() const;
188
189 inline void SetBuildLambdaTable(G4bool val);
190
191 inline G4PhysicsTable* LambdaTable() const;
192
193 // access particle type
194 inline const G4ParticleDefinition* Particle() const;
195
196 //------------------------------------------------------------------------
197 // Specific methods to set, access, modify models
198 //------------------------------------------------------------------------
199
200protected:
201 // Select model in run time
202 inline G4VEmModel* SelectModel(G4double kinEnergy);
203
204public:
205 // Select model in run time
206 inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
207 size_t& idxRegion) const;
208
209 // Add model for region, smaller value of order defines which
210 // model will be selected for a given energy interval
211 void AddEmModel(G4int order, G4VEmModel*, const G4Region* region = 0);
212
213 // Assign a model to a process
214 void SetModel(G4VMscModel*, G4int index = 1);
215
216 // return the assigned model
217 G4VMscModel* Model(G4int index = 1);
218
219 // Access to models by index
220 G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
221
222 //------------------------------------------------------------------------
223 // Get/Set parameters for simulation of multiple scattering
224 //------------------------------------------------------------------------
225
226 inline G4bool LateralDisplasmentFlag() const;
227 inline void SetLateralDisplasmentFlag(G4bool val);
228
229 inline G4double Skin() const;
230 inline void SetSkin(G4double val);
231
232 inline G4double RangeFactor() const;
233 inline void SetRangeFactor(G4double val);
234
235 inline G4double GeomFactor() const;
236 inline void SetGeomFactor(G4double val);
237
238 inline G4double PolarAngleLimit() const;
239 inline void SetPolarAngleLimit(G4double val);
240
241 inline G4MscStepLimitType StepLimitType() const;
242 inline void SetStepLimitType(G4MscStepLimitType val);
243
244 //------------------------------------------------------------------------
245 // Run time methods
246 //------------------------------------------------------------------------
247
248protected:
249
250 // This method is not used for tracking, it returns mean free path value
251 G4double GetMeanFreePath(const G4Track& track,
252 G4double,
253 G4ForceCondition* condition);
254
255 // This method is not used for tracking, it returns step limit
256 G4double GetContinuousStepLimit(const G4Track& track,
257 G4double previousStepSize,
258 G4double currentMinimalStep,
259 G4double& currentSafety);
260
261 // This method returns inversed transport cross section
262 inline G4double GetLambda(const G4ParticleDefinition* p,
263 G4double& kineticEnergy);
264
265 // defines current material in run time
266 inline void DefineMaterial(const G4MaterialCutsCouple* couple);
267
268 inline const G4MaterialCutsCouple* CurrentMaterialCutsCouple() const;
269
270private:
271
272 // hide assignment operator
273 G4VMultipleScattering(G4VMultipleScattering &);
274 G4VMultipleScattering & operator=(const G4VMultipleScattering &right);
275
276 // ======== Parameters of the class fixed at construction =========
277
278 G4EmModelManager* modelManager;
279 G4bool buildLambdaTable;
280
281 // ======== Parameters of the class fixed at initialisation =======
282
283 std::vector<G4VMscModel*> mscModels;
284
285 G4PhysicsTable* theLambdaTable;
286 const G4ParticleDefinition* firstParticle;
287
288 G4MscStepLimitType stepLimit;
289
290 G4double minKinEnergy;
291 G4double maxKinEnergy;
292 G4double skin;
293 G4double facrange;
294 G4double facgeom;
295 G4double polarAngleLimit;
296
297 G4int nBins;
298
299 G4bool latDisplasment;
300 G4bool isIon;
301
302 // ======== Cashed values - may be state dependent ================
303
304protected:
305
306 G4GPILSelection valueGPILSelectionMSC;
307 G4ParticleChangeForMSC fParticleChange;
308
309private:
310
311 G4VMscModel* currentModel;
312
313 // cache
314 const G4ParticleDefinition* currentParticle;
315 const G4MaterialCutsCouple* currentCouple;
316 size_t currentMaterialIndex;
317
318};
319
320//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
321//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
322
323inline G4double G4VMultipleScattering::ContinuousStepLimit(
324 const G4Track& track,
325 G4double previousStepSize,
326 G4double currentMinimalStep,
327 G4double& currentSafety)
328{
329 return GetContinuousStepLimit(track,previousStepSize,currentMinimalStep,
330 currentSafety);
331}
332
333//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
334
335inline void G4VMultipleScattering::SetBinning(G4int nbins)
336{
337 nBins = nbins;
338}
339
340//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
341
342inline G4int G4VMultipleScattering::Binning() const
343{
344 return nBins;
345}
346
347//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
348
349inline void G4VMultipleScattering::SetMinKinEnergy(G4double e)
350{
351 minKinEnergy = e;
352}
353
354//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
355
356inline G4double G4VMultipleScattering::MinKinEnergy() const
357{
358 return minKinEnergy;
359}
360
361//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
362
363inline void G4VMultipleScattering::SetMaxKinEnergy(G4double e)
364{
365 maxKinEnergy = e;
366}
367
368//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
369
370inline G4double G4VMultipleScattering::MaxKinEnergy() const
371{
372 return maxKinEnergy;
373}
374
375//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
376
377inline void G4VMultipleScattering::SetBuildLambdaTable(G4bool val)
378{
379 buildLambdaTable = val;
380}
381
382//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
383
384inline G4PhysicsTable* G4VMultipleScattering::LambdaTable() const
385{
386 return theLambdaTable;
387}
388
389//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
390
391inline const G4ParticleDefinition* G4VMultipleScattering::Particle() const
392{
393 return currentParticle;
394}
395
396//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
397
398inline G4VEmModel* G4VMultipleScattering::SelectModel(G4double kinEnergy)
399{
400 return modelManager->SelectModel(kinEnergy, currentMaterialIndex);
401}
402
403//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
404
405inline G4VEmModel* G4VMultipleScattering::SelectModelForMaterial(
406 G4double kinEnergy, size_t& idxRegion) const
407{
408 return modelManager->SelectModel(kinEnergy, idxRegion);
409}
410
411//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
412
413inline G4bool G4VMultipleScattering::LateralDisplasmentFlag() const
414{
415 return latDisplasment;
416}
417
418//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
419
420inline void G4VMultipleScattering::SetLateralDisplasmentFlag(G4bool val)
421{
422 latDisplasment = val;
423}
424
425//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
426
427inline G4double G4VMultipleScattering::Skin() const
428{
429 return skin;
430}
431
432//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
433
434inline void G4VMultipleScattering::SetSkin(G4double val)
435{
436 if(val < 1.0) skin = 0.0;
437 else skin = val;
438}
439
440//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
441
442inline G4double G4VMultipleScattering::RangeFactor() const
443{
444 return facrange;
445}
446
447//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
448
449inline void G4VMultipleScattering::SetRangeFactor(G4double val)
450{
451 if(val > 0.0) facrange = val;
452}
453
454//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
455
456inline G4double G4VMultipleScattering::GeomFactor() const
457{
458 return facgeom;
459}
460
461//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
462
463inline void G4VMultipleScattering::SetGeomFactor(G4double val)
464{
465 if(val > 0.0) facgeom = val;
466}
467
468//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
469
470inline G4double G4VMultipleScattering::PolarAngleLimit() const
471{
472 return polarAngleLimit;
473}
474
475//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
476
477inline void G4VMultipleScattering::SetPolarAngleLimit(G4double val)
478{
479 if(val < 0.0) polarAngleLimit = 0.0;
480 else if(val > pi) polarAngleLimit = pi;
481 else polarAngleLimit = val;
482}
483
484//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
485
486inline G4MscStepLimitType G4VMultipleScattering::StepLimitType() const
487{
488 return stepLimit;
489}
490
491//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
492
493inline void G4VMultipleScattering::SetStepLimitType(G4MscStepLimitType val)
494{
495 stepLimit = val;
496 if(val == fMinimal) facrange = 0.2;
497}
498
499//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
500
501inline
502G4double G4VMultipleScattering::GetLambda(const G4ParticleDefinition* p,
503 G4double& e)
504{
505 G4double x;
506 if(theLambdaTable) {
507 x = ((*theLambdaTable)[currentMaterialIndex])->Value(e);
508 } else {
509 x = currentModel->CrossSection(currentCouple,p,e);
510 }
511 if(x > DBL_MIN) x = 1./x;
512 else x = DBL_MAX;
513 return x;
514}
515
516//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
517
518inline
519void G4VMultipleScattering::DefineMaterial(const G4MaterialCutsCouple* couple)
520{
521 if(couple != currentCouple) {
522 currentCouple = couple;
523 currentMaterialIndex = couple->GetIndex();
524 }
525}
526
527//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
528
529inline const G4MaterialCutsCouple*
530G4VMultipleScattering::CurrentMaterialCutsCouple() const
531{
532 return currentCouple;
533}
534
535//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
536//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
537
538// Follwoing methods are virtual, they are inlined because they applied at
539// each simulation step and some compilers may inline these methods
540
541inline G4double
542G4VMultipleScattering::PostStepGetPhysicalInteractionLength(
543 const G4Track&, G4double, G4ForceCondition* condition)
544{
545 *condition = Forced;
546 return DBL_MAX;
547}
548
549//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
550
551inline G4VParticleChange*
552G4VMultipleScattering::AlongStepDoIt(const G4Track& track, const G4Step& step)
553{
554 if(currentModel->IsActive(track.GetKineticEnergy())) {
555 fParticleChange.ProposeTrueStepLength(currentModel->ComputeTrueStepLength(step.GetStepLength()));
556 }
557 return &fParticleChange;
558}
559
560//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
561
562inline G4VParticleChange*
563G4VMultipleScattering::PostStepDoIt(const G4Track& track, const G4Step& step)
564{
565 fParticleChange.Initialize(track);
566 if(currentModel->IsActive(track.GetKineticEnergy())) {
567 currentModel->SampleScattering(track.GetDynamicParticle(),
568 step.GetPostStepPoint()->GetSafety());
569 }
570 return &fParticleChange;
571}
572
573//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
574
575#endif
Note: See TracBrowser for help on using the repository browser.