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

Last change on this file since 991 was 991, checked in by garnier, 17 years ago

update

File size: 18.9 KB
RevLine 
[819]1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
[991]26// $Id: G4VMultipleScattering.hh,v 1.54 2008/07/31 13:01:26 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02 $
[819]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)
[961]65// 15-07-08 Reorder class members for further multi-thread development (VI)
[819]66//
67
68// -------------------------------------------------------------------
69//
70
71#ifndef G4VMultipleScattering_h
72#define G4VMultipleScattering_h 1
73
74#include "G4VContinuousDiscreteProcess.hh"
75#include "globals.hh"
76#include "G4Material.hh"
77#include "G4MaterialCutsCouple.hh"
78#include "G4ParticleChangeForMSC.hh"
79#include "G4Track.hh"
80#include "G4Step.hh"
81#include "G4EmModelManager.hh"
82#include "G4VEmModel.hh"
83#include "G4MscStepLimitType.hh"
84
85class G4ParticleDefinition;
86class G4DataVector;
87class G4PhysicsTable;
88class G4PhysicsVector;
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91
92class G4VMultipleScattering : public G4VContinuousDiscreteProcess
93{
94public:
95
96 G4VMultipleScattering(const G4String& name = "msc",
[961]97 G4ProcessType type = fElectromagnetic);
[819]98
99 virtual ~G4VMultipleScattering();
100
101 //------------------------------------------------------------------------
102 // Virtual methods to be implemented for the concrete model
103 //------------------------------------------------------------------------
104
105 virtual G4bool IsApplicable(const G4ParticleDefinition& p) = 0;
106
107 virtual void PrintInfo() = 0;
108
109protected:
110
111 virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
112
113public:
114
115 //------------------------------------------------------------------------
116 // Generic methods common to all ContinuousDiscrete processes
117 //------------------------------------------------------------------------
118
119 // Initialise for build of tables
120 void PreparePhysicsTable(const G4ParticleDefinition&);
121
122 // Build physics table during initialisation
123 void BuildPhysicsTable(const G4ParticleDefinition&);
124
125 // Print out of generic class parameters
126 void PrintInfoDefinition();
127
128 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&);
129
130 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
131
132 // Store PhysicsTable in a file.
133 // Return false in case of failure at I/O
134 G4bool StorePhysicsTable(const G4ParticleDefinition*,
135 const G4String& directory,
[961]136 G4bool ascii = false);
[819]137
138 // Retrieve Physics from a file.
139 // (return true if the Physics Table can be build by using file)
140 // (return false if the process has no functionality or in case of failure)
141 // File name should is constructed as processName+particleName and the
142 // should be placed under the directory specifed by the argument.
143 G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
144 const G4String& directory,
[961]145 G4bool ascii);
[819]146
[991]147 //------------------------------------------------------------------------
148 // Specific methods for msc processes
149 //------------------------------------------------------------------------
150
[819]151 // The function overloads the corresponding function of the base
152 // class.It limits the step near to boundaries only
153 // and invokes the method GetMscContinuousStepLimit at every step.
[991]154 virtual G4double AlongStepGetPhysicalInteractionLength(
[819]155 const G4Track&,
[961]156 G4double previousStepSize,
157 G4double currentMinimalStep,
158 G4double& currentSafety,
159 G4GPILSelection* selection);
[819]160
161 // The function overloads the corresponding function of the base
162 // class.
163 G4double PostStepGetPhysicalInteractionLength(
164 const G4Track&,
165 G4double previousStepSize,
166 G4ForceCondition* condition);
167
168 // This method does not used for tracking, it is intended only for tests
169 inline G4double ContinuousStepLimit(const G4Track& track,
170 G4double previousStepSize,
171 G4double currentMinimalStep,
172 G4double& currentSafety);
173
174 //------------------------------------------------------------------------
175 // Specific methods to build and access Physics Tables
176 //------------------------------------------------------------------------
177
178 // Build empty Physics Vector
179 G4PhysicsVector* PhysicsVector(const G4MaterialCutsCouple*);
180
181 inline void SetBinning(G4int nbins);
182 inline G4int Binning() const;
183
184 inline void SetMinKinEnergy(G4double e);
185 inline G4double MinKinEnergy() const;
186
187 inline void SetMaxKinEnergy(G4double e);
188 inline G4double MaxKinEnergy() const;
189
190 inline void SetBuildLambdaTable(G4bool val);
191
192 inline G4PhysicsTable* LambdaTable() const;
193
[991]194 //------------------------------------------------------------------------
195 // Define and access particle type
196 //------------------------------------------------------------------------
197
[819]198 inline const G4ParticleDefinition* Particle() const;
[991]199 inline void SetParticle(const G4ParticleDefinition*);
[819]200
201 //------------------------------------------------------------------------
202 // Specific methods to set, access, modify models
203 //------------------------------------------------------------------------
204
[991]205 inline void AddEmModel(G4int, G4VEmModel*, const G4Region* region = 0);
[819]206
207 inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
208 size_t& idxRegion) const;
209
[991]210 // Access to models
[961]211 inline G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false);
212
[819]213 //------------------------------------------------------------------------
[991]214 // Set parameters for simulation of multiple scattering
[819]215 //------------------------------------------------------------------------
216
217 inline void SetLateralDisplasmentFlag(G4bool val);
218
219 inline void SetSkin(G4double val);
220
221 inline void SetRangeFactor(G4double val);
222
223 inline void SetGeomFactor(G4double val);
224
[961]225 inline void SetPolarAngleLimit(G4double val);
226
[819]227 inline void SetStepLimitType(G4MscStepLimitType val);
228
229protected:
230
231 // This method is used for tracking, it returns mean free path value
232 G4double GetMeanFreePath(const G4Track& track,
233 G4double,
234 G4ForceCondition* condition);
235
[991]236 //------------------------------------------------------------------------
237 // Run time methods
238 //------------------------------------------------------------------------
239
[819]240 // This method is not used for tracking, it returns step limit
241 G4double GetContinuousStepLimit(const G4Track& track,
242 G4double previousStepSize,
243 G4double currentMinimalStep,
244 G4double& currentSafety);
245
[961]246 inline G4double GetLambda(const G4ParticleDefinition* p,
247 G4double& kineticEnergy);
[819]248
249 // This method is used for tracking, it returns step limit
250 inline G4double GetMscContinuousStepLimit(const G4Track& track,
[961]251 G4double scaledKinEnergy,
[819]252 G4double currentMinimalStep,
253 G4double& currentSafety);
254
[991]255 inline G4VEmModel* SelectModel(G4double kinEnergy);
256 // Select concrete model
257
258 inline const G4MaterialCutsCouple* CurrentMaterialCutsCouple() const;
259
[961]260 // define current material
261 inline void DefineMaterial(const G4MaterialCutsCouple* couple);
[819]262
[991]263 //------------------------------------------------------------------------
264 // Access parameters of multiple scattering
265 //------------------------------------------------------------------------
[819]266
[961]267 inline G4ParticleChangeForMSC* GetParticleChange();
[819]268
[991]269 inline G4double Skin() const;
[819]270
[991]271 inline G4double RangeFactor() const;
272
273 inline G4double GeomFactor() const;
274
275 inline G4double PolarAngleLimit() const;
276
277 inline G4MscStepLimitType StepLimitType() const;
278
279 inline G4bool LateralDisplasmentFlag() const;
280
[961]281private:
[819]282
[961]283 // hide assignment operator
[991]284
[961]285 G4VMultipleScattering(G4VMultipleScattering &);
286 G4VMultipleScattering & operator=(const G4VMultipleScattering &right);
[819]287
[961]288 // ======== Parameters of the class fixed at construction =========
[819]289
[961]290 G4EmModelManager* modelManager;
291 G4bool buildLambdaTable;
[819]292
[961]293 // ======== Parameters of the class fixed at initialisation =======
[819]294
[961]295 G4PhysicsTable* theLambdaTable;
296 const G4ParticleDefinition* firstParticle;
[819]297
[961]298 G4MscStepLimitType stepLimit;
[819]299
[961]300 G4double minKinEnergy;
301 G4double maxKinEnergy;
302 G4double skin;
303 G4double facrange;
304 G4double facgeom;
305 G4double polarAngleLimit;
[819]306
[961]307 G4int nBins;
[819]308
[961]309 G4bool latDisplasment;
310
311 // ======== Cashed values - may be state dependent ================
312
[819]313protected:
314
315 G4GPILSelection valueGPILSelectionMSC;
316 G4ParticleChangeForMSC fParticleChange;
317
318private:
319
320 G4VEmModel* currentModel;
321
322 // cache
323 const G4ParticleDefinition* currentParticle;
324 const G4MaterialCutsCouple* currentCouple;
325 size_t currentMaterialIndex;
326
[961]327};
[819]328
[961]329//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
330//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[819]331
[991]332inline
333void G4VMultipleScattering::DefineMaterial(const G4MaterialCutsCouple* couple)
[961]334{
[991]335 if(couple != currentCouple) {
336 currentCouple = couple;
337 currentMaterialIndex = couple->GetIndex();
338 }
[961]339}
[819]340
[961]341//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
[819]342
[991]343inline G4double G4VMultipleScattering::GetMscContinuousStepLimit(
344 const G4Track& track,
345 G4double scaledKinEnergy,
346 G4double currentMinimalStep,
347 G4double&)
[961]348{
[991]349 G4double x = currentMinimalStep;
350 DefineMaterial(track.GetMaterialCutsCouple());
351 currentModel = SelectModel(scaledKinEnergy);
352 if(x > 0.0 && scaledKinEnergy > 0.0) {
353 G4double tPathLength =
354 currentModel->ComputeTruePathLengthLimit(track, theLambdaTable, x);
355 if (tPathLength < x) valueGPILSelectionMSC = CandidateForSelection;
356 x = currentModel->ComputeGeomPathLength(tPathLength);
357 // G4cout << "tPathLength= " << tPathLength
358 // << " stepLimit= " << x
359 // << " currentMinimalStep= " << currentMinimalStep<< G4endl;
360 }
361 return x;
[961]362}
363
[819]364//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
365
[991]366inline G4double G4VMultipleScattering::ContinuousStepLimit(
367 const G4Track& track,
368 G4double previousStepSize,
369 G4double currentMinimalStep,
370 G4double& currentSafety)
[819]371{
[991]372 return GetMscContinuousStepLimit(track,previousStepSize,currentMinimalStep,
373 currentSafety);
[819]374}
375
376//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
377
[991]378inline
379G4double G4VMultipleScattering::GetLambda(const G4ParticleDefinition* p,
380 G4double& e)
[819]381{
[991]382 G4double x;
383 if(theLambdaTable) {
384 G4bool b;
385 x = ((*theLambdaTable)[currentMaterialIndex])->GetValue(e, b);
386 } else {
387 x = currentModel->CrossSection(currentCouple,p,e);
388 }
389 if(x > DBL_MIN) x = 1./x;
390 else x = DBL_MAX;
391 return x;
[819]392}
393
394//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
395
[991]396inline G4VEmModel* G4VMultipleScattering::SelectModel(G4double kinEnergy)
[819]397{
[991]398 return modelManager->SelectModel(kinEnergy, currentMaterialIndex);
[819]399}
400
401//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
402
[991]403inline G4VEmModel* G4VMultipleScattering::SelectModelForMaterial(
404 G4double kinEnergy, size_t& idxRegion) const
[819]405{
[991]406 return modelManager->SelectModel(kinEnergy, idxRegion);
[819]407}
408
409//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
410
[991]411inline void G4VMultipleScattering::SetBinning(G4int nbins)
[819]412{
[991]413 nBins = nbins;
[819]414}
415
416//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
417
[991]418inline G4int G4VMultipleScattering::Binning() const
[819]419{
[991]420 return nBins;
[819]421}
422
423//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
424
[991]425inline void G4VMultipleScattering::SetMinKinEnergy(G4double e)
[819]426{
[991]427 minKinEnergy = e;
[819]428}
429
430//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
431
[991]432inline G4double G4VMultipleScattering::MinKinEnergy() const
[819]433{
[991]434 return minKinEnergy;
[819]435}
436
437//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
438
[991]439inline void G4VMultipleScattering::SetMaxKinEnergy(G4double e)
[819]440{
[991]441 maxKinEnergy = e;
[819]442}
443
444//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
445
[991]446inline G4double G4VMultipleScattering::MaxKinEnergy() const
[819]447{
[991]448 return maxKinEnergy;
[819]449}
450
451//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
452
[991]453inline G4bool G4VMultipleScattering::LateralDisplasmentFlag() const
[819]454{
[991]455 return latDisplasment;
[819]456}
457
458//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
459
[991]460inline void G4VMultipleScattering::SetLateralDisplasmentFlag(G4bool val)
[819]461{
[991]462 latDisplasment = val;
[819]463}
464
465//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
466
[991]467inline G4ParticleChangeForMSC* G4VMultipleScattering::GetParticleChange()
[819]468{
[991]469 return &fParticleChange;
[819]470}
471
472//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
473
474inline G4double G4VMultipleScattering::Skin() const
475{
476 return skin;
477}
478
479//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
480
481inline void G4VMultipleScattering::SetSkin(G4double val)
482{
[961]483 if(val < 1.0) skin = 0.0;
484 else skin = val;
[819]485}
486
487//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
488
489inline G4double G4VMultipleScattering::RangeFactor() const
490{
491 return facrange;
492}
493
494//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
495
496inline void G4VMultipleScattering::SetRangeFactor(G4double val)
497{
498 if(val > 0.0) facrange = val;
499}
500
501//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
502
503inline G4double G4VMultipleScattering::GeomFactor() const
504{
505 return facgeom;
506}
507
508//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
509
510inline void G4VMultipleScattering::SetGeomFactor(G4double val)
511{
512 if(val > 0.0) facgeom = val;
513}
514
515//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
516
[961]517inline G4double G4VMultipleScattering::PolarAngleLimit() const
518{
519 return polarAngleLimit;
520}
521
522//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
523
524inline void G4VMultipleScattering::SetPolarAngleLimit(G4double val)
525{
526 if(val < 0.0) polarAngleLimit = 0.0;
527 else if(val > pi) polarAngleLimit = pi;
528 else polarAngleLimit = val;
529}
530
531//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
532
[819]533inline G4MscStepLimitType G4VMultipleScattering::StepLimitType() const
534{
535 return stepLimit;
536}
537
538//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
539
540inline void G4VMultipleScattering::SetStepLimitType(G4MscStepLimitType val)
541{
542 stepLimit = val;
[961]543 if(val == fMinimal) facrange = 0.2;
[819]544}
545
546//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
547
[991]548inline void G4VMultipleScattering::SetBuildLambdaTable(G4bool val)
[819]549{
[991]550 buildLambdaTable = val;
[819]551}
552
553//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
554
[991]555inline const G4ParticleDefinition* G4VMultipleScattering::Particle() const
[819]556{
[991]557 return currentParticle;
[819]558}
559
560//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
561
[991]562inline G4PhysicsTable* G4VMultipleScattering::LambdaTable() const
[819]563{
[991]564 return theLambdaTable;
[819]565}
566
567//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
568
[991]569inline
570const G4MaterialCutsCouple* G4VMultipleScattering::CurrentMaterialCutsCouple() const
[819]571{
572 return currentCouple;
573}
574
575//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
576
[991]577inline void G4VMultipleScattering::AddEmModel(G4int order, G4VEmModel* p,
578 const G4Region* region)
[819]579{
[991]580 G4VEmFluctuationModel* fm = 0;
581 modelManager->AddEmModel(order, p, fm, region);
582 if(p) p->SetParticleChange(pParticleChange);
[819]583}
584
585//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
586
[991]587inline
588G4VEmModel* G4VMultipleScattering::GetModelByIndex(G4int idx, G4bool ver)
589{
590 return modelManager->GetModel(idx, ver);
591}
592
593//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
594
[819]595#endif
Note: See TracBrowser for help on using the repository browser.