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

#include <G4BiasingProcessInterface.hh>

+ Inheritance diagram for G4BiasingProcessInterface:

Public Member Functions

 G4BiasingProcessInterface (G4String name="biasWrapper(0)")
 
 G4BiasingProcessInterface (G4VProcess *wrappedProcess, G4bool wrappedIsAtRest, G4bool wrappedIsAlongStep, G4bool wrappedIsPostStep, G4String useThisName="")
 
 ~G4BiasingProcessInterface ()
 
G4VProcessGetWrappedProcess () const
 
G4VBiasingOperatorGetCurrentBiasingOperator () const
 
G4VBiasingOperatorGetPreviousBiasingOperator () const
 
G4VBiasingOperationGetCurrentNonPhysicsBiasingOperation () const
 
G4VBiasingOperationGetPreviousNonPhysicsBiasingOperation () const
 
G4VBiasingOperationGetCurrentOccurenceBiasingOperation () const
 
G4VBiasingOperationGetPreviousOccurenceBiasingOperation () const
 
G4VBiasingOperationGetCurrentFinalStateBiasingOperation () const
 
G4VBiasingOperationGetPreviousFinalStateBiasingOperation () const
 
const std::vector< const G4BiasingProcessInterface * > & GetBiasingProcessInterfaces () const
 
const std::vector< const G4BiasingProcessInterface * > & GetPhysicsBiasingProcessInterfaces () const
 
const std::vector< const G4BiasingProcessInterface * > & GetNonPhysicsBiasingProcessInterfaces () const
 
const G4BiasingProcessSharedDataGetSharedData () const
 
G4bool GetIsFirstPostStepGPILInterface (G4bool physOnly=true) const
 
G4bool GetIsLastPostStepGPILInterface (G4bool physOnly=true) const
 
G4bool GetIsFirstPostStepDoItInterface (G4bool physOnly=true) const
 
G4bool GetIsLastPostStepDoItInterface (G4bool physOnly=true) const
 
G4bool IsFirstPostStepGPILInterface (G4bool physOnly=true) const
 
G4bool IsLastPostStepGPILInterface (G4bool physOnly=true) const
 
G4bool IsFirstPostStepDoItInterface (G4bool physOnly=true) const
 
G4bool IsLastPostStepDoItInterface (G4bool physOnly=true) const
 
G4bool GetWrappedProcessIsAtRest () const
 
G4bool GetWrappedProcessIsAlong () const
 
G4bool GetWrappedProcessIsPost () const
 
G4double GetPreviousStepSize () const
 
G4double GetCurrentMinimumStep () const
 
G4double GetProposedSafety () const
 
void SetProposedSafety (G4double sft)
 
G4double GetPostStepGPIL () const
 
G4double GetAlongStepGPIL () const
 
void StartTracking (G4Track *track)
 
void EndTracking ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &step)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &step)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &pd)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &pd)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &pd)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *pd, const G4String &s, G4bool f)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *pd, const G4String &s, G4bool f)
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &pd)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &pd)
 
- 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)
 
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
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
const G4VProcessGetMasterProcess () const
 

Static Public Member Functions

static const G4BiasingProcessSharedDataGetSharedData (const G4ProcessManager *)
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Additional Inherited Members

- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- 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
 

Detailed Description

Definition at line 69 of file G4BiasingProcessInterface.hh.

Constructor & Destructor Documentation

◆ G4BiasingProcessInterface() [1/2]

G4BiasingProcessInterface::G4BiasingProcessInterface ( G4String name = "biasWrapper(0)")

Definition at line 42 of file G4BiasingProcessInterface.cc.

43 : G4VProcess ( name ),
44 fCurrentTrack ( nullptr ),
45 fPreviousStepSize (-1.0), fCurrentMinimumStep( -1.0 ), fProposedSafety ( -1.0),
46 fOccurenceBiasingOperation( nullptr ), fFinalStateBiasingOperation( nullptr ), fNonPhysicsBiasingOperation( nullptr ),
47 fPreviousOccurenceBiasingOperation( nullptr ), fPreviousFinalStateBiasingOperation( nullptr ), fPreviousNonPhysicsBiasingOperation( nullptr ),
48 fResetWrappedProcessInteractionLength( true ),
49 fWrappedProcess ( nullptr ),
50 fIsPhysicsBasedBiasing ( false ),
51 fWrappedProcessIsAtRest ( false ),
52 fWrappedProcessIsAlong ( false ),
53 fWrappedProcessIsPost ( false ),
54 fWrappedProcessPostStepGPIL ( -1.0 ),
55 fBiasingPostStepGPIL ( -1.0 ),
56 fWrappedProcessInteractionLength ( -1.0 ),
57 fWrappedProcessForceCondition ( NotForced ),
58 fBiasingForceCondition ( NotForced ),
59 fWrappedProcessAlongStepGPIL ( -1.0 ),
60 fBiasingAlongStepGPIL ( -1.0 ),
61 fWrappedProcessGPILSelection ( NotCandidateForSelection ),
62 fBiasingGPILSelection ( NotCandidateForSelection ),
63 fBiasingInteractionLaw ( nullptr ),
64 fPreviousBiasingInteractionLaw ( nullptr ),
65 fPhysicalInteractionLaw ( nullptr ),
66 fOccurenceBiasingParticleChange ( nullptr ),
67 fDummyParticleChange ( nullptr ),
68 fIamFirstGPIL ( false ),
69 fProcessManager ( nullptr ),
70 fSharedData ( nullptr )
71{
72 for (G4int i = 0 ; i < 8 ; i++) fFirstLastFlags[i] = false;
73 fResetInteractionLaws.Put( true );
74 fCommonStart .Put( true );
75 fCommonEnd .Put( true );
76 fDoCommonConfigure .Put( true );
77}
@ NotForced
@ NotCandidateForSelection
int G4int
Definition G4Types.hh:85
void Put(const value_type &val) const
Definition G4Cache.hh:321

◆ G4BiasingProcessInterface() [2/2]

G4BiasingProcessInterface::G4BiasingProcessInterface ( G4VProcess * wrappedProcess,
G4bool wrappedIsAtRest,
G4bool wrappedIsAlongStep,
G4bool wrappedIsPostStep,
G4String useThisName = "" )

Definition at line 80 of file G4BiasingProcessInterface.cc.

83 : G4VProcess( useThisName != "" ? useThisName : "biasWrapper("+wrappedProcess->GetProcessName()+")",
84 wrappedProcess->GetProcessType()),
85 fCurrentTrack ( nullptr ),
86 fPreviousStepSize (-1.0), fCurrentMinimumStep( -1.0 ), fProposedSafety ( -1.0),
87 fOccurenceBiasingOperation( nullptr ), fFinalStateBiasingOperation( nullptr ), fNonPhysicsBiasingOperation( nullptr ),
88 fPreviousOccurenceBiasingOperation( nullptr ), fPreviousFinalStateBiasingOperation( nullptr ), fPreviousNonPhysicsBiasingOperation( nullptr ),
89 fResetWrappedProcessInteractionLength ( false ),
90 fWrappedProcess ( wrappedProcess ),
91 fIsPhysicsBasedBiasing ( true ),
92 fWrappedProcessIsAtRest ( wrappedIsAtRest ),
93 fWrappedProcessIsAlong ( wrappedIsAlongStep ),
94 fWrappedProcessIsPost ( wrappedIsPostStep ),
95 fWrappedProcessPostStepGPIL ( -1.0 ),
96 fBiasingPostStepGPIL ( -1.0 ),
97 fWrappedProcessInteractionLength ( -1.0 ),
98 fWrappedProcessForceCondition ( NotForced ),
99 fBiasingForceCondition ( NotForced ),
100 fWrappedProcessAlongStepGPIL ( -1.0 ),
101 fBiasingAlongStepGPIL ( -1.0 ),
102 fWrappedProcessGPILSelection ( NotCandidateForSelection ),
103 fBiasingGPILSelection ( NotCandidateForSelection ),
104 fBiasingInteractionLaw ( nullptr ),
105 fPreviousBiasingInteractionLaw ( nullptr ),
106 fPhysicalInteractionLaw ( nullptr ),
107 fOccurenceBiasingParticleChange ( nullptr ),
108 fIamFirstGPIL ( false ),
109 fProcessManager ( nullptr ),
110 fSharedData ( nullptr )
111{
112 for (G4int i = 0 ; i < 8 ; i++) fFirstLastFlags[i] = false;
113 fResetInteractionLaws.Put( true );
114 fCommonStart.Put(true);
115 fCommonEnd.Put(true);
116 fDoCommonConfigure.Put(true);
117
118 SetProcessSubType(fWrappedProcess->GetProcessSubType());
119
120 // -- create physical interaction law:
121 fPhysicalInteractionLaw = new G4InteractionLawPhysical("PhysicalInteractionLawFor("+GetProcessName()+")");
122 // -- instantiate particle change wrapper for occurrence biaising:
123 fOccurenceBiasingParticleChange = new G4ParticleChangeForOccurenceBiasing("biasingPCfor"+GetProcessName());
124 // -- instantiate a "do nothing" particle change:
125 fDummyParticleChange = new G4ParticleChangeForNothing();
126}
G4ProcessType GetProcessType() const
void SetProcessSubType(G4int)
G4int GetProcessSubType() const
const G4String & GetProcessName() const

