Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TransportationWithMsc Class Reference

#include <G4TransportationWithMsc.hh>

+ Inheritance diagram for G4TransportationWithMsc:

Public Types

enum class  ScatteringType { MultipleScattering , SingleScattering }
 

Public Member Functions

 G4TransportationWithMsc (ScatteringType type, G4int verbosity=0)
 
 ~G4TransportationWithMsc () override
 
void SetMultipleSteps (G4bool val)
 
G4bool MultipleSteps () const
 
void AddMscModel (G4VMscModel *mscModel, G4int order=0, const G4Region *region=nullptr)
 
void AddSSModel (G4VEmModel *model, G4int order=0, const G4Region *region=nullptr)
 
void PreparePhysicsTable (const G4ParticleDefinition &part) override
 
void BuildPhysicsTable (const G4ParticleDefinition &part) override
 
void StartTracking (G4Track *track) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection) override
 
- Public Member Functions inherited from G4Transportation
 G4Transportation (G4int verbosityLevel=1, const G4String &aName="Transportation")
 
 ~G4Transportation ()
 
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
 
G4bool FieldExertedForce ()
 
G4PropagatorInFieldGetPropagatorInField ()
 
void SetPropagatorInField (G4PropagatorInField *pFieldPropagator)
 
G4double GetThresholdWarningEnergy () const
 
G4double GetThresholdImportantEnergy () const
 
G4int GetThresholdTrials () const
 
void SetThresholdWarningEnergy (G4double newEnWarn)
 
void SetThresholdImportantEnergy (G4double newEnImp)
 
void SetThresholdTrials (G4int newMaxTrials)
 
void SetHighLooperThresholds ()
 
void SetLowLooperThresholds ()
 
void PushThresholdsToLogger ()
 
void ReportLooperThresholds ()
 
G4double GetMaxEnergyKilled () const
 
G4double GetSumEnergyKilled () const
 
void ResetKilledStatistics (G4int report=1)
 
void EnableShortStepOptimisation (G4bool optimise=true)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual void ProcessDescription (std::ostream &outFile) const
 
void PrintStatistics (std::ostream &outStr) const
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4Transportation
static G4bool EnableMagneticMoment (G4bool useMoment=true)
 
static G4bool EnableGravity (G4bool useGravity=true)
 
static void SetSilenceLooperWarnings (G4bool val)
 
static G4bool GetSilenceLooperWarnings ()
 
static G4bool EnableUseMagneticMoment (G4bool useMoment=true)
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4Transportation
void SetTouchableInformation (const G4TouchableHandle &touchable)
 
void ReportMissingLogger (const char *methodName)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4Transportation
G4NavigatorfLinearNavigator
 
G4PropagatorInFieldfFieldPropagator
 
G4ThreeVector fTransportEndPosition = G4ThreeVector( 0.0, 0.0, 0.0 )
 
G4ThreeVector fTransportEndMomentumDir = G4ThreeVector( 0.0, 0.0, 0.0 )
 
G4double fTransportEndKineticEnergy = 0.0
 
G4ThreeVector fTransportEndSpin = G4ThreeVector( 0.0, 0.0, 0.0 )
 
G4bool fMomentumChanged = true
 
G4bool fEndGlobalTimeComputed = false
 
G4double fCandidateEndGlobalTime = 0.0
 
G4bool fAnyFieldExists = false
 
G4bool fParticleIsLooping = false
 
G4bool fNewTrack = true
 
G4bool fFirstStepInVolume = true
 
G4bool fLastStepInVolume = false
 
G4bool fGeometryLimitedStep = true
 
G4bool fFieldExertedForce = false
 
G4TouchableHandle fCurrentTouchableHandle
 
G4ThreeVector fPreviousSftOrigin
 
G4double fPreviousSafety
 
G4ParticleChangeForTransport fParticleChange
 
G4double fEndPointDistance
 
G4double fThreshold_Warning_Energy = 1.0 * CLHEP::keV
 
G4double fThreshold_Important_Energy = 1.0 * CLHEP::MeV
 
G4int fThresholdTrials = 10
 
