Geant4 11.2.2
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.

Constructor & Destructor Documentation

◆ G4TransportationWithMsc()

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

Definition at line 68 of file G4TransportationWithMsc.cc.

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

◆ ~G4TransportationWithMsc()

G4TransportationWithMsc::~G4TransportationWithMsc ( )
override

Definition at line 93 of file G4TransportationWithMsc.cc.

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

Member Function Documentation

◆ AddMscModel()

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

Definition at line 108 of file G4TransportationWithMsc.cc.

110{
112 G4Exception("G4TransportationWithMsc::AddMscModel", "em0051", FatalException,
113 "not allowed unless type == MultipleScattering");
114 }
115
116 fModelManager->AddEmModel(order, mscModel, nullptr, region);
117 mscModel->SetParticleChange(fParticleChangeForMSC);
118}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fm, const G4Region *r)
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 122 of file G4TransportationWithMsc.cc.

123{
125 G4Exception("G4TransportationWithMsc::AddSSModel", "em0051", FatalException,
126 "not allowed unless type == SingleScattering");
127 }
128
129 fModelManager->AddEmModel(order, model, nullptr, region);
130 model->SetPolarAngleLimit(0.0);
131 model->SetParticleChange(fParticleChangeForSS);
132}
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 276 of file G4TransportationWithMsc.cc.

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

◆ BuildPhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 186 of file G4TransportationWithMsc.cc.

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

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

◆ 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 254 of file G4TransportationWithMsc.cc.

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

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