◆ ~G4BiasingProcessInterface()

G4BiasingProcessInterface::~G4BiasingProcessInterface ( )

Definition at line 130 of file G4BiasingProcessInterface.cc.

131{
132 if ( fPhysicalInteractionLaw != 0 ) delete fPhysicalInteractionLaw;
133 if ( fOccurenceBiasingParticleChange ) delete fOccurenceBiasingParticleChange;
134 if ( fDummyParticleChange ) delete fDummyParticleChange;
135}

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4BiasingProcessInterface::AlongStepDoIt ( const G4Track & track,
const G4Step & step )
virtual

Implements G4VProcess.

Definition at line 648 of file G4BiasingProcessInterface.cc.

650{
651
652 // ---------------------------------------
653 // -- case outside of volume with biasing:
654 // ---------------------------------------
655 if ( fSharedData->fCurrentBiasingOperator == 0 )
656 {
657 if ( fWrappedProcessIsAlong ) return fWrappedProcess->AlongStepDoIt(track, step);
658 else
659 {
660 fDummyParticleChange->Initialize( track );
661 return fDummyParticleChange;
662 }
663 }
664
665 // -----------------------------------
666 // -- case inside volume with biasing:
667 // -----------------------------------
668 if ( fWrappedProcessIsAlong ) fOccurenceBiasingParticleChange->SetWrappedParticleChange( fWrappedProcess->AlongStepDoIt(track, step) );
669 else
670 {
671 fOccurenceBiasingParticleChange->SetWrappedParticleChange ( 0 );
672 fOccurenceBiasingParticleChange->ProposeTrackStatus( track.GetTrackStatus() );
673 }
674 G4double weightForNonInteraction (1.0);
675 if ( fBiasingInteractionLaw != 0 )
676 {
677 weightForNonInteraction =
678 fPhysicalInteractionLaw->ComputeNonInteractionProbabilityAt(step.GetStepLength()) /
679 fBiasingInteractionLaw ->ComputeNonInteractionProbabilityAt(step.GetStepLength());
680
681 fOccurenceBiasingOperation->AlongMoveBy( this, &step, weightForNonInteraction );
682
683 if ( weightForNonInteraction <= 0. )
684 {
686 ed << " Negative non interaction weight : w_NI = " << weightForNonInteraction <<
687 " p_NI(phys) = " << fPhysicalInteractionLaw->ComputeNonInteractionProbabilityAt(step.GetStepLength()) <<
688 " p_NI(bias) = " << fBiasingInteractionLaw ->ComputeNonInteractionProbabilityAt(step.GetStepLength()) <<
689 " step length = " << step.GetStepLength() <<
690 " biasing interaction law = `" << fBiasingInteractionLaw->GetName() << "'" << G4endl;
691 G4Exception(" G4BiasingProcessInterface::AlongStepDoIt(...)",
692 "BIAS.GEN.04",
694 ed);
695 }
696
697 }
698
699 fOccurenceBiasingParticleChange->SetOccurenceWeightForNonInteraction( weightForNonInteraction );
700
701 return fOccurenceBiasingParticleChange;
702
703}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
virtual G4double ComputeNonInteractionProbabilityAt(G4double length) const
virtual void Initialize(const G4Track &track)
G4double GetStepLength() const
G4TrackStatus GetTrackStatus() const
const G4String & GetName() const
virtual G4double ComputeNonInteractionProbabilityAt(G4double length) const =0
virtual void AlongMoveBy(const G4BiasingProcessInterface *, const G4Step *, G4double)
void ProposeTrackStatus(G4TrackStatus status)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)=0

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 569 of file G4BiasingProcessInterface.cc.

574{
575
576 // -- for helper methods:
577 fCurrentMinimumStep = currentMinimumStep;
578 fProposedSafety = proposedSafety;
579
580
581 // -- initialization default case:
582 fWrappedProcessAlongStepGPIL = DBL_MAX;
583 *selection = NotCandidateForSelection;
584 // ---------------------------------------
585 // -- case outside of volume with biasing:
586 // ---------------------------------------
587 if ( fSharedData->fCurrentBiasingOperator == 0 )
588 {
589 if ( fWrappedProcessIsAlong ) fWrappedProcessAlongStepGPIL =
590 fWrappedProcess->AlongStepGetPhysicalInteractionLength(track,
591 previousStepSize,
592 currentMinimumStep,
593 proposedSafety,
594 selection);
595 return fWrappedProcessAlongStepGPIL;
596 }
597
598 // --------------------------------------------------------------------
599 // -- non-physics based biasing: no along operation expected (for now):
600 // --------------------------------------------------------------------
601 if ( !fIsPhysicsBasedBiasing ) return fWrappedProcessAlongStepGPIL;
602
603 // ----------------------
604 // -- physics-based case:
605 // ----------------------
606 if ( fOccurenceBiasingOperation == 0 )
607 {
608 if ( fWrappedProcessIsAlong ) fWrappedProcessAlongStepGPIL =
609 fWrappedProcess->AlongStepGetPhysicalInteractionLength(track,
610 previousStepSize,
611 currentMinimumStep,
612 proposedSafety,
613 selection);
614 return fWrappedProcessAlongStepGPIL;
615 }
616
617
618 // -----------------------------------------------------------
619 // -- From here we have an valid occurrence biasing operation:
620 // -----------------------------------------------------------
621 // -- Give operation opportunity to shorten step proposed by physics process:
622 fBiasingAlongStepGPIL = fOccurenceBiasingOperation->ProposeAlongStepLimit( this );
623 G4double minimumStep = fBiasingAlongStepGPIL < currentMinimumStep ? fBiasingAlongStepGPIL : currentMinimumStep ;
624 // -- wrapped process is called with minimum step ( <= currentMinimumStep passed ) : an along process can not
625 // -- have its operation stretched over what it expects:
626 if ( fWrappedProcessIsAlong )
627 {
628 fWrappedProcessAlongStepGPIL = fWrappedProcess->AlongStepGetPhysicalInteractionLength(track,
629 previousStepSize,
630 minimumStep,
631 proposedSafety,
632 selection);
633 fWrappedProcessGPILSelection = *selection;
634 fBiasingGPILSelection = fOccurenceBiasingOperation->ProposeGPILSelection( fWrappedProcessGPILSelection );
635 }
636 else
637 {
638 fBiasingGPILSelection = fOccurenceBiasingOperation->ProposeGPILSelection( NotCandidateForSelection );
639 fWrappedProcessAlongStepGPIL = fBiasingAlongStepGPIL;
640 }
641
642 *selection = fBiasingGPILSelection;
643
644 return fWrappedProcessAlongStepGPIL;
645
646}
virtual G4GPILSelection ProposeGPILSelection(const G4GPILSelection wrappedProcessSelection)
virtual G4double ProposeAlongStepLimit(const G4BiasingProcessInterface *)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
#define DBL_MAX
Definition templates.hh:62