G4int fAbandonUnstableTrials = 0
 
G4int fNoLooperTrials = 0
 
G4double fSumEnergyKilled = 0.0
 
G4double fSumEnerSqKilled = 0.0
 
G4double fMaxEnergyKilled = -1.0
 
G4int fMaxEnergyKilledPDG = 0
 
unsigned long fNumLoopersKilled = 0
 
G4double fSumEnergyKilled_NonElectron = 0.0
 
G4double fSumEnerSqKilled_NonElectron = 0.0
 
G4double fMaxEnergyKilled_NonElectron = -1.0
 
G4int fMaxEnergyKilled_NonElecPDG = 0
 
unsigned long fNumLoopersKilled_NonElectron = 0
 
G4double fSumEnergySaved = 0.0
 
G4double fMaxEnergySaved = -1.0
 
G4double fSumEnergyUnstableSaved = 0.0
 
G4bool fShortStepOptimisation
 
G4SafetyHelperfpSafetyHelper
 
G4TransportationLoggerfpLogger
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 
- Static Protected Attributes inherited from G4Transportation
static G4bool fUseMagneticMoment =false
 
static G4bool fUseGravity = false
 
static G4bool fSilenceLooperWarnings = false
 

Detailed Description

Definition at line 57 of file G4TransportationWithMsc.hh.

Member Enumeration Documentation

◆ ScatteringType

Enumerator
MultipleScattering 
SingleScattering 

Definition at line 60 of file G4TransportationWithMsc.hh.

61 {
62 MultipleScattering,
63 SingleScattering,
64 };

Constructor & Destructor Documentation

◆ G4TransportationWithMsc()

G4TransportationWithMsc::G4TransportationWithMsc ( ScatteringType type,
G4int verbosity = 0 )
explicit

Definition at line 69 of file G4TransportationWithMsc.cc.

70 : G4Transportation(verbosity, "TransportationWithMsc"), fType(type)
71{
74
75 fEmManager = G4LossTableManager::Instance();
76 fModelManager = new G4EmModelManager;
77
79 fParticleChangeForMSC = new G4ParticleChangeForMSC;
80 }
81 else if (type == ScatteringType::SingleScattering) {
82 fParticleChangeForSS = new G4ParticleChangeForGamma;
83 fSecondariesSS = new std::vector<G4DynamicParticle*>;
84 }
85
86 G4ThreeVector zero;
87 fSubStepDynamicParticle = new G4DynamicParticle(G4Electron::Definition(), zero);
88 fSubStepTrack = new G4Track(fSubStepDynamicParticle, 0, zero);
89 fSubStep = new G4Step;
90 fSubStepTrack->SetStep(fSubStep);
91}
CLHEP::Hep3Vector G4ThreeVector
int G4int
Definition G4Types.hh:85
static G4Electron * Definition()
Definition G4Electron.cc:45
static G4LossTableManager * Instance()
G4Transportation(G4int verbosityLevel=1, const G4String &aName="Transportation")
void SetVerboseLevel(G4int value)
void SetProcessSubType(G4int)

Referenced by BuildPhysicsTable().

◆ ~G4TransportationWithMsc()

G4TransportationWithMsc::~G4TransportationWithMsc ( )
override

Definition at line 95 of file G4TransportationWithMsc.cc.

96{
97 delete fModelManager;
98 delete fParticleChangeForMSC;
99 delete fEmData;
100 delete fParticleChangeForSS;
101 delete fSecondariesSS;
102
103 // fSubStepDynamicParticle is owned and also deleted by fSubStepTrack!
104 delete fSubStepTrack;
105 delete fSubStep;
106}

Member Function Documentation

◆ AddMscModel()

void G4TransportationWithMsc::AddMscModel ( G4VMscModel * mscModel,
G4int order = 0,
const G4Region * region = nullptr )

Definition at line 110 of file G4TransportationWithMsc.cc.

112{
114 G4Exception("G4TransportationWithMsc::AddMscModel", "em0051", FatalException,
115 "not allowed unless type == MultipleScattering");
116 }
117
118 fModelManager->AddEmModel(order, mscModel, nullptr, region);
119 mscModel->SetParticleChange(fParticleChangeForMSC);
120}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)

