Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VMultipleScattering.cc
Go to the documentation of this file.
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$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4VMultipleScattering
34//
35// Author: Vladimir Ivanchenko on base of Laszlo Urban code
36//
37// Creation date: 25.03.2003
38//
39// Modifications:
40//
41// 13.04.03 Change printout (V.Ivanchenko)
42// 04-06-03 Fix compilation warnings (V.Ivanchenko)
43// 16-07-03 Use G4VMscModel interface (V.Ivanchenko)
44// 03-11-03 Fix initialisation problem in RetrievePhysicsTable (V.Ivanchenko)
45// 04-11-03 Update PrintInfoDefinition (V.Ivanchenko)
46// 01-03-04 SampleCosineTheta signature changed
47// 22-04-04 SampleCosineTheta signature changed back to original
48// 27-08-04 Add InitialiseForRun method (V.Ivanchneko)
49// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
50// 11-03-05 Shift verbose level by 1 (V.Ivantchenko)
51// 15-04-05 optimize internal interface (V.Ivanchenko)
52// 15-04-05 remove boundary flag (V.Ivanchenko)
53// 27-10-05 introduce virtual function MscStepLimitation() (V.Ivanchenko)
54// 12-04-07 Add verbosity at destruction (V.Ivanchenko)
55// 27-10-07 Virtual functions moved to source (V.Ivanchenko)
56// 11-03-08 Set skin value does not effect step limit type (V.Ivanchenko)
57// 24-06-09 Removed hidden bin in G4PhysicsVector (V.Ivanchenko)
58//
59// Class Description:
60//
61// It is the generic process of multiple scattering it includes common
62// part of calculations for all charged particles
63
64// -------------------------------------------------------------------
65//
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
68
71#include "G4SystemOfUnits.hh"
72#include "G4LossTableManager.hh"
74#include "G4Step.hh"
77#include "G4UnitsTable.hh"
79#include "G4Electron.hh"
80#include "G4GenericIon.hh"
82#include "G4SafetyHelper.hh"
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85
88 numberOfModels(0),
89 firstParticle(0),
90 currParticle(0),
91 stepLimit(fUseSafety),
92 skin(1.0),
93 facrange(0.04),
94 facgeom(2.5),
95 latDisplasment(true),
96 isIon(false)
97{
100 if("ionmsc" == name) { firstParticle = G4GenericIon::GenericIon(); }
101
102 geomMin = 1.e-6*CLHEP::mm;
103 lowestKinEnergy = 1*eV;
104
105 // default limit on polar angle
106 polarAngleLimit = 0.0;
107
108 physStepLimit = gPathLength = tPathLength = 0.0;
109 fIonisation = 0;
110
112 safetyHelper = 0;
113 fPositionChanged = false;
114 isActive = false;
115
116 modelManager = new G4EmModelManager();
117 emManager = G4LossTableManager::Instance();
118 emManager->Register(this);
119
120 warn = 0;
121}
122
123//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
124
126{
127 if(1 < verboseLevel) {
128 G4cout << "G4VMultipleScattering destruct " << GetProcessName()
129 << G4endl;
130 }
131 delete modelManager;
132 emManager->DeRegister(this);
133}
134
135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
136
138 const G4Region* region)
139{
140 G4VEmFluctuationModel* fm = 0;
141 modelManager->AddEmModel(order, p, fm, region);
143}
144
145//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
146
148{
149 ++warn;
150 if(warn < 10) {
151 G4cout << "### G4VMultipleScattering::SetModel is obsolete method "
152 << "and will be removed for the next release." << G4endl;
153 G4cout << " Please, use SetEmModel" << G4endl;
154 }
155 G4int n = mscModels.size();
156 if(index >= n) { for(G4int i=n; i<=index; ++i) {mscModels.push_back(0);} }
157 mscModels[index] = p;
158}
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
161
163{
164 ++warn;
165 if(warn < 10) {
166 G4cout << "### G4VMultipleScattering::Model is obsolete method "
167 << "and will be removed for the next release." << G4endl;
168 G4cout << " Please, use EmModel" << G4endl;
169 }
170 G4VMscModel* p = 0;
171 if(index >= 0 && index < G4int(mscModels.size())) { p = mscModels[index]; }
172 return p;
173}
174
175//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
176
178{
179 G4int n = mscModels.size();
180 if(index >= n) { for(G4int i=n; i<=index; ++i) { mscModels.push_back(0); } }
181 mscModels[index] = p;
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
185
187{
188 G4VMscModel* p = 0;
189 if(index >= 0 && index < G4int(mscModels.size())) { p = mscModels[index]; }
190 return p;
191}
192
193//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
194
197{
198 return modelManager->GetModel(idx, ver);
199}
200
201//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
202
203void
205{
206 if(!firstParticle) { firstParticle = &part; }
207 if(part.GetParticleType() == "nucleus") {
210 SetRangeFactor(0.2);
211 if(&part == G4GenericIon::GenericIon()) { firstParticle = &part; }
212 isIon = true;
213 }
214
215 emManager->PreparePhysicsTable(&part, this);
216 currParticle = 0;
217
218 if(1 < verboseLevel) {
219 G4cout << "### G4VMultipleScattering::PrepearPhysicsTable() for "
220 << GetProcessName()
221 << " and particle " << part.GetParticleName()
222 << " local particle " << firstParticle->GetParticleName()
223 << " isIon= " << isIon
224 << G4endl;
225 }
226
227 if(firstParticle == &part) {
228
229 InitialiseProcess(firstParticle);
230
231 // initialisation of models
232 numberOfModels = modelManager->NumberOfModels();
233 for(G4int i=0; i<numberOfModels; ++i) {
234 G4VMscModel* msc = static_cast<G4VMscModel*>(modelManager->GetModel(i));
235 msc->SetIonisation(0, firstParticle);
236 if(0 == i) { currentModel = msc; }
237 if(isIon) {
239 msc->SetLateralDisplasmentFlag(false);
240 msc->SetRangeFactor(0.2);
241 } else {
244 msc->SetSkin(Skin());
247 }
248 msc->SetPolarAngleLimit(polarAngleLimit);
249 G4double emax =
250 std::min(msc->HighEnergyLimit(),emManager->MaxKinEnergy());
251 msc->SetHighEnergyLimit(emax);
252 }
253
254 modelManager->Initialise(firstParticle, G4Electron::Electron(),
255 10.0, verboseLevel);
256
257 if(!safetyHelper) {
259 ->GetSafetyHelper();
260 safetyHelper->InitialiseHelper();
261 }
262 }
263}
264
265//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
266
268{
269 G4String num = part.GetParticleName();
270 if(1 < verboseLevel) {
271 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
272 << GetProcessName()
273 << " and particle " << num
274 << G4endl;
275 }
276
277 emManager->BuildPhysicsTable(firstParticle);
278
279 // explicitly defined printout by particle name
280 if(1 < verboseLevel ||
281 (0 < verboseLevel && (num == "e-" ||
282 num == "e+" || num == "mu+" ||
283 num == "mu-" || num == "proton"||
284 num == "pi+" || num == "pi-" ||
285 num == "kaon+" || num == "kaon-" ||
286 num == "alpha" || num == "anti_proton" ||
287 num == "GenericIon")))
288 {
290 << ": for " << num
291 << " SubType= " << GetProcessSubType()
292 << G4endl;
293 PrintInfo();
294 modelManager->DumpModelList(verboseLevel);
295 }
296
297 if(1 < verboseLevel) {
298 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() done for "
299 << GetProcessName()
300 << " and particle " << num
301 << G4endl;
302 }
303}
304
305//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
306
308{
309 if (0 < verboseLevel) {
311 << ": for " << firstParticle->GetParticleName()
312 << " SubType= " << GetProcessSubType()
313 << G4endl;
314 PrintInfo();
315 modelManager->DumpModelList(verboseLevel);
316 }
317}
318
319//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
320
322{
323 G4VEnergyLossProcess* eloss = 0;
324 if(track->GetParticleDefinition() != currParticle) {
325 currParticle = track->GetParticleDefinition();
326 fIonisation = emManager->GetEnergyLossProcess(currParticle);
327 eloss = fIonisation;
328 }
329 // one model
330 if(1 == numberOfModels) {
331 currentModel->StartTracking(track);
332 if(eloss) { currentModel->SetIonisation(fIonisation, currParticle); }
333
334 // many models
335 } else {
336 for(G4int i=0; i<numberOfModels; ++i) {
337 G4VMscModel* msc = static_cast<G4VMscModel*>(modelManager->GetModel(i));
338 msc->StartTracking(track);
339 if(eloss) { msc->SetIonisation(fIonisation, currParticle); }
340 }
341 }
342}
343
344//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
345
347 const G4Track& track,
348 G4double,
349 G4double currentMinimalStep,
350 G4double&,
351 G4GPILSelection* selection)
352{
353 // get Step limit proposed by the process
354 *selection = NotCandidateForSelection;
355 physStepLimit = gPathLength = tPathLength = currentMinimalStep;
356
357 G4double ekin = track.GetKineticEnergy();
358 // isIon flag is used only to select a model
359 if(isIon) {
360 ekin *= proton_mass_c2/track.GetParticleDefinition()->GetPDGMass();
361 }
362
363 // select new model
364 if(1 < numberOfModels) {
365 currentModel = static_cast<G4VMscModel*>(
367 }
368
369 // step limit
370 if(currentModel->IsActive(ekin) && gPathLength >= geomMin
371 && ekin >= lowestKinEnergy) {
372 isActive = true;
373 tPathLength = currentModel->ComputeTruePathLengthLimit(track, gPathLength);
374 if (tPathLength < physStepLimit) {
375 *selection = CandidateForSelection;
376 }
377 } else { isActive = false; }
378 /*
379 if(currParticle->GetPDGMass() > GeV)
380 G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
381 << " gPathLength= " << gPathLength
382 << " tPathLength= " << tPathLength
383 << " currentMinimalStep= " << currentMinimalStep
384 << " isActive " << isActive << G4endl;
385 */
386 return gPathLength;
387}
388
389//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
390
394{
395 *condition = Forced;
396 //*condition = NotForced;
397 return DBL_MAX;
398}
399
400//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
401
404{
406 fNewPosition = step.GetPostStepPoint()->GetPosition();
407 fParticleChange.ProposePosition(fNewPosition);
408 fPositionChanged = false;
409
410 G4double geomLength = step.GetStepLength();
411
412 // very small step - no msc
413 if(!isActive) {
414 tPathLength = geomLength;
415
416 // sample msc
417 } else {
418 G4double range =
419 currentModel->GetRange(currParticle,track.GetKineticEnergy(),
420 track.GetMaterialCutsCouple());
421
422 G4double trueLength = currentModel->ComputeTrueStepLength(geomLength);
423
424
425 // protection against wrong t->g->t conversion
426 // if(trueLength > tPathLength)
427 /*
428 if(currParticle->GetPDGMass() > GeV)
429 G4cout << "G4VMsc::AlongStepDoIt: GeomLength= "
430 << geomLength
431 << " trueLenght= " << trueLength
432 << " tPathLength= " << tPathLength
433 << " dr= " << range - trueLength
434 << " ekin= " << track.GetKineticEnergy() << G4endl;
435 */
436 if (trueLength <= physStepLimit) {
437 tPathLength = trueLength;
438 } else {
439 tPathLength = physStepLimit - 0.5*geomMin;
440 }
441
442 // do not sample scattering at the last or at a small step
443 if(tPathLength + geomMin < range && tPathLength > geomMin) {
444
445 G4double preSafety = step.GetPreStepPoint()->GetSafety();
446 G4double postSafety= preSafety - geomLength;
447 G4bool safetyRecomputed = false;
448 if( postSafety < geomMin ) {
449 safetyRecomputed = true;
450 postSafety = currentModel->ComputeSafety(fNewPosition,0.0);
451 }
452 G4ThreeVector displacement =
453 currentModel->SampleScattering(track.GetDynamicParticle(),postSafety);
454
455 G4double r2 = displacement.mag2();
456
457 //G4cout << "R= " << sqrt(r2) << " postSafety= " << postSafety << G4endl;
458
459 // make correction for displacement
460 if(r2 > 0.0) {
461
462 fPositionChanged = true;
463 G4double fac = 1.0;
464
465 // displaced point is definitely within the volume
466 if(r2 > postSafety*postSafety) {
467 if(!safetyRecomputed) {
468 postSafety = currentModel->ComputeSafety(fNewPosition, 0.0);
469 }
470 // add a factor which ensure numerical stability
471 if(r2 > postSafety*postSafety) { fac = 0.99*postSafety/std::sqrt(r2); }
472 }
473 // compute new endpoint of the Step
474 fNewPosition += fac*displacement;
475 //safetyHelper->ReLocateWithinVolume(fNewPosition);
476 }
477 }
478 }
480 //fParticleChange.ProposePosition(fNewPosition);
481 return &fParticleChange;
482}
483
484//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
485
488{
490 if(fPositionChanged) {
491 safetyHelper->ReLocateWithinVolume(fNewPosition);
492 fParticleChange.ProposePosition(fNewPosition);
493 }
494 return &fParticleChange;
495}
496
497//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
498
500 const G4Track& track,
501 G4double previousStepSize,
502 G4double currentMinimalStep,
503 G4double& currentSafety)
504{
506 G4double x = AlongStepGetPhysicalInteractionLength(track,previousStepSize,
507 currentMinimalStep,
508 currentSafety, &selection);
509 return x;
510}
511
512//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
513
515 const G4Track& track,
516 G4double previousStepSize,
517 G4double currentMinimalStep,
518 G4double& currentSafety)
519{
520 return GetContinuousStepLimit(track,previousStepSize,currentMinimalStep,
521 currentSafety);
522}
523
524//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
525
528{
529 *condition = Forced;
530 return DBL_MAX;
531}
532
533//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
534
535G4bool
537 const G4String& directory,
538 G4bool ascii)
539{
540 G4bool yes = true;
541 if(part != firstParticle) { return yes; }
542 G4int nmod = modelManager->NumberOfModels();
543 const G4String ss[4] = {"1","2","3","4"};
544 for(G4int i=0; i<nmod; ++i) {
545 G4VEmModel* msc = modelManager->GetModel(i);
546 yes = true;
547 G4PhysicsTable* table = msc->GetCrossSectionTable();
548 if (table) {
549 G4int j = std::min(i,3);
550 G4String name =
551 GetPhysicsTableFileName(part,directory,"LambdaMod"+ss[j],ascii);
552 yes = table->StorePhysicsTable(name,ascii);
553
554 if ( yes ) {
555 if ( verboseLevel>0 ) {
556 G4cout << "Physics table are stored for " << part->GetParticleName()
557 << " and process " << GetProcessName()
558 << " with a name <" << name << "> " << G4endl;
559 }
560 } else {
561 G4cout << "Fail to store Physics Table for " << part->GetParticleName()
562 << " and process " << GetProcessName()
563 << " in the directory <" << directory
564 << "> " << G4endl;
565 }
566 }
567 }
568 return yes;
569}
570
571//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
572
573G4bool
575 const G4String&,
576 G4bool)
577{
578 return true;
579}
580
581//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
582
584{
585 for(G4int i=0; i<numberOfModels; ++i) {
586 G4VMscModel* msc = static_cast<G4VMscModel*>(modelManager->GetModel(i));
587 msc->SetIonisation(p, firstParticle);
588 }
589}
590
591//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
592
@ fMultipleScattering
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ Forced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
@ fMinimal
@ fUseSafety
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
double mag2() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *, const G4Region *)
G4int NumberOfModels() const
void DumpModelList(G4int verb)
G4VEmModel * GetModel(G4int, G4bool ver=false)
const G4DataVector * Initialise(const G4ParticleDefinition *, const G4ParticleDefinition *, G4double, G4int)
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
static G4LossTableManager * Instance()
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
G4double MaxKinEnergy() const
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
void DeRegister(G4VEnergyLossProcess *p)
void Register(G4VEnergyLossProcess *p)
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
virtual void Initialize(const G4Track &)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void ProposePosition(const G4ThreeVector &finalPosition)
const G4String & GetParticleType() const
const G4String & GetParticleName() const
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
void ReLocateWithinVolume(const G4ThreeVector &pGlobalPoint)
void InitialiseHelper()
G4double GetSafety() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
Definition: G4Step.hh:78
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4ParticleDefinition * GetParticleDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
G4bool IsActive(G4double kinEnergy)
Definition: G4VEmModel.hh:613
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:620
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:318
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:660
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:211
virtual G4ThreeVector & SampleScattering(const G4DynamicParticle *, G4double safety)
Definition: G4VMscModel.cc:131
void SetSkin(G4double)
Definition: G4VMscModel.hh:203
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:210
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:324
void SetLateralDisplasmentFlag(G4bool val)
Definition: G4VMscModel.hh:196
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:288
G4double ComputeSafety(const G4ThreeVector &position, G4double limit)
Definition: G4VMscModel.hh:238
void SetGeomFactor(G4double)
Definition: G4VMscModel.hh:217
virtual G4double ComputeTrueStepLength(G4double geomPathLength)
Definition: G4VMscModel.cc:152
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)
Definition: G4VMscModel.cc:138
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:224
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
void SetIonisation(G4VEnergyLossProcess *)
virtual void PrintInfo()=0
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
G4VMscModel * EmModel(G4int index=1)
G4double GetMeanFreePath(const G4Track &track, G4double, G4ForceCondition *condition)
G4ParticleChangeForMSC fParticleChange
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection)
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
G4VMultipleScattering(const G4String &name="msc", G4ProcessType type=fElectromagnetic)
void PreparePhysicsTable(const G4ParticleDefinition &)
void SetRangeFactor(G4double val)
void SetModel(G4VMscModel *, G4int index=1)
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *condition)
void SetStepLimitType(G4MscStepLimitType val)
void SetEmModel(G4VMscModel *, G4int index=1)
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
G4VMscModel * Model(G4int index=1)
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
void BuildPhysicsTable(const G4ParticleDefinition &)
void SetLateralDisplasmentFlag(G4bool val)
G4bool LateralDisplasmentFlag() const
G4VEmModel * SelectModel(G4double kinEnergy, size_t idx)
G4MscStepLimitType StepLimitType() const
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
void ProposeTrueStepLength(G4double truePathLength)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:408
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4int GetProcessSubType() const
Definition: G4VProcess.hh:397
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:214
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
#define DBL_MAX
Definition: templates.hh:83