◆ AtRestDoIt()

G4VParticleChange * G4BiasingProcessInterface::AtRestDoIt ( const G4Track & track,
const G4Step & step )
virtual

Implements G4VProcess.

Definition at line 711 of file G4BiasingProcessInterface.cc.

713{
714 return fWrappedProcess->AtRestDoIt(track, step);
715}
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)=0

◆ AtRestGetPhysicalInteractionLength()

G4double G4BiasingProcessInterface::AtRestGetPhysicalInteractionLength ( const G4Track & track,
G4ForceCondition * condition )
virtual

Implements G4VProcess.

Definition at line 706 of file G4BiasingProcessInterface.cc.

708{
709 return fWrappedProcess->AtRestGetPhysicalInteractionLength(track, condition);
710}
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)=0

◆ BuildPhysicsTable()

void G4BiasingProcessInterface::BuildPhysicsTable ( const G4ParticleDefinition & pd)
virtual

Reimplemented from G4VProcess.

Definition at line 741 of file G4BiasingProcessInterface.cc.

742{
743 // -- Sequential mode : called second (after PreparePhysicsTable(..))
744 // -- MT mode : called second (after PreparePhysicsTable(..)) by master thread.
745 // -- Corresponding process instance not used then by tracking.
746 // -- PreparePhysicsTable(...) has been called first for all processes,
747 // -- so the first/last flags and G4BiasingProcessInterface vector of processes have
748 // -- been properly setup, fIamFirstGPIL is valid.
749 if ( fWrappedProcess != 0 )
750 {
751 fWrappedProcess->BuildPhysicsTable(pd);
752 }
753
754 if ( fIamFirstGPIL )
755 {
756 // -- Re-order vector of processes to match that of the GPIL
757 // -- (made for fIamFirstGPIL, but important is to have it made once):
758 ReorderBiasingVectorAsGPIL();
759 // -- Let operators to configure themselves for the master thread or for sequential mode.
760 // -- Intended here is in particular the registration to physics model catalog.
761 // -- The fDoCommonConfigure is to ensure that this Configure is made by only one process (othewise each first process makes the call):
762 if ( fDoCommonConfigure.Get() )
763 {
764 for ( size_t optr = 0 ; optr < ( G4VBiasingOperator::GetBiasingOperators() ).size() ; optr ++)
765 ( G4VBiasingOperator::GetBiasingOperators() )[optr]->Configure( );
766 fDoCommonConfigure.Put(false);
767 }
768
769 }
770}
value_type & Get() const
Definition G4Cache.hh:315
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators()
virtual void BuildPhysicsTable(const G4ParticleDefinition &)

◆ BuildWorkerPhysicsTable()

void G4BiasingProcessInterface::BuildWorkerPhysicsTable ( const G4ParticleDefinition & pd)
virtual

Reimplemented from G4VProcess.

Definition at line 838 of file G4BiasingProcessInterface.cc.

839{
840 // -- Sequential mode : not called
841 // -- MT mode : called after PrepareWorkerPhysicsTable(..)
842 // -- PrepareWorkerPhysicsTable(...) has been called first for all processes,
843 // -- so the first/last flags and G4BiasingProcessInterface vector of processes have
844 // -- been properly setup, fIamFirstGPIL is valid.
845 if ( fWrappedProcess != 0 )
846 {
847 fWrappedProcess->BuildWorkerPhysicsTable(pd);
848 }
849
850 if ( fIamFirstGPIL )
851 {
852 // -- Re-order vector of processes to match that of the GPIL
853 // -- (made for fIamFirstGPIL, but important is to have it made once):
854 ReorderBiasingVectorAsGPIL();
855 // -- Let operators to configure themselves for the worker thread, if needed.
856 // -- Registration to physics model catalog **IS NOT** to be made here, but in Configure().
857 // -- The fDoCommonConfigure is to ensure that this Configure is made by only one process (othewise each first process makes the call):
858 if ( fDoCommonConfigure.Get() )
859 {
860 for ( size_t optr = 0 ; optr < ( G4VBiasingOperator::GetBiasingOperators() ).size() ; optr ++)
861 ( G4VBiasingOperator::GetBiasingOperators() )[optr]->ConfigureForWorker( );
862 fDoCommonConfigure.Put(false);
863 }
864 }
865}
virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition &part)

◆ EndTracking()

void G4BiasingProcessInterface::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 184 of file G4BiasingProcessInterface.cc.

185{
186 if ( fIsPhysicsBasedBiasing ) fWrappedProcess->EndTracking();
187 if ( fSharedData->fCurrentBiasingOperator) (fSharedData->fCurrentBiasingOperator)->ExitingBiasing( fCurrentTrack, this );
188 fBiasingInteractionLaw = 0;
189
190 // -- Inform operators of end of tracking:
191 if ( fCommonEnd.Get() )
192 {
193 fCommonEnd .Put( false );// = false;
194 fCommonStart.Put( true );// = true;
195
196 for ( size_t optr = 0 ; optr < ( G4VBiasingOperator::GetBiasingOperators() ).size() ; optr ++)
197 ( G4VBiasingOperator::GetBiasingOperators() )[optr]->EndTracking( );
198
199 // -- §§ for above loop, do as in StartTracking.
200 }
201}
virtual void EndTracking()

◆ GetAlongStepGPIL()

G4double G4BiasingProcessInterface::GetAlongStepGPIL ( ) const
inline

Definition at line 153 of file G4BiasingProcessInterface.hh.

153{ return fWrappedProcessAlongStepGPIL; }

Referenced by G4BOptnForceCommonTruncatedExp::ApplyFinalStateBiasing().

◆ GetBiasingProcessInterfaces()

const std::vector< const G4BiasingProcessInterface * > & G4BiasingProcessInterface::GetBiasingProcessInterfaces ( ) const
inline

Definition at line 106 of file G4BiasingProcessInterface.hh.

107 {return fSharedData-> fPublicBiasingProcessInterfaces;}

◆ GetCurrentBiasingOperator()

G4VBiasingOperator * G4BiasingProcessInterface::GetCurrentBiasingOperator ( ) const
inline

Definition at line 92 of file G4BiasingProcessInterface.hh.

92{return fSharedData-> fCurrentBiasingOperator; }

◆ GetCurrentFinalStateBiasingOperation()

G4VBiasingOperation * G4BiasingProcessInterface::GetCurrentFinalStateBiasingOperation ( ) const
inline

Definition at line 99 of file G4BiasingProcessInterface.hh.

99{ return fFinalStateBiasingOperation; }

◆ GetCurrentMinimumStep()

G4double G4BiasingProcessInterface::GetCurrentMinimumStep ( ) const
inline

Definition at line 148 of file G4BiasingProcessInterface.hh.

148{ return fCurrentMinimumStep;}

◆ GetCurrentNonPhysicsBiasingOperation()

G4VBiasingOperation * G4BiasingProcessInterface::GetCurrentNonPhysicsBiasingOperation ( ) const
inline

Definition at line 95 of file G4BiasingProcessInterface.hh.

95{ return fNonPhysicsBiasingOperation; }

◆ GetCurrentOccurenceBiasingOperation()

G4VBiasingOperation * G4BiasingProcessInterface::GetCurrentOccurenceBiasingOperation ( ) const
inline

Definition at line 97 of file G4BiasingProcessInterface.hh.

97{ return fOccurenceBiasingOperation; }

◆ GetIsFirstPostStepDoItInterface()