Referenced by G4EmBuilder::ConstructElectronMscProcess(), and G4EmConfigurator::PrepareModels().

◆ AddSSModel()

void G4TransportationWithMsc::AddSSModel ( G4VEmModel * model,
G4int order = 0,
const G4Region * region = nullptr )

Definition at line 124 of file G4TransportationWithMsc.cc.

125{
127 G4Exception("G4TransportationWithMsc::AddSSModel", "em0051", FatalException,
128 "not allowed unless type == SingleScattering");
129 }
130
131 fModelManager->AddEmModel(order, model, nullptr, region);
132 model->SetPolarAngleLimit(0.0);
133 model->SetParticleChange(fParticleChangeForSS);
134}
void SetPolarAngleLimit(G4double)

Referenced by G4EmBuilder::ConstructElectronSSProcess().

◆ AlongStepGetPhysicalInteractionLength()

G4double G4TransportationWithMsc::AlongStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & proposedSafety,
G4GPILSelection * selection )
overridevirtual

Reimplemented from G4Transportation.

Definition at line 277 of file G4TransportationWithMsc.cc.

282{
283 *selection = NotCandidateForSelection;
284
285 const G4double physStepLimit = currentMinimumStep;
286
287 switch (fType) {
289 // Select the MSC model for the current kinetic energy.
290 G4VMscModel* mscModel = nullptr;
291 const G4double ekin = track.GetKineticEnergy();
292 const auto* couple = track.GetMaterialCutsCouple();
293 const auto* particleDefinition = track.GetParticleDefinition();
294 if (physStepLimit > kGeomMin) {
295 G4double ekinForSelection = ekin;
296 G4double pdgMass = particleDefinition->GetPDGMass();
297 if (pdgMass > CLHEP::GeV) {
298 ekinForSelection *= proton_mass_c2 / pdgMass;
299 }
300
301 if (ekinForSelection >= kLowestKinEnergy) {
302 mscModel = static_cast<G4VMscModel*>(
303 fModelManager->SelectModel(ekinForSelection, couple->GetIndex()));
304 if (mscModel == nullptr) {
305 G4Exception("G4TransportationWithMsc::AlongStepGPIL", "em0052", FatalException,
306 "no MSC model found");
307 }
308 if (!mscModel->IsActive(ekinForSelection)) {
309 mscModel = nullptr;
310 }
311 }
312 }
313
314 // Call the MSC model to potentially limit the step and convert to
315 // geometric path length.
316 if (mscModel != nullptr) {
317 mscModel->SetCurrentCouple(couple);
318
319 // Use the provided track for the first step.
320 const G4Track* currentTrackPtr = &track;
321
322 G4double currentSafety = proposedSafety;
323 G4double currentEnergy = ekin;
324
325 G4double stepLimitLeft = physStepLimit;
326 G4double totalGeometryStepLength = 0, totalTruePathLength = 0;
327 G4bool firstStep = true, continueStepping = fMultipleSteps;
328
329 do {
330 G4double gPathLength = stepLimitLeft;
331 G4double tPathLength =
332 mscModel->ComputeTruePathLengthLimit(*currentTrackPtr, gPathLength);
333 G4bool mscLimitsStep = (tPathLength < stepLimitLeft);
334 if (!fMultipleSteps && mscLimitsStep) {
335 // MSC limits the step.
336 *selection = CandidateForSelection;
337 }
338
339 if (!firstStep) {
340 // Move the navigator to where the previous step ended.
341 fLinearNavigator->LocateGlobalPointWithinVolume(fTransportEndPosition);
342 }
343
344 G4GPILSelection transportSelection;
346 *currentTrackPtr, previousStepSize, gPathLength, currentSafety, &transportSelection);
347 if (geometryStepLength < gPathLength) {
348 // Transportation limits the step, ie the track hit a boundary.
349 *selection = CandidateForSelection;
350 continueStepping = false;
351 }
352 if (fTransportEndKineticEnergy != currentEnergy) {
353 // Field propagation changed the energy, it's not possible to
354 // estimate the continuous energy loss and continue stepping.
355 continueStepping = false;
356 }
357
358 if (firstStep) {
359 proposedSafety = currentSafety;
360 }
361 totalGeometryStepLength += geometryStepLength;
362
363 // Sample MSC direction change and displacement.
364 const G4double range = mscModel->GetRange(particleDefinition, currentEnergy, couple);
365
366 tPathLength = mscModel->ComputeTrueStepLength(geometryStepLength);
367
368 // Protect against wrong t->g->t conversion.
369 tPathLength = std::min(tPathLength, stepLimitLeft);
370
371 totalTruePathLength += tPathLength;
372 if (*selection != CandidateForSelection && !mscLimitsStep) {
373 // If neither MSC nor transportation limits the step, we got the
374 // distance we want - make sure we exit the loop.
375 continueStepping = false;
376 }
377 else if (tPathLength >= range) {
378 // The particle will stop, exit the loop.
379 continueStepping = false;
380 }
381 else {
382 stepLimitLeft -= tPathLength;
383 }
384
385 // Do not sample scattering at the last or at a small step.
386 if (tPathLength < range && tPathLength > kGeomMin) {
387 static constexpr G4double minSafety = 1.20 * CLHEP::nm;
388 static constexpr G4double sFact = 0.99;
389
390 // The call to SampleScattering() *may* directly fill in the changed
391 // direction into fParticleChangeForMSC, so we have to:
392 // 1) Make sure the momentum direction is initialized.
393 fParticleChangeForMSC->ProposeMomentumDirection(fTransportEndMomentumDir);
394 // 2) Call SampleScattering(), which *may* change it.
395 const G4ThreeVector displacement =
396 mscModel->SampleScattering(fTransportEndMomentumDir, minSafety);
397 // 3) Get the changed direction.
398 fTransportEndMomentumDir = *fParticleChangeForMSC->GetProposedMomentumDirection();
399
400 const G4double r2 = displacement.mag2();
401 if (r2 > kMinDisplacement2) {
402 G4bool positionChanged = true;
403 G4double dispR = std::sqrt(r2);
404 G4double postSafety =
405 sFact * fpSafetyHelper->ComputeSafety(fTransportEndPosition, dispR);
406
407 // Far away from geometry boundary
408 if (postSafety > 0.0 && dispR <= postSafety) {
409 fTransportEndPosition += displacement;
410
411 // Near the boundary
412 }
413 else {
414 // displaced point is definitely within the volume
415 if (dispR < postSafety) {
416 fTransportEndPosition += displacement;
417
418 // reduced displacement
419 }
420 else if (postSafety > kGeomMin) {
421 fTransportEndPosition += displacement * (postSafety / dispR);
422
423 // very small postSafety
424 }
425 else {
426 positionChanged = false;
427 }
428 }
429 if (positionChanged) {
430 fpSafetyHelper->ReLocateWithinVolume(fTransportEndPosition);
431 }
432 }
433 }
434
435 if (continueStepping) {
436 // Update safety according to the geometry distance.
437 if (currentSafety < fEndPointDistance) {
438 currentSafety = 0;
439 }
440 else {
441 currentSafety -= fEndPointDistance;
442 }
443
444 // Update the kinetic energy according to the continuous loss.
445 currentEnergy = mscModel->GetEnergy(particleDefinition, range - tPathLength, couple);
446
447 // From now on, use the track that we can update below.
448 currentTrackPtr = fSubStepTrack;
449
450 fSubStepDynamicParticle->SetKineticEnergy(currentEnergy);
451 fSubStepDynamicParticle->SetMomentumDirection(fTransportEndMomentumDir);
452 fSubStepTrack->SetPosition(fTransportEndPosition);
453
454 G4StepPoint& subPreStepPoint = *fSubStep->GetPreStepPoint();
455 subPreStepPoint.SetMaterialCutsCouple(couple);
456 subPreStepPoint.SetPosition(fTransportEndPosition);
457 subPreStepPoint.SetSafety(currentSafety);
458 subPreStepPoint.SetStepStatus(fAlongStepDoItProc);
459 }
460 firstStep = false;
461 } while (continueStepping);
462
463 // Note: currentEnergy is only updated if continueStepping is true.
464 // In case field propagation changed the energy, this flag is
465 // immediately set to false and currentEnergy is still equal to the
466 // initial kinetic energy stored in ekin.
467 if (currentEnergy != ekin) {
468 // If field propagation didn't change the energy and we potentially
469 // did multiple steps, reset the energy that G4Transportation will
470 // propose to not subtract the energy loss twice.
472 // Also ask for the range again with the initial energy so it is
473 // correctly cached in the G4VEnergyLossProcess.
474 // FIXME: Asking for a range should never change the cached values!
475 (void)mscModel->GetRange(particleDefinition, ekin, couple);
476 }
477
478 fParticleChange.ProposeTrueStepLength(totalTruePathLength);
479
480 // Inform G4Transportation that the momentum might have changed due
481 // to scattering. We do this unconditionally to avoid the situation
482 // where the last step is done without MSC and G4Transportation reset
483 // the flag, for example when running without field.
484 fMomentumChanged = true;
485
486 return totalGeometryStepLength;
487 }
488 break;
489 }
490
492 // Select the model for the current kinetic energy.
493 const G4double ekin = track.GetKineticEnergy();
494 const auto* couple = track.GetMaterialCutsCouple();
495 const auto* particleDefinition = track.GetParticleDefinition();
496
497 G4double ekinForSelection = ekin;
498 G4double pdgMass = particleDefinition->GetPDGMass();
499 if (pdgMass > CLHEP::GeV) {
500 ekinForSelection *= proton_mass_c2 / pdgMass;
501 }
502
503 G4VEmModel* currentModel = fModelManager->SelectModel(ekinForSelection, couple->GetIndex());
504 if (currentModel == nullptr) {
505 G4Exception("G4TransportationWithMsc::AlongStepGPIL", "em0052", FatalException,
506 "no scattering model found");
507 }
508 if (!currentModel->IsActive(ekinForSelection)) {
509 currentModel = nullptr;
510 }
511
512 if (currentModel != nullptr) {
513 currentModel->SetCurrentCouple(couple);
514 G4int coupleIndex = couple->GetIndex();
515
516 // Compute mean free path.
518 G4double lambda = ((*fLambdaTable)[coupleIndex])->LogVectorValue(ekin, logEkin);
519 if (lambda > 0.0) {
520 // Assume that the mean free path and dE/dx are constant along the
521 // step, which is a valid approximation for most cases.
522 G4double meanFreePath = 1.0 / lambda;
523 G4double dedx = fIonisation->GetDEDX(ekin, couple);
524
525 G4double currentSafety = proposedSafety;
526 G4double currentEnergy = ekin;
527
528 // Use the provided track for the first step.
529 const G4Track* currentTrackPtr = &track;
530
531 G4double stepLimitLeft = physStepLimit;
532 G4double totalStepLength = 0;
533 G4bool firstStep = true, continueStepping = fMultipleSteps;
534
535 do {
536 G4double interactionLength = meanFreePath * -G4Log(G4UniformRand());
537
538 G4bool ssLimitsStep = (interactionLength < stepLimitLeft);
539 G4double gPathLength = stepLimitLeft;
540 if (ssLimitsStep) {
541 if (!fMultipleSteps) {
542 // Scattering limits the step.
543 *selection = CandidateForSelection;
544 }
545 gPathLength = interactionLength;
546 }
547
548 if (!firstStep) {
549 // Move the navigator to where the previous step ended.
550 fLinearNavigator->LocateGlobalPointWithinVolume(fTransportEndPosition);
551 }
552
553 G4GPILSelection transportSelection;
555 *currentTrackPtr, previousStepSize, gPathLength, currentSafety, &transportSelection);
556 if (geometryStepLength < gPathLength) {
557 // Transportation limits the step, ie the track hit a boundary.
558 *selection = CandidateForSelection;
559 ssLimitsStep = false;
560 continueStepping = false;
561 }
562 if (fTransportEndKineticEnergy != currentEnergy) {
563 // Field propagation changed the energy, it's not possible to
564 // estimate the continuous energy loss and continue stepping.
565 continueStepping = false;
566 }
567
568 if (firstStep) {
569 proposedSafety = currentSafety;
570 }
571 totalStepLength += geometryStepLength;
572
573 if (*selection != CandidateForSelection && !ssLimitsStep) {
574 // If neither scattering nor transportation limits the step, we
575 // got the distance we want - make sure we exit the loop.
576 continueStepping = false;
577 }
578 else {
579 stepLimitLeft -= geometryStepLength;
580 }
581
582 // Update the kinetic energy according to the continuous loss.
583 G4double energyAfterLinearLoss =
584 fTransportEndKineticEnergy - geometryStepLength * dedx;
585
586 if (ssLimitsStep) {
587 fSubStepDynamicParticle->SetKineticEnergy(energyAfterLinearLoss);
588
589 // The call to SampleSecondaries() directly fills in the changed
590 // direction into fParticleChangeForSS, so we have to:
591 // 1) Set the momentum direction in dynamic particle.
592 fSubStepDynamicParticle->SetMomentumDirection(fTransportEndMomentumDir);
593 // 2) Call SampleSecondaries(), which changes the direction.
594 currentModel->SampleSecondaries(fSecondariesSS, couple, fSubStepDynamicParticle,
595 (*fCuts)[coupleIndex]);
596 // 3) Get the changed direction.
597 fTransportEndMomentumDir = fParticleChangeForSS->GetProposedMomentumDirection();
598
599 // Check that the model neither created secondaries nor proposed
600 // a local energy deposit because this process does not know how
601 // to handle these cases.
602 if (fSecondariesSS->size() > 0) {
603 G4Exception("G4TransportationWithMsc::AlongStepGPIL", "em0053", FatalException,
604 "scattering model created secondaries");
605 }
606 if (fParticleChangeForSS->GetLocalEnergyDeposit() > 0) {
607 G4Exception("G4TransportationWithMsc::AlongStepGPIL", "em0053", FatalException,
608 "scattering model proposed energy deposit");
609 }
610 }
611
612 if (continueStepping) {
613 // Update safety according to the geometry distance.
614 if (currentSafety < fEndPointDistance) {
615 currentSafety = 0;
616 }
617 else {
618 currentSafety -= fEndPointDistance;
619 }
620
621 // Update the energy taking continuous loss into account.
622 currentEnergy = energyAfterLinearLoss;
623
624 // From now on, use the track that we can update below.
625 currentTrackPtr = fSubStepTrack;
626
627 fSubStepDynamicParticle->SetKineticEnergy(currentEnergy);
628 fSubStepDynamicParticle->SetMomentumDirection(fTransportEndMomentumDir);
629 fSubStepTrack->SetPosition(fTransportEndPosition);
630
631 G4StepPoint& subPreStepPoint = *fSubStep->GetPreStepPoint();
632 subPreStepPoint.SetMaterialCutsCouple(couple);
633 subPreStepPoint.SetPosition(fTransportEndPosition);
634 subPreStepPoint.SetSafety(currentSafety);
635 subPreStepPoint.SetStepStatus(fAlongStepDoItProc);
636 }
637 firstStep = false;
638 } while (continueStepping);
639
640 // Note: currentEnergy is only updated if continueStepping is true.
641 // In case field propagation changed the energy, this flag is
642 // immediately set to false and currentEnergy is still equal to the
643 // initial kinetic energy stored in ekin.
644 if (currentEnergy != ekin) {
645 // If field propagation didn't change the energy and we potentially
646 // did multiple steps, reset the energy that G4Transportation will
647 // propose to not subtract the energy loss twice.
649 }
650
651 fParticleChange.ProposeTrueStepLength(totalStepLength);
652
653 // Inform G4Transportation that the momentum might have changed due
654 // to scattering, even if there is no field.
655 fMomentumChanged = true;
656
657 return totalStepLength;
658 }
659 }
660 break;
661 }
662 }
663
664 // If we get here, no scattering has happened.
666 track, previousStepSize, currentMinimumStep, proposedSafety, selection);
667}
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
G4double G4Log(G4double x)
Definition G4Log.hh:227
@ fAlongStepDoItProc
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
#define G4UniformRand()
Definition Randomize.hh:52
double mag2() const
G4double GetLogKineticEnergy() const
void SetSafety(const G4double aValue)
void SetStepStatus(const G4StepStatus aValue)
void SetMaterialCutsCouple(const G4MaterialCutsCouple *)
void SetPosition(const G4ThreeVector &aValue)
const G4ParticleDefinition * GetParticleDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
void SetKineticEnergy(const G4double aValue)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4double fTransportEndKineticEnergy
G4ThreeVector fTransportEndPosition
G4ParticleChangeForTransport fParticleChange
G4Navigator * fLinearNavigator
G4ThreeVector fTransportEndMomentumDir
G4SafetyHelper * fpSafetyHelper
void SetCurrentCouple(const G4MaterialCutsCouple *)
G4bool IsActive(G4double kinEnergy) const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)=0
virtual G4double ComputeTrueStepLength(G4double geomPathLength)=0
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)=0

