Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4VMultipleScattering
33//
34// Author: Vladimir Ivanchenko on base of Laszlo Urban code
35//
36// Creation date: 25.03.2003
37//
38// Modifications:
39//
40// 16-07-03 Use G4VMscModel interface (V.Ivanchenko)
41// 03-11-03 Fix initialisation problem in RetrievePhysicsTable (V.Ivanchenko)
42// 04-11-03 Update PrintInfoDefinition (V.Ivanchenko)
43// 01-03-04 SampleCosineTheta signature changed
44// 22-04-04 SampleCosineTheta signature changed back to original
45// 27-08-04 Add InitialiseForRun method (V.Ivanchneko)
46// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
47// 11-03-05 Shift verbose level by 1 (V.Ivantchenko)
48// 15-04-05 optimize internal interface (V.Ivanchenko)
49// 15-04-05 remove boundary flag (V.Ivanchenko)
50// 27-10-05 introduce virtual function MscStepLimitation() (V.Ivanchenko)
51// 12-04-07 Add verbosity at destruction (V.Ivanchenko)
52// 27-10-07 Virtual functions moved to source (V.Ivanchenko)
53// 11-03-08 Set skin value does not effect step limit type (V.Ivanchenko)
54// 24-06-09 Removed hidden bin in G4PhysicsVector (V.Ivanchenko)
55// 04-06-13 Adoptation to MT mode (V.Ivanchenko)
56//
57
58// -------------------------------------------------------------------
59//
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
65#include "G4SystemOfUnits.hh"
66#include "G4LossTableManager.hh"
68#include "G4Step.hh"
71#include "G4UnitsTable.hh"
73#include "G4Electron.hh"
74#include "G4GenericIon.hh"
76#include "G4SafetyHelper.hh"
77#include "G4ParticleTable.hh"
78#include "G4ProcessVector.hh"
79#include "G4ProcessManager.hh"
80#include "G4LossTableBuilder.hh"
81#include <iostream>
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84
87 fNewPosition(0.,0.,0.),
88 fNewDirection(0.,0.,1.)
89{
93
95
96 geomMin = 0.05*CLHEP::nm;
98
100
103 mscModels.reserve(2);
104 emManager->Register(this);
105}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
108
110{
111 delete modelManager;
112 emManager->DeRegister(this);
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116
118 const G4Region* region)
119{
120 if(nullptr == ptr) { return; }
121 G4VEmFluctuationModel* fm = nullptr;
122 modelManager->AddEmModel(order, ptr, fm, region);
124}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127
129{
130 if(nullptr == ptr) { return; }
131 if(!mscModels.empty()) {
132 for(auto & msc : mscModels) { if(msc == ptr) { return; } }
133 }
134 mscModels.push_back(ptr);
135}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
138
139void
141{
142 if(1 < verboseLevel) {
143 G4cout << "### G4VMultipleScattering::PrepearPhysicsTable() for "
144 << GetProcessName()
145 << " and particle " << part.GetParticleName()
146 << G4endl;
147 }
148 G4bool master = emManager->IsMaster();
149
150 if(nullptr == firstParticle) { firstParticle = &part; }
151 if(part.GetPDGMass() > CLHEP::GeV) {
152 // flag declears that mass scaling is applied
153 isIon = true;
154 }
155
156 emManager->PreparePhysicsTable(&part, this, master);
157 currParticle = nullptr;
158
159 if(1 < verboseLevel) {
160 G4cout << "### G4VMultipleScattering::PrepearPhysicsTable() for "
161 << GetProcessName()
162 << " and particle " << part.GetParticleName()
163 << " local particle " << firstParticle->GetParticleName()
164 << " isIon: " << isIon << " isMaster: " << master
165 << G4endl;
166 }
167
168 if(firstParticle == &part) {
169
170 // initialise process
172
173 // heavy particles
174 if(part.GetPDGMass() > CLHEP::MeV) {
178 } else {
182 }
183 if(master) { SetVerboseLevel(theParameters->Verbose()); }
185
186 // initialisation of models
188 /*
189 std::cout << "### G4VMultipleScattering::PreparePhysicsTable() for "
190 << GetProcessName()
191 << " and particle " << part.GetParticleName()
192 << " Nmodels= " << mscModels.size() << " " << this << std::endl;
193 */
196
197 for(G4int i=0; i<numberOfModels; ++i) {
199 if(nullptr == msc) { continue; }
200 if(nullptr == currentModel) { currentModel = msc; }
201 msc->SetIonisation(nullptr, firstParticle);
202 msc->SetMasterThread(master);
204 G4double emax =
208 }
209
211 1.0, verboseLevel);
212
213 if(nullptr == safetyHelper) {
215 ->GetSafetyHelper();
217 }
218 }
219}
220
221//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
222
224{
225 const G4String& num = part.GetParticleName();
226 G4bool master = emManager->IsMaster();
227 if(1 < verboseLevel) {
228 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
229 << GetProcessName()
230 << " and particle " << num << " isIon: " << isIon
231 << " IsMaster: " << master << G4endl;
232 }
233 const G4VMultipleScattering* masterProcess =
234 static_cast<const G4VMultipleScattering*>(GetMasterProcess());
235
236 if(firstParticle == &part) {
237 /*
238 std::cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
239 << GetProcessName() << " and particle " << num
240 << " IsMaster= " << G4LossTableManager::Instance()->IsMaster()
241 << " " << this << std::endl;
242 */
244
245 if(!master) {
246 // initialisation of models
247 /*
248 std::cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
249 << GetProcessName() << " and particle " << num
250 << " Nmod= " << mscModels.size() << " NOT master" << std::endl;
251 */
252 baseMat = masterProcess->UseBaseMaterial();
253 for(G4int i=0; i<numberOfModels; ++i) {
255 if(nullptr == msc) { continue; }
256 G4VMscModel* msc0 = masterProcess->GetModelByIndex(i);
258 msc->SetCrossSectionTable(msc0->GetCrossSectionTable(), false);
259 msc->InitialiseLocal(firstParticle, msc0);
260 }
261 }
262 }
263 // protection against double printout
264 if(theParameters->IsPrintLocked()) { return; }
265
266 // explicitly defined printout by particle name
267 if(1 < verboseLevel ||
268 (0 < verboseLevel && (num == "e-" ||
269 num == "e+" || num == "mu+" ||
270 num == "mu-" || num == "proton"||
271 num == "pi+" || num == "pi-" ||
272 num == "kaon+" || num == "kaon-" ||
273 num == "alpha" || num == "anti_proton" ||
274 num == "GenericIon" || num == "alpha+" ||
275 num == "alpha++" )))
276 {
277 StreamInfo(G4cout, part);
278 }
279
280 if(1 < verboseLevel) {
281 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() done for "
282 << GetProcessName()
283 << " and particle " << num << G4endl;
284 }
285}
286
287//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
288
289void G4VMultipleScattering::StreamInfo(std::ostream& outFile,
290 const G4ParticleDefinition& part, G4bool rst) const
291{
292 G4String indent = (rst ? " " : "");
293 outFile << G4endl << indent << GetProcessName() << ": ";
294 if (!rst) outFile << " for " << part.GetParticleName();
295 outFile << " SubType= " << GetProcessSubType() << G4endl;
297}
298
299//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
300
302{
303 G4VEnergyLossProcess* eloss = nullptr;
304 if(track->GetParticleDefinition() != currParticle) {
307 eloss = fIonisation;
308 }
309 /*
310 G4cout << "G4VMultipleScattering::StartTracking Nmod= " << numberOfModels
311 << " " << currParticle->GetParticleName()
312 << " E(MeV)= " << track->GetKineticEnergy()
313 << " Ion= " << fIonisation << " IsMaster= "
314 << G4LossTableManager::Instance()->IsMaster()
315 << G4endl;
316 */
317 for(G4int i=0; i<numberOfModels; ++i) {
319 /*
320 G4cout << "Next model " << msc
321 << " Emin= " << msc->LowEnergyLimit()
322 << " Emax= " << msc->HighEnergyLimit()
323 << " Eact= " << msc->LowEnergyActivationLimit() << G4endl;
324 */
325 msc->StartTracking(track);
326 if(nullptr != eloss) {
327 msc->SetIonisation(eloss, currParticle);
328 }
329 }
330}
331
332//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
333
335 const G4Track& track,
336 G4double,
337 G4double currentMinimalStep,
338 G4double&,
339 G4GPILSelection* selection)
340{
341 // get Step limit proposed by the process
342 *selection = NotCandidateForSelection;
343 physStepLimit = gPathLength = tPathLength = currentMinimalStep;
344
345 G4double ekin = track.GetKineticEnergy();
346 /*
347 G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
348 << " " << currParticle->GetParticleName()
349 << " currMod " << currentModel
350 << G4endl;
351 */
352 // isIon flag is used only to select a model
353 if(isIon) {
355 }
356 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
357
358 // select new model
359 if(1 < numberOfModels) {
360 currentModel = static_cast<G4VMscModel*>(SelectModel(ekin,couple->GetIndex()));
361 }
363 // msc is active is model is active, energy above the limit,
364 // and step size is above the limit;
365 // if it is active msc may limit the step
367 && ekin >= lowestKinEnergy) {
368 isActive = true;
369 tPathLength =
371 if (tPathLength < physStepLimit) {
372 *selection = CandidateForSelection;
373 }
374 } else { isActive = false; }
375
376 //if(currParticle->GetPDGMass() > GeV)
377 /*
378 G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
379 << " " << currParticle->GetParticleName()
380 << " gPathLength= " << gPathLength
381 << " tPathLength= " << tPathLength
382 << " currentMinimalStep= " << currentMinimalStep
383 << " isActive " << isActive << G4endl;
384 */
385 return gPathLength;
386}
387
388//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
389
393{
395 return DBL_MAX;
396}
397
398//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
399
402{
407 fPositionChanged = false;
408
409 G4double geomLength = step.GetStepLength();
410
411 // very small step - no msc
412 if(!isActive) {
413 tPathLength = geomLength;
414
415 // sample msc
416 } else {
417 G4double range =
419 track.GetMaterialCutsCouple());
420
422
423 /*
424 if(currParticle->GetPDGMass() > 0.9*GeV)
425 G4cout << "G4VMsc::AlongStepDoIt: GeomLength= "
426 << geomLength
427 << " tPathLength= " << tPathLength
428 << " physStepLimit= " << physStepLimit
429 << " dr= " << range - tPathLength
430 << " ekin= " << track.GetKineticEnergy() << G4endl;
431 */
432 // protection against wrong t->g->t conversion
434
435 // do not sample scattering at the last or at a small step
436 if(tPathLength < range && tPathLength > geomMin) {
437
438 static const G4double minSafety = 1.20*CLHEP::nm;
439 static const G4double sFact = 0.99;
440
442 step.GetPostStepPoint()->GetMomentumDirection(),minSafety);
443
444 G4double r2 = displacement.mag2();
445 //G4cout << " R= " << sqrt(r2) << " Rmin= " << sqrt(minDisplacement2)
446 // << " flag= " << fDispBeyondSafety << G4endl;
447 if(r2 > minDisplacement2) {
448
449 fPositionChanged = true;
450 G4double dispR = std::sqrt(r2);
451 G4double postSafety =
453 //G4cout<<" R= "<< dispR<<" postSafety= "<<postSafety<<G4endl;
454
455 // far away from geometry boundary
456 if(postSafety > 0.0 && dispR <= postSafety) {
457 fNewPosition += displacement;
458
459 //near the boundary
460 } else {
461 // displaced point is definitely within the volume
462 //G4cout<<" R= "<<dispR<<" postSafety= "<<postSafety<<G4endl;
463 if(dispR < postSafety) {
464 fNewPosition += displacement;
465
466 // reduced displacement
467 } else if(postSafety > geomMin) {
468 fNewPosition += displacement*(postSafety/dispR);
469
470 // very small postSafety
471 } else {
472 fPositionChanged = false;
473 }
474 }
475 if(fPositionChanged) {
478 }
479 }
480 }
481 }
483 return &fParticleChange;
484}
485
486//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
487
490{
492 return &fParticleChange;
493}
494
495//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
496
498 const G4Track& track,
499 G4double previousStepSize,
500 G4double currentMinimalStep,
501 G4double& currentSafety)
502{
504 G4double x = AlongStepGetPhysicalInteractionLength(track,previousStepSize,
505 currentMinimalStep,
506 currentSafety,
507 &selection);
508 return x;
509}
510
511//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
512
514 const G4Track& track,
515 G4double previousStepSize,
516 G4double currentMinimalStep,
517 G4double& currentSafety)
518{
519 return GetContinuousStepLimit(track,previousStepSize,currentMinimalStep,
520 currentSafety);
521}
522
523//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
524
527{
528 *condition = Forced;
529 return DBL_MAX;
530}
531
532//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
533
534G4bool
536 const G4String& directory,
537 G4bool ascii)
538{
539 G4bool yes = true;
540 if(part != firstParticle) { return yes; }
541 const G4VMultipleScattering* masterProcess =
542 static_cast<const G4VMultipleScattering*>(GetMasterProcess());
543 if(nullptr != masterProcess && masterProcess != this) { return yes; }
544
546 static const G4String ss[4] = {"1","2","3","4"};
547 for(G4int i=0; i<nmod; ++i) {
549 if(nullptr == msc) { continue; }
550 yes = true;
551 G4PhysicsTable* table = msc->GetCrossSectionTable();
552 if (nullptr != table) {
553 G4int j = std::min(i,3);
554 G4String name =
555 GetPhysicsTableFileName(part,directory,"LambdaMod"+ss[j],ascii);
556 yes = table->StorePhysicsTable(name,ascii);
557
558 if ( yes ) {
559 if ( verboseLevel>0 ) {
560 G4cout << "Physics table are stored for "
561 << part->GetParticleName()
562 << " and process " << GetProcessName()
563 << " with a name <" << name << "> " << G4endl;
564 }
565 } else {
566 G4cout << "Fail to store Physics Table for "
567 << part->GetParticleName()
568 << " and process " << GetProcessName()
569 << " in the directory <" << directory
570 << "> " << G4endl;
571 }
572 }
573 }
574 return yes;
575}
576
577//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
578
579G4bool
581 const G4String&,
582 G4bool)
583{
584 return true;
585}
586
587//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
588
590{
591 for(auto & msc : mscModels) {
592 if(nullptr != msc) { msc->SetIonisation(p, firstParticle); }
593 }
594}
595
596//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
597
598void G4VMultipleScattering::ProcessDescription(std::ostream& outFile) const
599{
600 if(firstParticle) {
601 StreamInfo(outFile, *firstParticle, true);
602 }
603}
604
605//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
606
static const G4double emax
@ fMultipleScattering
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
@ Forced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double mag2() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fm, const G4Region *r)
void DumpModelList(std::ostream &out, G4int verb)
G4int NumberOfModels() const
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4double, G4int verb)
G4VEmModel * GetModel(G4int idx, G4bool ver=false)
G4bool IsPrintLocked() const
static G4EmParameters * Instance()
G4double MscMuHadRangeFactor() const
G4double MscThetaLimit() const
G4MscStepLimitType MscMuHadStepLimitType() const
G4int Verbose() const
G4int WorkerVerbose() const
G4MscStepLimitType MscStepLimitType() const
G4double MaxKinEnergy() const
G4bool LateralDisplacement() const
G4bool MuHadLateralDisplacement() const
G4double MscRangeFactor() const
static G4LossTableManager * Instance()
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
G4LossTableBuilder * GetTableBuilder()
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
void DeRegister(G4VEnergyLossProcess *p)
void Register(G4VEnergyLossProcess *p)
virtual void Initialize(const G4Track &)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void ProposePosition(const G4ThreeVector &finalPosition)
const G4String & GetParticleName() const
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
void ReLocateWithinVolume(const G4ThreeVector &pGlobalPoint)
void InitialiseHelper()
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
Definition: G4Step.hh:62
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:455
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:802
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:739
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:472
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:447
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:655
G4bool IsActive(G4double kinEnergy) const
Definition: G4VEmModel.hh:795
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:870
void SetUseBaseMaterials(G4bool val)
Definition: G4VEmModel.hh:753
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:270
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:204
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)=0
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:309
virtual G4double ComputeTrueStepLength(G4double geomPathLength)=0
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:189
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)=0
void SetIonisation(G4VEnergyLossProcess *)
G4EmModelManager * modelManager
G4double GetMeanFreePath(const G4Track &track, G4double, G4ForceCondition *condition) override
G4ParticleChangeForMSC fParticleChange
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=nullptr)
G4LossTableManager * emManager
G4VEnergyLossProcess * fIonisation
void StartTracking(G4Track *) override
void PreparePhysicsTable(const G4ParticleDefinition &) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
G4VMultipleScattering(const G4String &name="msc", G4ProcessType type=fElectromagnetic)
G4VMscModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
void ProcessDescription(std::ostream &outFile) const override
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
const G4ParticleDefinition * firstParticle
void SetEmModel(G4VMscModel *, G4int idx=0)
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection) override
void StreamInfo(std::ostream &outFile, const G4ParticleDefinition &, G4bool rst=false) const
G4MscStepLimitType stepLimit
std::vector< G4VMscModel * > mscModels
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *condition) override
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4VEmModel * SelectModel(G4double kinEnergy, size_t idx)
const G4ParticleDefinition * currParticle
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety) override
void ProposeTrueStepLength(G4double truePathLength)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:518
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
G4int GetProcessSubType() const
Definition: G4VProcess.hh:400
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:182
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
static constexpr double GeV
static constexpr double MeV
static constexpr double eV
static constexpr double nm
Definition: SystemOfUnits.h:93
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
const char * name(G4int ptype)
float proton_mass_c2
Definition: hepunit.py:274
#define DBL_MAX
Definition: templates.hh:62