G4bool G4BiasingProcessInterface::GetIsFirstPostStepDoItInterface ( G4bool physOnly = true) const

Definition at line 902 of file G4BiasingProcessInterface.cc.

903{
904 G4int iPhys = ( physOnly ) ? 1 : 0;
905 return fFirstLastFlags[IdxFirstLast( 1, 0, iPhys)];
906}

Referenced by G4BOptnForceFreeFlight::ApplyFinalStateBiasing().

◆ GetIsFirstPostStepGPILInterface()

G4bool G4BiasingProcessInterface::GetIsFirstPostStepGPILInterface ( G4bool physOnly = true) const

Definition at line 888 of file G4BiasingProcessInterface.cc.

889{
890 G4int iPhys = ( physOnly ) ? 1 : 0;
891 return fFirstLastFlags[IdxFirstLast( 1, 1, iPhys)];
892}

◆ GetIsLastPostStepDoItInterface()

G4bool G4BiasingProcessInterface::GetIsLastPostStepDoItInterface ( G4bool physOnly = true) const

Definition at line 909 of file G4BiasingProcessInterface.cc.

910{
911 G4int iPhys = ( physOnly ) ? 1 : 0;
912 return fFirstLastFlags[IdxFirstLast( 0, 0, iPhys)];
913}

◆ GetIsLastPostStepGPILInterface()

G4bool G4BiasingProcessInterface::GetIsLastPostStepGPILInterface ( G4bool physOnly = true) const

Definition at line 895 of file G4BiasingProcessInterface.cc.

896{
897 G4int iPhys = ( physOnly ) ? 1 : 0;
898 return fFirstLastFlags[IdxFirstLast( 0, 1, iPhys)];
899}

◆ GetNonPhysicsBiasingProcessInterfaces()

const std::vector< const G4BiasingProcessInterface * > & G4BiasingProcessInterface::GetNonPhysicsBiasingProcessInterfaces ( ) const
inline

Definition at line 110 of file G4BiasingProcessInterface.hh.

111 {return fSharedData->fPublicNonPhysicsBiasingProcessInterfaces;}

◆ GetPhysicsBiasingProcessInterfaces()

const std::vector< const G4BiasingProcessInterface * > & G4BiasingProcessInterface::GetPhysicsBiasingProcessInterfaces ( ) const
inline

Definition at line 108 of file G4BiasingProcessInterface.hh.

109 {return fSharedData-> fPublicPhysicsBiasingProcessInterfaces;}

◆ GetPostStepGPIL()

G4double G4BiasingProcessInterface::GetPostStepGPIL ( ) const
inline

Definition at line 152 of file G4BiasingProcessInterface.hh.

152{ return fBiasingPostStepGPIL; }

Referenced by G4BOptnForceCommonTruncatedExp::ApplyFinalStateBiasing().

◆ GetPreviousBiasingOperator()

G4VBiasingOperator * G4BiasingProcessInterface::GetPreviousBiasingOperator ( ) const
inline

Definition at line 93 of file G4BiasingProcessInterface.hh.

93{return fSharedData->fPreviousBiasingOperator; }

◆ GetPreviousFinalStateBiasingOperation()

G4VBiasingOperation * G4BiasingProcessInterface::GetPreviousFinalStateBiasingOperation ( ) const
inline

Definition at line 100 of file G4BiasingProcessInterface.hh.

100{ return fPreviousFinalStateBiasingOperation; }

◆ GetPreviousNonPhysicsBiasingOperation()

G4VBiasingOperation * G4BiasingProcessInterface::GetPreviousNonPhysicsBiasingOperation ( ) const
inline

Definition at line 96 of file G4BiasingProcessInterface.hh.

96{ return fPreviousNonPhysicsBiasingOperation; }

◆ GetPreviousOccurenceBiasingOperation()

G4VBiasingOperation * G4BiasingProcessInterface::GetPreviousOccurenceBiasingOperation ( ) const
inline

Definition at line 98 of file G4BiasingProcessInterface.hh.

98{ return fPreviousOccurenceBiasingOperation; }

◆ GetPreviousStepSize()

G4double G4BiasingProcessInterface::GetPreviousStepSize ( ) const
inline

Definition at line 147 of file G4BiasingProcessInterface.hh.

147{ return fPreviousStepSize;}

◆ GetProcessManager()

const G4ProcessManager * G4BiasingProcessInterface::GetProcessManager ( )
virtual

Reimplemented from G4VProcess.

Definition at line 831 of file G4BiasingProcessInterface.cc.

832{
833 if ( fWrappedProcess != 0 ) return fWrappedProcess->GetProcessManager();
834 else return G4VProcess::GetProcessManager();
835}
virtual const G4ProcessManager * GetProcessManager()

◆ GetProposedSafety()

G4double G4BiasingProcessInterface::GetProposedSafety ( ) const
inline

Definition at line 149 of file G4BiasingProcessInterface.hh.

149{ return fProposedSafety;}

◆ GetSharedData() [1/2]

const G4BiasingProcessSharedData * G4BiasingProcessInterface::GetSharedData ( ) const
inline

Definition at line 114 of file G4BiasingProcessInterface.hh.

114{ return fSharedData; }

Referenced by G4BOptrForceCollision::ConfigureForWorker(), and G4ChannelingOptrChangeCrossSection::StartRun().

◆ GetSharedData() [2/2]

const G4BiasingProcessSharedData * G4BiasingProcessInterface::GetSharedData ( const G4ProcessManager * mgr)
static

Definition at line 138 of file G4BiasingProcessInterface.cc.

139{
141 G4BiasingProcessSharedData* >::const_iterator itr = G4BiasingProcessSharedData::fSharedDataMap.Find( mgr );
142 if ( itr != G4BiasingProcessSharedData::fSharedDataMap.End( ) )
143 {
144 return (*itr).second;
145 }
146 else return 0;
147}
iterator Find(const key_type &k)
Definition G4Cache.hh:458

◆ GetWrappedProcess()

◆ GetWrappedProcessIsAlong()

G4bool G4BiasingProcessInterface::GetWrappedProcessIsAlong ( ) const
inline

Definition at line 142 of file G4BiasingProcessInterface.hh.

142{ return fWrappedProcessIsAlong; }

◆ GetWrappedProcessIsAtRest()

G4bool G4BiasingProcessInterface::GetWrappedProcessIsAtRest ( ) const
inline

Definition at line 141 of file G4BiasingProcessInterface.hh.

141{ return fWrappedProcessIsAtRest; }

◆ GetWrappedProcessIsPost()

G4bool G4BiasingProcessInterface::GetWrappedProcessIsPost ( ) const
inline

Definition at line 143 of file G4BiasingProcessInterface.hh.

143{ return fWrappedProcessIsPost; }

◆ IsApplicable()

G4bool G4BiasingProcessInterface::IsApplicable ( const G4ParticleDefinition & pd)
virtual

Reimplemented from G4VProcess.

Definition at line 718 of file G4BiasingProcessInterface.cc.

719{
720 if ( fWrappedProcess != 0 ) return fWrappedProcess->IsApplicable(pd);
721 else return true;
722}
virtual G4bool IsApplicable(const G4ParticleDefinition &)

◆ IsFirstPostStepDoItInterface()

G4bool G4BiasingProcessInterface::IsFirstPostStepDoItInterface ( G4bool physOnly = true) const

Definition at line 976 of file G4BiasingProcessInterface.cc.