◆ BuildPhysicsTable()

void G4TransportationWithMsc::BuildPhysicsTable ( const G4ParticleDefinition & part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 187 of file G4TransportationWithMsc.cc.

188{
189 if (fFirstParticle == &part) {
190 fEmManager->BuildPhysicsTable(fFirstParticle);
191
192 if (fEmManager->IsMaster()) {
194 const auto* theParameters = G4EmParameters::Instance();
195 G4LossTableBuilder* bld = fEmManager->GetTableBuilder();
196 const G4ProductionCutsTable* theCoupleTable =
198 std::size_t numOfCouples = theCoupleTable->GetTableSize();
199
200 G4double emin = theParameters->MinKinEnergy();
201 G4double emax = theParameters->MaxKinEnergy();
202
203 G4double scale = emax / emin;
204 G4int nbin = theParameters->NumberOfBinsPerDecade() * G4lrint(std::log10(scale));
205 scale = nbin / G4Log(scale);
206
207 G4int bin = G4lrint(scale * G4Log(emax / emin));
208 bin = std::max(bin, 5);
209
210 for (std::size_t i = 0; i < numOfCouples; ++i) {
211 if (!bld->GetFlag(i)) continue;
212
213 // Create physics vector and fill it
214 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple((G4int)i);
215
216 auto* aVector = new G4PhysicsLogVector(emin, emax, bin, /*splineFlag*/ true);
217 fModelManager->FillLambdaVector(aVector, couple, /*startNull*/ false);
218 aVector->FillSecondDerivatives();
219 G4PhysicsTableHelper::SetPhysicsVector(fLambdaTable, i, aVector);
220 }
221 }
222 }
223 else {
224 const auto masterProcess = static_cast<const G4TransportationWithMsc*>(GetMasterProcess());
225
226 // Initialisation of models.
227 const G4int numberOfModels = fModelManager->NumberOfModels();
229 for (G4int i = 0; i < numberOfModels; ++i) {
230 auto msc = static_cast<G4VMscModel*>(fModelManager->GetModel(i));
231 auto msc0 = static_cast<G4VMscModel*>(masterProcess->fModelManager->GetModel(i));
232 msc->SetCrossSectionTable(msc0->GetCrossSectionTable(), false);
233 msc->InitialiseLocal(fFirstParticle, msc0);
234 }
235 }
236 else if (fType == ScatteringType::SingleScattering) {
237 this->fLambdaTable = masterProcess->fLambdaTable;
238 }
239 }
240 }
241
242 if (!G4EmParameters::Instance()->IsPrintLocked() && verboseLevel > 0) {
243 G4cout << G4endl;
244 G4cout << GetProcessName() << ": for " << part.GetParticleName();
245 if (fMultipleSteps) {
246 G4cout << " (multipleSteps: 1)";
247 }
248 G4cout << G4endl;
249 fModelManager->DumpModelList(G4cout, verboseLevel);
250 }
251}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4EmParameters * Instance()
static G4bool GetFlag(std::size_t idx)
const G4String & GetParticleName() const
static void SetPhysicsVector(G4PhysicsTable *physTable, std::size_t idx, G4PhysicsVector *vec)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4TransportationWithMsc(ScatteringType type, G4int verbosity=0)
const G4VProcess * GetMasterProcess() const
G4int verboseLevel
const G4String & GetProcessName() const
int G4lrint(double ad)
Definition templates.hh:134