977{
978 G4bool isFirst = true;
979 const G4ProcessVector* pv = fProcessManager->GetPostStepProcessVector(typeDoIt);
980 G4int thisIdx(-1);
981 for (G4int i = 0; i < (G4int)pv->size(); ++i )
982 if ( (*pv)(i) == this ) { thisIdx = i; break; }
983 if ( thisIdx < 0 ) return false; // -- to ignore pure along processes
984 for ( std::size_t i = 0; i < (fSharedData->fBiasingProcessInterfaces).size(); ++i )
985 {
986 if ( (fSharedData->fBiasingProcessInterfaces)[i]->fIsPhysicsBasedBiasing || !physOnly )
987 {
988 G4int thatIdx(-1);
989 for (G4int j = 0; j < (G4int)pv->size(); ++j )
990 if ( (*pv)(j) == (fSharedData->fBiasingProcessInterfaces)[i] )
991 { thatIdx = j; break; }
992 if ( thatIdx >= 0 ) // -- to ignore pure along processes
993 {
994 if ( thisIdx > thatIdx )
995 {
996 isFirst = false;
997 break;
998 }
999 }
1000 }
1001 }
1002 return isFirst;
1003}
@ typeDoIt
bool G4bool
Definition G4Types.hh:86
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const

◆ IsFirstPostStepGPILInterface()

G4bool G4BiasingProcessInterface::IsFirstPostStepGPILInterface ( G4bool physOnly = true) const

Definition at line 916 of file G4BiasingProcessInterface.cc.

917{
918 G4bool isFirst = true;
919 const G4ProcessVector* pv = fProcessManager->GetPostStepProcessVector(typeGPIL);
920 G4int thisIdx(-1);
921 for (G4int i = 0; i < (G4int)pv->size(); ++i )
922 if ( (*pv)(i) == this ) { thisIdx = i; break; }
923 if ( thisIdx < 0 ) return false; // -- to ignore pure along processes
924 for ( std::size_t i = 0; i < (fSharedData->fBiasingProcessInterfaces).size(); ++i )
925 {
926 if ( (fSharedData->fBiasingProcessInterfaces)[i]->fIsPhysicsBasedBiasing || !physOnly )
927 {
928 G4int thatIdx(-1);
929 for (G4int j = 0; j < (G4int)pv->size(); ++j )
930 if ( (*pv)(j) == (fSharedData->fBiasingProcessInterfaces)[i] )
931 { thatIdx = j; break; }
932 if ( thatIdx >= 0 ) // -- to ignore pure along processes
933 {
934 if ( thisIdx > thatIdx )
935 {
936 isFirst = false;
937 break;
938 }
939 }
940 }
941 }
942 return isFirst;
943}
@ typeGPIL

◆ IsLastPostStepDoItInterface()

G4bool G4BiasingProcessInterface::IsLastPostStepDoItInterface ( G4bool physOnly = true) const

Definition at line 1006 of file G4BiasingProcessInterface.cc.

1007{
1008 G4bool isLast = true;
1009 const G4ProcessVector* pv = fProcessManager->GetPostStepProcessVector(typeDoIt);
1010 G4int thisIdx(-1);
1011 for (G4int i = 0; i < (G4int)pv->size(); ++i )
1012 if ( (*pv)(i) == this ) { thisIdx = i; break; }
1013 if ( thisIdx < 0 ) return false; // -- to ignore pure along processes
1014 for ( std::size_t i = 0; i < (fSharedData->fBiasingProcessInterfaces).size(); ++i )
1015 {
1016 if ( (fSharedData->fBiasingProcessInterfaces)[i]->fIsPhysicsBasedBiasing || !physOnly )
1017 {
1018 G4int thatIdx(-1);
1019 for (G4int j = 0; j < (G4int)pv->size(); ++j )
1020 if ( (*pv)(j) == (fSharedData->fBiasingProcessInterfaces)[i] )
1021 { thatIdx = j; break; }
1022 if ( thatIdx >= 0 ) // -- to ignore pure along processes
1023 {
1024 if ( thisIdx < thatIdx )
1025 {
1026 isLast = false;
1027 break;
1028 }
1029 }
1030 }
1031 }
1032 return isLast;
1033}

◆ IsLastPostStepGPILInterface()

G4bool G4BiasingProcessInterface::IsLastPostStepGPILInterface ( G4bool physOnly = true) const

Definition at line 946 of file G4BiasingProcessInterface.cc.

947{
948 G4bool isLast = true;
949 const G4ProcessVector* pv = fProcessManager->GetPostStepProcessVector(typeGPIL);
950 G4int thisIdx(-1);
951 for (G4int i = 0; i < (G4int)pv->size(); ++i )
952 if ( (*pv)(i) == this ) { thisIdx = i; break; }
953 if ( thisIdx < 0 ) return false; // -- to ignore pure along processes
954 for ( std::size_t i = 0; i < (fSharedData->fBiasingProcessInterfaces).size(); ++i )
955 {
956 if ( (fSharedData->fBiasingProcessInterfaces)[i]->fIsPhysicsBasedBiasing || !physOnly )
957 {
958 G4int thatIdx(-1);
959 for (G4int j = 0; j < (G4int)pv->size(); ++j )
960 if ( (*pv)(j) == (fSharedData->fBiasingProcessInterfaces)[i] )
961 { thatIdx = j; break; }
962 if ( thatIdx >= 0 ) // -- to ignore pure along processes
963 {
964 if ( thisIdx < thatIdx )
965 {
966 isLast = false;
967 break;
968 }
969 }
970 }
971 }
972 return isLast;
973}

◆ PostStepDoIt()

G4VParticleChange * G4BiasingProcessInterface::PostStepDoIt ( const G4Track & track,
const G4Step & step )
virtual

Implements G4VProcess.

Definition at line 447 of file G4BiasingProcessInterface.cc.

449{
450 // ---------------------------------------
451 // -- case outside of volume with biasing:
452 // ---------------------------------------
453 if ( fSharedData->fCurrentBiasingOperator == 0 ) return fWrappedProcess->PostStepDoIt(track, step);
454
455 // ----------------------------
456 // -- non-physics biasing case:
457 // ----------------------------
458 if ( !fIsPhysicsBasedBiasing )
459 {
460 G4VParticleChange* particleChange = fNonPhysicsBiasingOperation->GenerateBiasingFinalState( &track, &step );
461 (fSharedData->fCurrentBiasingOperator)->ReportOperationApplied( this, BAC_NonPhysics, fNonPhysicsBiasingOperation, particleChange );
462 return particleChange;
463 }
464
465 // -- physics biasing case:
466 // ------------------------
467 // -- It proceeds with the following logic:
468 // -- 1) Obtain the final state
469 // -- This final state may be analog or biased.
470 // -- The biased final state is obtained through a biasing operator
471 // -- returned by the operator.
472 // -- 2) The biased final state may be asked to be "force as it is"
473 // -- in what case the particle change is returned as is to the
474 // -- stepping.
475 // -- In all other cases (analog final state or biased final but
476 // -- not forced) the final state weight may be modified by the
477 // -- occurrence biasing, if such an occurrence biasing is at play.
478
479 // -- Get final state, biased or analog:
480 G4VParticleChange* finalStateParticleChange;
482 fFinalStateBiasingOperation = (fSharedData->fCurrentBiasingOperator)->GetProposedFinalStateBiasingOperation( &track, this );
483 // -- Flag below is to force the biased generated particle change to be returned "as is" to the stepping, disregarding there
484 // -- was or not a occurrence biasing that would apply. Weight relevance under full responsibility of the biasing operation.
485 G4bool forceBiasedFinalState = false;
486 if ( fFinalStateBiasingOperation != 0 )
487 {
488 finalStateParticleChange = fFinalStateBiasingOperation->ApplyFinalStateBiasing( this, &track, &step, forceBiasedFinalState );
489 BAC = BAC_FinalState;
490 }
491 else
492 {
493 finalStateParticleChange = fWrappedProcess->PostStepDoIt(track, step);
494 BAC = BAC_None ;
495 }
496
497 // -- if no occurrence biasing operation, we're done:
498 if ( fOccurenceBiasingOperation == 0 )
499 {
500 (fSharedData->fCurrentBiasingOperator)->ReportOperationApplied( this, BAC, fFinalStateBiasingOperation, finalStateParticleChange );
501 return finalStateParticleChange;
502 }
503
504 // -- if biased final state has been asked to be forced, we're done:
505 if ( forceBiasedFinalState )
506 {
507 (fSharedData->fCurrentBiasingOperator)->ReportOperationApplied( this, BAC, fFinalStateBiasingOperation, finalStateParticleChange );
508 return finalStateParticleChange;
509 }
510
511
512 // -- If occurrence biasing, applies the occurrence biasing weight correction on top of final state (biased or not):
513 G4double weightForInteraction = 1.0;
514 if ( !fBiasingInteractionLaw->IsSingular() ) weightForInteraction =
515 fPhysicalInteractionLaw->ComputeEffectiveCrossSectionAt(step.GetStepLength()) /
516 fBiasingInteractionLaw ->ComputeEffectiveCrossSectionAt(step.GetStepLength());
517 else
518 {
519 // -- at this point effective XS can only be infinite, if not, there is a logic problem
520 if ( !fBiasingInteractionLaw->IsEffectiveCrossSectionInfinite() )
521 {
523 ed << "Internal inconsistency in cross-section handling. Please report !" << G4endl;
524 G4Exception(" G4BiasingProcessInterface::PostStepDoIt(...)",
525 "BIAS.GEN.02",
527 ed);
528 // -- if XS is infinite, weight is zero (and will stay zero), but we'll do differently.
529 // -- Should foresee in addition something to remember that in case of singular
530 // -- distribution, weight can only be partly calculated
531 }
532 }
533
534 if ( weightForInteraction <= 0. )
535 {
537 ed << " Negative interaction weight : w_I = "
538 << weightForInteraction <<
539 " XS_I(phys) = " << fBiasingInteractionLaw ->ComputeEffectiveCrossSectionAt(step.GetStepLength()) <<
540 " XS_I(bias) = " << fPhysicalInteractionLaw->ComputeEffectiveCrossSectionAt(step.GetStepLength()) <<
541 " step length = " << step.GetStepLength() <<
542 " Interaction law = `" << fBiasingInteractionLaw << "'" <<
543 G4endl;
544 G4Exception(" G4BiasingProcessInterface::PostStepDoIt(...)",
545 "BIAS.GEN.03",
547 ed);
548
549 }
550
551 (fSharedData->fCurrentBiasingOperator)->ReportOperationApplied( this, BAC,
552 fOccurenceBiasingOperation, weightForInteraction,
553 fFinalStateBiasingOperation, finalStateParticleChange );
554
555
556 fOccurenceBiasingParticleChange->SetOccurenceWeightForInteraction( weightForInteraction );
557 fOccurenceBiasingParticleChange->SetSecondaryWeightByProcess( true );
558 fOccurenceBiasingParticleChange->SetWrappedParticleChange( finalStateParticleChange );
559 fOccurenceBiasingParticleChange->ProposeTrackStatus( finalStateParticleChange->GetTrackStatus() );
560 fOccurenceBiasingParticleChange->StealSecondaries(); // -- this also makes weightForInteraction applied to secondaries stolen
561
562 // -- finish:
563 return fOccurenceBiasingParticleChange;
564
565}
G4BiasingAppliedCase
@ BAC_FinalState
@ BAC_NonPhysics
virtual G4double ComputeEffectiveCrossSectionAt(G4double length) const
virtual G4bool IsEffectiveCrossSectionInfinite() const
virtual G4bool IsSingular() const
virtual G4double ComputeEffectiveCrossSectionAt(G4double length) const =0
virtual G4VParticleChange * ApplyFinalStateBiasing(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)=0
virtual G4VParticleChange * GenerateBiasingFinalState(const G4Track *, const G4Step *)=0
void SetSecondaryWeightByProcess(G4bool)
G4TrackStatus GetTrackStatus() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0

◆ PostStepGetPhysicalInteractionLength()

G4double G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
virtual

Implements G4VProcess.

Definition at line 205 of file G4BiasingProcessInterface.cc.