◆ MultipleSteps()

G4bool G4TransportationWithMsc::MultipleSteps ( ) const
inline

Definition at line 123 of file G4TransportationWithMsc.hh.

124{
125 return fMultipleSteps;
126}

◆ PreparePhysicsTable()

void G4TransportationWithMsc::PreparePhysicsTable ( const G4ParticleDefinition & part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 138 of file G4TransportationWithMsc.cc.

139{
140 if (nullptr == fFirstParticle) {
141 fFirstParticle = &part;
142 G4VMultipleScattering* ptr = nullptr;
143 auto emConfigurator = fEmManager->EmConfigurator();
144 emConfigurator->PrepareModels(&part, ptr, this);
145 }
146
147 if (fFirstParticle == &part) {
148 G4bool master = fEmManager->IsMaster();
149 G4LossTableBuilder* bld = fEmManager->GetTableBuilder();
150 G4bool baseMat = bld->GetBaseMaterialFlag();
151 const auto* theParameters = G4EmParameters::Instance();
152
153 if (master) {
154 SetVerboseLevel(theParameters->Verbose());
155 }
156 else {
157 SetVerboseLevel(theParameters->WorkerVerbose());
158 }
159
160 const G4int numberOfModels = fModelManager->NumberOfModels();
162 for (G4int i = 0; i < numberOfModels; ++i) {
163 auto msc = static_cast<G4VMscModel*>(fModelManager->GetModel(i));
164 msc->SetPolarAngleLimit(theParameters->MscThetaLimit());
165 G4double emax = std::min(msc->HighEnergyLimit(), theParameters->MaxKinEnergy());
166 msc->SetHighEnergyLimit(emax);
167 msc->SetUseBaseMaterials(baseMat);
168 }
169 }
170 else if (fType == ScatteringType::SingleScattering) {
171 if (master) {
172 if (fEmData == nullptr) {
173 fEmData = new G4EmDataHandler(2);
174 }
175
176 fLambdaTable = fEmData->MakeTable(0);
177 bld->InitialiseBaseMaterials(fLambdaTable);
178 }
179 }
180
181 fCuts = fModelManager->Initialise(fFirstParticle, G4Electron::Electron(), verboseLevel);
182 }
183}
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4bool GetBaseMaterialFlag()
void InitialiseBaseMaterials(const G4PhysicsTable *table=nullptr)