208{
209
210 // ---------------------------------------------------------------------------------------------------
211 // -- The "biasing process master" takes care of updating the biasing operator, and for all biasing
212 // -- processes it invokes the PostStepGPIL of physical wrapped processes (anticipate stepping manager
213 // -- call ! ) to make all cross-sections updated with current step, and hence available before the
214 // -- first call to the biasing operator.
215 // ---------------------------------------------------------------------------------------------------
216 if ( fIamFirstGPIL )
217 {
218 // -- Update previous biasing operator, and assume the operator stays the same by
219 // -- default and that it is not left at the beginning of this step. These
220 // -- assumptions might be wrong if there is a volume change (in paralllel or
221 // -- mass geometries) in what case the flags will be updated.
222 fSharedData->fPreviousBiasingOperator = fSharedData->fCurrentBiasingOperator;
223 fSharedData->fIsNewOperator = false;
224 fSharedData->fLeavingPreviousOperator = false;
225 // -- If new volume, either in mass or parallel geometries, get possible new biasing operator:
226 // -------------------------------------------------------------------------------------------
227 // -- Get biasing operator in parallel geometries:
228 G4bool firstStepInParallelVolume = false;
229 if ( fSharedData->fParallelGeometriesLimiterProcess )
230 {
231 G4VBiasingOperator* newParallelOperator( nullptr );
232 G4bool firstStep = ( track.GetCurrentStepNumber() == 1 );
233 size_t iParallel = 0;
234 for ( auto wasLimiting : fSharedData->fParallelGeometriesLimiterProcess->GetWasLimiting() )
235 {
236 if ( firstStep || wasLimiting )
237 {
238 firstStepInParallelVolume = true;
239
240 auto tmpParallelOperator = G4VBiasingOperator::GetBiasingOperator( (fSharedData->fParallelGeometriesLimiterProcess->GetCurrentVolumes()[iParallel])
241 ->GetLogicalVolume() );
242 if ( newParallelOperator )
243 {
244 if ( tmpParallelOperator )
245 {
247 ed << " Several biasing operators are defined at the same place in parallel geometries ! Found:\n";
248 ed << " - `" << newParallelOperator->GetName() << "' and \n";
249 ed << " - `" << tmpParallelOperator->GetName() << "'.\n";
250 ed << " Keeping `" << newParallelOperator->GetName() << "'. Behavior not guaranteed ! Please consider having only one operator at a place. " << G4endl;
251 G4Exception(" G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength(...)",
252 "BIAS.GEN.30",
254 ed);
255 }
256 }
257 else newParallelOperator = tmpParallelOperator;
258 }
259 iParallel++;
260 }
261 fSharedData->fParallelGeometryOperator = newParallelOperator;
262 } // -- end of " if ( fSharedData->fParallelGeometriesLimiterProcess )"
263
264 // -- Get biasing operator in mass geometry:
265 // -- [§§ Note : bug with this first step ? Does not work if previous step was concurrently limited with geometry. Might make use of safety at last point ?]
266 G4bool firstStepInVolume = ( (track.GetStep()->GetPreStepPoint()->GetStepStatus() == fGeomBoundary) || (track.GetCurrentStepNumber() == 1) );
267 // fSharedData->fIsNewOperator = false;
268 // fSharedData->fLeavingPreviousOperator = false;
269 if ( firstStepInVolume )
270 {
272 fSharedData->fMassGeometryOperator = newOperator;
273 if ( ( newOperator != nullptr ) && ( fSharedData->fParallelGeometryOperator != nullptr ) )
274 {
276 ed << " Biasing operators are defined at the same place in mass and parallel geometries ! Found:\n";
277 ed << " - `" << fSharedData->fParallelGeometryOperator->GetName() << "' in parallel geometry and \n";
278 ed << " - `" << newOperator->GetName() << "' in mass geometry.\n";
279 ed << " Keeping `" << fSharedData->fParallelGeometryOperator->GetName() << "'. Behavior not guaranteed ! Please consider having only one operator at a place. " << G4endl;
280 G4Exception(" G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength(...)",
281 "BIAS.GEN.31",
283 ed);
284 }
285 }
286
287 // -- conclude the operator selection, giving priority to parallel geometry (as told in exception message BIAS.GEN.30):
288 if ( firstStepInVolume || firstStepInParallelVolume )
289 {
290 G4VBiasingOperator* newOperator = fSharedData->fParallelGeometryOperator;
291 if ( newOperator == nullptr ) newOperator = fSharedData->fMassGeometryOperator;
292
293 fSharedData->fCurrentBiasingOperator = newOperator ;
294
295 if ( newOperator != fSharedData->fPreviousBiasingOperator )
296 {
297 fSharedData->fLeavingPreviousOperator = ( fSharedData->fPreviousBiasingOperator != nullptr ) ;
298 fSharedData->fIsNewOperator = ( newOperator != nullptr );
299 }
300 }
301
302
303 // -- calls to wrapped process PostStepGPIL's:
304 // -------------------------------------------
305 // -- Each physics wrapper process has its
306 // -- fWrappedProcessPostStepGPIL ,
307 // -- fWrappedProcessForceCondition ,
308 // -- fWrappedProcessInteractionLength
309 // -- updated.
310 if ( fSharedData->fCurrentBiasingOperator != nullptr )
311 {
312 for ( size_t i = 0 ; i < (fSharedData->fPhysicsBiasingProcessInterfaces).size(); i++ )
313 (fSharedData->fPhysicsBiasingProcessInterfaces)[i]->InvokeWrappedProcessPostStepGPIL( track, previousStepSize, condition );
314 }
315 } // -- end of "if ( fIamFirstGPIL )"
316
317
318
319 // -- Remember previous operator and proposed operations, if any, and reset:
320 // -------------------------------------------------------------------------
321 // -- remember only in case some biasing might be called
322 if ( ( fSharedData->fPreviousBiasingOperator != 0 ) ||
323 ( fSharedData->fCurrentBiasingOperator != 0 ) )
324 {
325 fPreviousOccurenceBiasingOperation = fOccurenceBiasingOperation;
326 fPreviousFinalStateBiasingOperation = fFinalStateBiasingOperation;
327 fPreviousNonPhysicsBiasingOperation = fNonPhysicsBiasingOperation;
328 fPreviousBiasingInteractionLaw = fBiasingInteractionLaw;
329 // -- reset:
330 fOccurenceBiasingOperation = 0;
331 fFinalStateBiasingOperation = 0;
332 fNonPhysicsBiasingOperation = 0;
333 fBiasingInteractionLaw = 0;
334 // -- Physics PostStep and AlongStep GPIL
335 // fWrappedProcessPostStepGPIL : updated by InvokeWrappedProcessPostStepGPIL(...) above
336 fBiasingPostStepGPIL = DBL_MAX;
337 // fWrappedProcessInteractionLength : updated by InvokeWrappedProcessPostStepGPIL(...) above; inverse of analog cross-section.
338 // fWrappedProcessForceCondition : updated by InvokeWrappedProcessPostStepGPIL(...) above
339 fBiasingForceCondition = NotForced;
340 fWrappedProcessAlongStepGPIL = DBL_MAX;
341 fBiasingAlongStepGPIL = DBL_MAX;
342 fWrappedProcessGPILSelection = NotCandidateForSelection;
343 fBiasingGPILSelection = NotCandidateForSelection;
344 // -- for helper:
345 fPreviousStepSize = previousStepSize;
346 }
347
348
349 // -- previous step size value; it is switched to zero if resetting a wrapped process:
350 // -- (same trick used than in InvokedWrappedProcessPostStepGPIL )
351 G4double usedPreviousStepSize = previousStepSize;
352
353 // ----------------------------------------------
354 // -- If leaving a biasing operator, let it know:
355 // ----------------------------------------------
356 if ( fSharedData->fLeavingPreviousOperator )
357 {
358 (fSharedData->fPreviousBiasingOperator)->ExitingBiasing( &track, this );
359 // -- if no further biasing operator, reset process behavior to standard tracking:
360 if ( fSharedData->fCurrentBiasingOperator == 0 )
361 {
362 ResetForUnbiasedTracking();
363 if ( fIsPhysicsBasedBiasing )
364 {
365 // -- if the physics process has been under occurrence biasing, reset it:
366 if ( fResetWrappedProcessInteractionLength )
367 {
368 fResetWrappedProcessInteractionLength = false;
369 fWrappedProcess->ResetNumberOfInteractionLengthLeft();
370 // -- We set "previous step size" as 0.0, to let the process believe this is first step:
371 usedPreviousStepSize = 0.0;
372 }
373 }
374 }
375 }
376
377
378 // --------------------------------------------------------------
379 // -- no operator : analog tracking if physics-based, or nothing:
380 // --------------------------------------------------------------
381 if ( fSharedData->fCurrentBiasingOperator == 0 )
382 {
383 // -- take note of the "usedPreviousStepSize" value:
384 if ( fIsPhysicsBasedBiasing ) return fWrappedProcess->PostStepGetPhysicalInteractionLength(track, usedPreviousStepSize, condition);
385 else
386 {
388 return DBL_MAX;
389 }
390 }
391
392 // --------------------------------------------------
393 // -- A biasing operator exists. Proceed with
394 // -- treating non-physics and physics biasing cases:
395 //---------------------------------------------------
396
397 // -- non-physics-based biasing case:
398 // ----------------------------------
399 if ( !fIsPhysicsBasedBiasing )
400 {
401 fNonPhysicsBiasingOperation = (fSharedData->fCurrentBiasingOperator)->GetProposedNonPhysicsBiasingOperation( &track, this );
402 if ( fNonPhysicsBiasingOperation == 0 )
403 {
405 return DBL_MAX;
406 }
407 return fNonPhysicsBiasingOperation->DistanceToApplyOperation(&track, previousStepSize, condition);
408 }
409
410
411 // -- Physics based biasing case:
412 // ------------------------------
413 // -- Ask for possible GPIL biasing operation:
414 fOccurenceBiasingOperation = (fSharedData->fCurrentBiasingOperator)->GetProposedOccurenceBiasingOperation( &track, this );
415
416
417 // -- no operation for occurrence biasing, analog GPIL returns the wrapped process GPIL and condition values
418 if ( fOccurenceBiasingOperation == 0 )
419 {
420 *condition = fWrappedProcessForceCondition;
421 return fWrappedProcessPostStepGPIL;
422 }
423
424 // -- A valid GPIL biasing operation has been proposed:
425 // -- 0) remember wrapped process will need to be reset on biasing exit, if particle survives:
426 fResetWrappedProcessInteractionLength = true;
427 // -- 1) update process interaction length for reference analog interaction law ( fWrappedProcessInteractionLength updated/collected above):
428 fPhysicalInteractionLaw->SetPhysicalCrossSection( 1.0 / fWrappedProcessInteractionLength );
429 // -- 2) Collect biasing interaction law:
430 // -- The interaction law pointer is collected as a const pointer to the interaction law object.
431 // -- This interaction law will be kept under control of the biasing operation, which is the only
432 // -- entity that will change the state of the biasing interaction law.
433 // -- The force condition for biasing is asked at the same time, passing the analog one as default:
434 fBiasingForceCondition = fWrappedProcessForceCondition;
435 fBiasingInteractionLaw = fOccurenceBiasingOperation->ProvideOccurenceBiasingInteractionLaw( this, fBiasingForceCondition );
436 // -- 3) Ask operation to sample the biasing interaction law:
437 fBiasingPostStepGPIL = fBiasingInteractionLaw->GetSampledInteractionLength();
438
439 // -- finish
440 *condition = fBiasingForceCondition;
441 return fBiasingPostStepGPIL;
442
443}
@ fGeomBoundary
void SetPhysicalCrossSection(G4double crossSection)
const std::vector< const G4VPhysicalVolume * > & GetCurrentVolumes() const
const std::vector< G4bool > & GetWasLimiting() const
G4VPhysicalVolume * GetVolume() const
G4int GetCurrentStepNumber() const
const G4Step * GetStep() const
virtual G4double DistanceToApplyOperation(const G4Track *, G4double, G4ForceCondition *)=0
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *, G4ForceCondition &)=0
const G4String GetName() const
static G4VBiasingOperator * GetBiasingOperator(const G4LogicalVolume *)
G4LogicalVolume * GetLogicalVolume() const
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
virtual void ResetNumberOfInteractionLengthLeft()
Definition G4VProcess.cc:80

◆ PreparePhysicsTable()

void G4BiasingProcessInterface::PreparePhysicsTable ( const G4ParticleDefinition & pd)
virtual

Reimplemented from G4VProcess.

Definition at line 773 of file G4BiasingProcessInterface.cc.

774{
775 // -- Sequential mode : called first (before BuildPhysicsTable(..))
776 // -- MT mode : called first (before BuildPhysicsTable(..)) by master thread.
777 // -- Corresponding process instance not used then by tracking.
778 // -- Let process finding its first/last position in the process manager:
779 SetUpFirstLastFlags();
780 if ( fWrappedProcess != 0 )
781 {
782 fWrappedProcess->PreparePhysicsTable(pd);
783 }
784}
virtual void PreparePhysicsTable(const G4ParticleDefinition &)

◆ PrepareWorkerPhysicsTable()

void G4BiasingProcessInterface::PrepareWorkerPhysicsTable ( const G4ParticleDefinition & pd)
virtual

Reimplemented from G4VProcess.

Definition at line 868 of file G4BiasingProcessInterface.cc.

869{
870 // -- Sequential mode : not called
871 // -- MT mode : called first, before BuildWorkerPhysicsTable(..)
872 // -- Let process finding its first/last position in the process manager:
873 SetUpFirstLastFlags();
874
875 if ( fWrappedProcess != 0 )
876 {
877 fWrappedProcess->PrepareWorkerPhysicsTable(pd);
878 }
879}
virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition &)

◆ ResetNumberOfInteractionLengthLeft()

void G4BiasingProcessInterface::ResetNumberOfInteractionLengthLeft ( )
virtual

Reimplemented from G4VProcess.

Definition at line 882 of file G4BiasingProcessInterface.cc.

883{
884 if ( fWrappedProcess != 0 ) fWrappedProcess->ResetNumberOfInteractionLengthLeft();
885}

◆ RetrievePhysicsTable()

G4bool G4BiasingProcessInterface::RetrievePhysicsTable ( const G4ParticleDefinition * pd,
const G4String & s,
G4bool f )
virtual

Reimplemented from G4VProcess.

Definition at line 794 of file G4BiasingProcessInterface.cc.

795{
796 if ( fWrappedProcess != 0 ) return fWrappedProcess->RetrievePhysicsTable(pd, s, f);
797 else return false;
798}
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)

◆ SetMasterProcess()

void G4BiasingProcessInterface::SetMasterProcess ( G4VProcess * masterP)
virtual

Reimplemented from G4VProcess.

Definition at line 725 of file G4BiasingProcessInterface.cc.

726{
727 // -- Master for this process:
729 // -- Master for wrapped process:
730 if ( fWrappedProcess != 0 )
731 {
732 const G4BiasingProcessInterface* thisWrapperMaster = (const G4BiasingProcessInterface *)GetMasterProcess();
733 // -- paranoia check: (?)
734 G4VProcess* wrappedMaster = 0;
735 wrappedMaster = thisWrapperMaster->GetWrappedProcess();
736 fWrappedProcess->SetMasterProcess( wrappedMaster );
737 }
738}
const G4VProcess * GetMasterProcess() const
virtual void SetMasterProcess(G4VProcess *masterP)

◆ SetProcessManager()

void G4BiasingProcessInterface::SetProcessManager ( const G4ProcessManager * mgr)
virtual

Reimplemented from G4VProcess.

Definition at line 801 of file G4BiasingProcessInterface.cc.

802{
803 if ( fWrappedProcess != 0 ) fWrappedProcess->SetProcessManager(mgr);
805
806 // -- initialize fSharedData pointer:
807 if ( G4BiasingProcessSharedData::fSharedDataMap.Find(mgr) == G4BiasingProcessSharedData::fSharedDataMap.End() )
808 {
809 fSharedData = new G4BiasingProcessSharedData( mgr );
810 G4BiasingProcessSharedData::fSharedDataMap[mgr] = fSharedData;
811 }
812 else fSharedData = G4BiasingProcessSharedData::fSharedDataMap[mgr] ;
813 // -- augment list of co-operating processes:
814 fSharedData-> fBiasingProcessInterfaces.push_back( this );
815 fSharedData-> fPublicBiasingProcessInterfaces.push_back( this );
816 if ( fIsPhysicsBasedBiasing )
817 {
818 fSharedData-> fPhysicsBiasingProcessInterfaces.push_back( this );
819 fSharedData-> fPublicPhysicsBiasingProcessInterfaces.push_back( this );
820 }
821 else
822 {
823 fSharedData-> fNonPhysicsBiasingProcessInterfaces.push_back( this );
824 fSharedData-> fPublicNonPhysicsBiasingProcessInterfaces.push_back( this );
825 }
826 // -- remember process manager:
827 fProcessManager = mgr;
828}
virtual void SetProcessManager(const G4ProcessManager *)

◆ SetProposedSafety()

void G4BiasingProcessInterface::SetProposedSafety ( G4double sft)
inline

Definition at line 150 of file G4BiasingProcessInterface.hh.

150{ fProposedSafety = sft;}

◆ StartTracking()

void G4BiasingProcessInterface::StartTracking ( G4Track * track)
virtual

Reimplemented from G4VProcess.

Definition at line 150 of file G4BiasingProcessInterface.cc.

151{
152 fCurrentTrack = track;
153 if ( fIsPhysicsBasedBiasing ) fWrappedProcess->StartTracking(fCurrentTrack);
154 fOccurenceBiasingOperation = 0;
155 fPreviousOccurenceBiasingOperation = 0;
156 fFinalStateBiasingOperation = 0;
157 fPreviousFinalStateBiasingOperation = 0;
158 fNonPhysicsBiasingOperation = 0;
159 fPreviousNonPhysicsBiasingOperation = 0;
160 fBiasingInteractionLaw = 0;
161 fPreviousBiasingInteractionLaw = 0;
162
163 fPreviousStepSize = -1.0;
164
165 fResetWrappedProcessInteractionLength = false;
166
167 if ( fCommonStart.Get() )
168 {
169 fCommonStart.Put( false );// = false;
170 fCommonEnd .Put( true );// = true;
171
172 fSharedData-> fCurrentBiasingOperator = 0;
173 fSharedData->fPreviousBiasingOperator = 0;
174
175 // -- §§ Add a "fSharedData->nStarting" here and outside bracket "fSharedData->nStarting++" and " if (fSharedData->nStarting) == fSharedData->(vector interface length)"
176 // -- §§ call to the loop "StartTracking" of operators"
177
178 for ( size_t optr = 0 ; optr < ( G4VBiasingOperator::GetBiasingOperators() ).size() ; optr ++)
179 ( G4VBiasingOperator::GetBiasingOperators() )[optr]->StartTracking( fCurrentTrack );
180 }
181}
virtual void StartTracking(const G4Track *)
virtual void StartTracking(G4Track *)
Definition G4VProcess.cc:87

◆ StorePhysicsTable()

G4bool G4BiasingProcessInterface::StorePhysicsTable ( const G4ParticleDefinition * pd,
const G4String & s,
G4bool f )
virtual

Reimplemented from G4VProcess.

Definition at line 787 of file G4BiasingProcessInterface.cc.

788{
789 if ( fWrappedProcess != 0 ) return fWrappedProcess->StorePhysicsTable(pd, s, f);
790 else return false;
791}
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)

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