◆ SetMultipleSteps()

void G4TransportationWithMsc::SetMultipleSteps ( G4bool val)
inline

Definition at line 116 of file G4TransportationWithMsc.hh.

117{
118 fMultipleSteps = val;
119}

Referenced by G4EmBuilder::ConstructElectronMscProcess(), and G4EmBuilder::ConstructElectronSSProcess().

◆ StartTracking()

void G4TransportationWithMsc::StartTracking ( G4Track * track)
overridevirtual

Reimplemented from G4Transportation.

Definition at line 255 of file G4TransportationWithMsc.cc.

256{
257 auto* currParticle = track->GetParticleDefinition();
258 fIonisation = fEmManager->GetEnergyLossProcess(currParticle);
259
260 fSubStepDynamicParticle->SetDefinition(currParticle);
261
262 const G4int numberOfModels = fModelManager->NumberOfModels();
264 for (G4int i = 0; i < numberOfModels; ++i) {
265 auto msc = static_cast<G4VMscModel*>(fModelManager->GetModel(i));
266 msc->StartTracking(track);
267 msc->SetIonisation(fIonisation, currParticle);
268 }
269 }
270
271 // Ensure that field propagation state is also cleared / prepared
273}
void StartTracking(G4Track *aTrack)

The documentation for this class was generated from the following files: