Geant4-11
Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
G4BOptrForceCollision Class Reference

#include <G4BOptrForceCollision.hh>

Inheritance diagram for G4BOptrForceCollision:
G4VBiasingOperator

Public Member Functions

void AttachTo (const G4LogicalVolume *)
 
virtual void Configure () final
 
virtual void ConfigureForWorker () final
 
virtual void EndTracking () final
 
virtual void ExitBiasing (const G4Track *, const G4BiasingProcessInterface *) final
 
void ExitingBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
 G4BOptrForceCollision (const G4ParticleDefinition *particleToForce, G4String name="ForceCollision")
 
 G4BOptrForceCollision (G4String particleToForce, G4String name="ForceCollision")
 
const G4String GetName () const
 
G4BiasingAppliedCase GetPreviousBiasingAppliedCase () const
 
const G4VBiasingOperationGetPreviousNonPhysicsAppliedOperation ()
 
G4VBiasingOperationGetProposedFinalStateBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
G4VBiasingOperationGetProposedNonPhysicsBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
G4VBiasingOperationGetProposedOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced) final
 
void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced) final
 
void ReportOperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 
void ReportOperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void StartRun () final
 
virtual void StartTracking (const G4Track *track) final
 
 ~G4BOptrForceCollision ()
 

Static Public Member Functions

static G4VBiasingOperatorGetBiasingOperator (const G4LogicalVolume *)
 
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators ()
 

Private Member Functions

virtual G4VBiasingOperationProposeFinalStateBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess) final
 
virtual G4VBiasingOperationProposeNonPhysicsBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess) final
 
virtual G4VBiasingOperationProposeOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess) final
 

Private Attributes

G4BOptnCloningfCloningOperation
 
const G4TrackfCurrentTrack
 
G4BOptrForceCollisionTrackDatafCurrentTrackData
 
std::map< const G4LogicalVolume *, G4intfDepthInTree
 
G4VBiasingOperationfFinalStateBiasingOperation
 
G4int fForceCollisionModelID
 
std::map< const G4BiasingProcessInterface *, G4BOptnForceFreeFlight * > fFreeFlightOperations
 
G4double fInitialTrackWeight
 
const G4String fName
 
G4VBiasingOperationfNonPhysicsBiasingOperation
 
G4VBiasingOperationfOccurenceBiasingOperation
 
const G4ParticleDefinitionfParticleToBias
 
const G4VBiasingOperationfPreviousAppliedFinalStateBiasingOperation
 
const G4VBiasingOperationfPreviousAppliedNonPhysicsBiasingOperation
 
const G4VBiasingOperationfPreviousAppliedOccurenceBiasingOperation
 
G4BiasingAppliedCase fPreviousBiasingAppliedCase
 
const G4VBiasingOperationfPreviousProposedFinalStateBiasingOperation
 
const G4VBiasingOperationfPreviousProposedNonPhysicsBiasingOperation
 
const G4VBiasingOperationfPreviousProposedOccurenceBiasingOperation
 
std::vector< const G4LogicalVolume * > fRootVolumes
 
G4bool fSetup
 
G4BOptnForceCommonTruncatedExpfSharedForceInteractionOperation
 

Static Private Attributes

static G4MapCache< const G4LogicalVolume *, G4VBiasingOperator * > fLogicalToSetupMap
 
static G4VectorCache< G4VBiasingOperator * > fOperators
 
static G4Cache< G4BiasingOperatorStateNotifier * > fStateNotifier
 

Detailed Description

Definition at line 58 of file G4BOptrForceCollision.hh.

Constructor & Destructor Documentation

◆ G4BOptrForceCollision() [1/2]

G4BOptrForceCollision::G4BOptrForceCollision ( G4String  particleToForce,
G4String  name = "ForceCollision" 
)

Definition at line 46 of file G4BOptrForceCollision.cc.

48 fForceCollisionModelID(G4PhysicsModelCatalog::GetModelID("model_GenBiasForceCollision")),
49 fCurrentTrack(nullptr),
50 fCurrentTrackData(nullptr),
52 fSetup(true)
53{
55 fCloningOperation = new G4BOptnCloning("Cloning");
57
58 if ( fParticleToBias == 0 )
59 {
61 ed << " Particle `" << particleName << "' not found !" << G4endl;
62 G4Exception(" G4BOptrForceCollision::G4BOptrForceCollision(...)",
63 "BIAS.GEN.07",
65 ed);
66 }
67}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4endl
Definition: G4ios.hh:57
G4BOptnCloning * fCloningOperation
G4BOptrForceCollisionTrackData * fCurrentTrackData
G4BOptnForceCommonTruncatedExp * fSharedForceInteractionOperation
const G4ParticleDefinition * fParticleToBias
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4int GetModelID(const G4int modelIndex)
G4VBiasingOperator(G4String name)
const char * name(G4int ptype)

References fCloningOperation, G4ParticleTable::FindParticle(), fParticleToBias, fSharedForceInteractionOperation, G4endl, G4Exception(), G4ParticleTable::GetParticleTable(), and JustWarning.

◆ G4BOptrForceCollision() [2/2]

G4BOptrForceCollision::G4BOptrForceCollision ( const G4ParticleDefinition particleToForce,
G4String  name = "ForceCollision" 
)

Definition at line 70 of file G4BOptrForceCollision.cc.

72 fForceCollisionModelID(G4PhysicsModelCatalog::GetModelID("model_GenBiasForceCollision")),
73 fCurrentTrack(nullptr),
74 fCurrentTrackData(nullptr),
76 fSetup(true)
77{
79 fCloningOperation = new G4BOptnCloning("Cloning");
80 fParticleToBias = particle;
81}

References fCloningOperation, fParticleToBias, and fSharedForceInteractionOperation.

◆ ~G4BOptrForceCollision()

G4BOptrForceCollision::~G4BOptrForceCollision ( )

Definition at line 84 of file G4BOptrForceCollision.cc.

85{
86 for ( std::map< const G4BiasingProcessInterface*, G4BOptnForceFreeFlight* >::iterator it = fFreeFlightOperations.begin() ;
87 it != fFreeFlightOperations.end() ;
88 it++ ) delete (*it).second;
90 delete fCloningOperation;
91}
std::map< const G4BiasingProcessInterface *, G4BOptnForceFreeFlight * > fFreeFlightOperations

References fCloningOperation, fFreeFlightOperations, and fSharedForceInteractionOperation.

Member Function Documentation

◆ AttachTo()

void G4VBiasingOperator::AttachTo ( const G4LogicalVolume logical)
inherited

Definition at line 59 of file G4VBiasingOperator.cc.

60{
62 it = fLogicalToSetupMap.Find(logical);
63 if ( it == fLogicalToSetupMap.End() ) fLogicalToSetupMap[logical] = this;
64 else if ( (*it).second != this )
65 {
67 ed << "Biasing operator `" << GetName()
68 << "' can not be attached to Logical volume `"
69 << logical->GetName() << "' which is already used by another operator !" << G4endl;
70 G4Exception("G4VBiasingOperator::AttachTo(...)",
71 "BIAS.MNG.01",
73 ed);
74 }
75}
const G4String & GetName() const
iterator Find(const key_type &k)
Definition: G4Cache.hh:459
iterator End()
Definition: G4Cache.hh:453
static G4MapCache< const G4LogicalVolume *, G4VBiasingOperator * > fLogicalToSetupMap
const G4String GetName() const

References G4MapCache< KEYTYPE, VALTYPE >::End(), G4MapCache< KEYTYPE, VALTYPE >::Find(), G4VBiasingOperator::fLogicalToSetupMap, G4endl, G4Exception(), G4LogicalVolume::GetName(), G4VBiasingOperator::GetName(), and JustWarning.

◆ Configure()

void G4BOptrForceCollision::Configure ( )
finalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 94 of file G4BOptrForceCollision.cc.

95{
96 // -- build free flight operations:
98}
virtual void ConfigureForWorker() final

References ConfigureForWorker().

◆ ConfigureForWorker()

void G4BOptrForceCollision::ConfigureForWorker ( )
finalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 101 of file G4BOptrForceCollision.cc.

102{
103 // -- start by remembering processes under biasing, and create needed biasing operations:
104 if ( fSetup )
105 {
106 const G4ProcessManager* processManager = fParticleToBias->GetProcessManager();
107 const G4BiasingProcessSharedData* interfaceProcessSharedData = G4BiasingProcessInterface::GetSharedData( processManager );
108 if ( interfaceProcessSharedData ) // -- sharedData tested, as is can happen a user attaches an operator
109 // -- to a volume but without defining BiasingProcessInterface processes.
110 {
111 for ( size_t i = 0 ; i < (interfaceProcessSharedData->GetPhysicsBiasingProcessInterfaces()).size(); i++ )
112 {
113 const G4BiasingProcessInterface* wrapperProcess =
114 (interfaceProcessSharedData->GetPhysicsBiasingProcessInterfaces())[i];
115 G4String operationName = "FreeFlight-"+wrapperProcess->GetWrappedProcess()->GetProcessName();
116 fFreeFlightOperations[wrapperProcess] = new G4BOptnForceFreeFlight(operationName);
117 }
118 }
119 fSetup = false;
120 }
121}
const G4BiasingProcessSharedData * GetSharedData() const
G4VProcess * GetWrappedProcess() const
const std::vector< const G4BiasingProcessInterface * > & GetPhysicsBiasingProcessInterfaces() const
G4ProcessManager * GetProcessManager() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

References fFreeFlightOperations, fParticleToBias, fSetup, G4BiasingProcessSharedData::GetPhysicsBiasingProcessInterfaces(), G4ParticleDefinition::GetProcessManager(), G4VProcess::GetProcessName(), G4BiasingProcessInterface::GetSharedData(), and G4BiasingProcessInterface::GetWrappedProcess().

Referenced by Configure().

◆ EndTracking()

void G4BOptrForceCollision::EndTracking ( )
finalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 310 of file G4BOptrForceCollision.cc.

311{
312 // -- check for consistency, operator should have cleaned the track:
313 if ( fCurrentTrackData != nullptr )
314 {
316 {
318 {
320 ed << "Current track deleted while under biasing by " << GetName() << ". Will result in inconsistencies.";
321 G4Exception(" G4BOptrForceCollision::EndTracking()",
322 "BIAS.GEN.18",
324 ed);
325 }
326 }
327 }
328}
@ fKillTrackAndSecondaries
@ fStopAndKill
G4TrackStatus GetTrackStatus() const

References fCurrentTrack, fCurrentTrackData, fKillTrackAndSecondaries, fStopAndKill, G4Exception(), G4VBiasingOperator::GetName(), G4Track::GetTrackStatus(), G4BOptrForceCollisionTrackData::IsFreeFromBiasing(), and JustWarning.

◆ ExitBiasing()

virtual void G4BOptrForceCollision::ExitBiasing ( const G4Track ,
const G4BiasingProcessInterface  
)
inlinefinalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 75 of file G4BOptrForceCollision.hh.

75{};

◆ ExitingBiasing()

void G4VBiasingOperator::ExitingBiasing ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
inherited

Definition at line 152 of file G4VBiasingOperator.cc.

153{
154 ExitBiasing( track, callingProcess );
155
156 // -- reset all data members:
167}
virtual void ExitBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
const G4VBiasingOperation * fPreviousProposedNonPhysicsBiasingOperation
G4VBiasingOperation * fNonPhysicsBiasingOperation
G4BiasingAppliedCase fPreviousBiasingAppliedCase
G4VBiasingOperation * fOccurenceBiasingOperation
G4VBiasingOperation * fFinalStateBiasingOperation
const G4VBiasingOperation * fPreviousProposedOccurenceBiasingOperation
const G4VBiasingOperation * fPreviousAppliedOccurenceBiasingOperation
const G4VBiasingOperation * fPreviousAppliedNonPhysicsBiasingOperation
const G4VBiasingOperation * fPreviousProposedFinalStateBiasingOperation
const G4VBiasingOperation * fPreviousAppliedFinalStateBiasingOperation

References BAC_None, G4VBiasingOperator::ExitBiasing(), G4VBiasingOperator::fFinalStateBiasingOperation, G4VBiasingOperator::fNonPhysicsBiasingOperation, G4VBiasingOperator::fOccurenceBiasingOperation, G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation, G4VBiasingOperator::fPreviousAppliedNonPhysicsBiasingOperation, G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation, G4VBiasingOperator::fPreviousBiasingAppliedCase, G4VBiasingOperator::fPreviousProposedFinalStateBiasingOperation, G4VBiasingOperator::fPreviousProposedNonPhysicsBiasingOperation, and G4VBiasingOperator::fPreviousProposedOccurenceBiasingOperation.

◆ GetBiasingOperator()

G4VBiasingOperator * G4VBiasingOperator::GetBiasingOperator ( const G4LogicalVolume logical)
staticinherited

◆ GetBiasingOperators()

static const std::vector< G4VBiasingOperator * > & G4VBiasingOperator::GetBiasingOperators ( )
inlinestaticinherited

◆ GetName()

const G4String G4VBiasingOperator::GetName ( ) const
inlineinherited

◆ GetPreviousBiasingAppliedCase()

G4BiasingAppliedCase G4VBiasingOperator::GetPreviousBiasingAppliedCase ( ) const
inlineinherited

◆ GetPreviousNonPhysicsAppliedOperation()

const G4VBiasingOperation * G4VBiasingOperator::GetPreviousNonPhysicsAppliedOperation ( )
inlineinherited

◆ GetProposedFinalStateBiasingOperation()

G4VBiasingOperation * G4VBiasingOperator::GetProposedFinalStateBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
inherited

Definition at line 92 of file G4VBiasingOperator.cc.

93{
96}
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

References G4VBiasingOperator::fFinalStateBiasingOperation, and G4VBiasingOperator::ProposeFinalStateBiasingOperation().

◆ GetProposedNonPhysicsBiasingOperation()

G4VBiasingOperation * G4VBiasingOperator::GetProposedNonPhysicsBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
inherited

Definition at line 98 of file G4VBiasingOperator.cc.

99{
102}
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

References G4VBiasingOperator::fNonPhysicsBiasingOperation, and G4VBiasingOperator::ProposeNonPhysicsBiasingOperation().

◆ GetProposedOccurenceBiasingOperation()

G4VBiasingOperation * G4VBiasingOperator::GetProposedOccurenceBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
inherited

Definition at line 86 of file G4VBiasingOperator.cc.

87{
90}
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

References G4VBiasingOperator::fOccurenceBiasingOperation, and G4VBiasingOperator::ProposeOccurenceBiasingOperation().

Referenced by G4ChannelingOptrMultiParticleChangeCrossSection::ProposeOccurenceBiasingOperation().

◆ OperationApplied() [1/2]

void G4BOptrForceCollision::OperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation occurenceOperationApplied,
G4double  weightForOccurenceInteraction,
G4VBiasingOperation finalStateOperationApplied,
const G4VParticleChange particleChangeProduced 
)
finalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 401 of file G4BOptrForceCollision.cc.

404{
405
407 {
408 if ( finalStateOperationApplied != fSharedForceInteractionOperation )
409 {
411 ed << " Internal inconsistency : please submit bug report. " << G4endl;
412 G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
413 "BIAS.GEN.20.5",
415 ed);
416 }
417 if ( fSharedForceInteractionOperation->GetInteractionOccured() ) fCurrentTrackData->Reset(); // -- off biasing for this track
418 }
419 else
420 {
422 ed << " Internal inconsistency : please submit bug report. " << G4endl;
423 G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
424 "BIAS.GEN.20.6",
426 ed);
427 }
428
429}

References fCurrentTrackData, G4BOptrForceCollisionTrackData::fForceCollisionState, fSharedForceInteractionOperation, G4endl, G4Exception(), G4BOptnForceCommonTruncatedExp::GetInteractionOccured(), JustWarning, G4BOptrForceCollisionTrackData::Reset(), and toBeForced.

◆ OperationApplied() [2/2]

void G4BOptrForceCollision::OperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation operationApplied,
const G4VParticleChange particleChangeProduced 
)
finalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 331 of file G4BOptrForceCollision.cc.

335{
336
337 if ( fCurrentTrackData == nullptr )
338 {
339 if ( BAC != BAC_None )
340 {
342 ed << " Internal inconsistency : please submit bug report. " << G4endl;
343 G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
344 "BIAS.GEN.20.1",
346 ed);
347 }
348 return;
349 }
350
352 {
354 auto cloneData = new G4BOptrForceCollisionTrackData( this );
355 cloneData->fForceCollisionState = ForceCollisionState::toBeForced;
357 }
359 {
360 if ( fFreeFlightOperations[callingProcess]->OperationComplete() ) fCurrentTrackData->Reset(); // -- off biasing for this track
361 }
363 {
364 if ( operationApplied != fSharedForceInteractionOperation )
365 {
367 ed << " Internal inconsistency : please submit bug report. " << G4endl;
368 G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
369 "BIAS.GEN.20.2",
371 ed);
372 }
374 {
375 if ( operationApplied != fSharedForceInteractionOperation )
376 {
378 ed << " Internal inconsistency : please submit bug report. " << G4endl;
379 G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
380 "BIAS.GEN.20.3",
382 ed);
383 }
384 }
385 }
386 else
387 {
389 {
391 ed << " Internal inconsistency : please submit bug report. " << G4endl;
392 G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
393 "BIAS.GEN.20.4",
395 ed);
396 }
397 }
398}
G4Track * GetCloneTrack() const
void SetAuxiliaryTrackInformation(G4int id, G4VAuxiliaryTrackInformation *info) const
Definition: G4Track.cc:205

References BAC_None, fCloningOperation, fCurrentTrackData, fForceCollisionModelID, G4BOptrForceCollisionTrackData::fForceCollisionState, fFreeFlightOperations, free, fSharedForceInteractionOperation, G4endl, G4Exception(), G4BOptnCloning::GetCloneTrack(), G4BOptnForceCommonTruncatedExp::GetInteractionOccured(), JustWarning, G4BOptrForceCollisionTrackData::Reset(), G4Track::SetAuxiliaryTrackInformation(), toBeCloned, toBeForced, and toBeFreeFlight.

◆ ProposeFinalStateBiasingOperation()

G4VBiasingOperation * G4BOptrForceCollision::ProposeFinalStateBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
finalprivatevirtual

Implements G4VBiasingOperator.

Definition at line 294 of file G4BOptrForceCollision.cc.

295{
296 // -- Propose at final state generation the same operation which was proposed at GPIL level,
297 // -- (which is either the force free flight or the force interaction operation).
298 // -- count on the process interface to collect this operation:
299 return callingProcess->GetCurrentOccurenceBiasingOperation();
300}
G4VBiasingOperation * GetCurrentOccurenceBiasingOperation() const

References G4BiasingProcessInterface::GetCurrentOccurenceBiasingOperation().

◆ ProposeNonPhysicsBiasingOperation()

G4VBiasingOperation * G4BOptrForceCollision::ProposeNonPhysicsBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
finalprivatevirtual

Implements G4VBiasingOperator.

Definition at line 251 of file G4BOptrForceCollision.cc.

253{
254 if ( track->GetDefinition() != fParticleToBias ) return nullptr;
255
256 // -- Apply biasing scheme only to tracks entering the volume.
257 // -- A "cloning" is done:
258 // -- - the primary will be forced flight under a zero weight up to volume exit,
259 // -- where the weight will be restored with proper weight for free flight
260 // -- - the clone will be forced to interact in the volume.
261 if ( track->GetStep()->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ) // -- §§§ extent to case of a track shoot on the boundary ?
262 {
263 // -- check that track is free of undergoing biasing scheme ( no biasing data, or no active active )
264 // -- Get possible track data:
266 if ( fCurrentTrackData != nullptr )
267 {
269 {
270 // -- takes "ownership" (some track data created before, left free, reuse it):
272 }
273 else
274 {
275 // §§§ Would something be really wrong in this case ? Could this be that a process made a zero step ?
276 }
277 }
278 else
279 {
282 }
286 return fCloningOperation;
287 }
288
289 // --
290 return nullptr;
291}
@ fGeomBoundary
Definition: G4StepStatus.hh:43
void SetCloneWeights(G4double clone1Weight, G4double clone2Weight)
const G4BOptrForceCollision * fForceCollisionOperator
G4StepStatus GetStepStatus() const
G4StepPoint * GetPreStepPoint() const
G4double GetWeight() const
G4VAuxiliaryTrackInformation * GetAuxiliaryTrackInformation(G4int id) const
Definition: G4Track.cc:225
G4ParticleDefinition * GetDefinition() const
const G4Step * GetStep() const

References fCloningOperation, fCurrentTrackData, fForceCollisionModelID, G4BOptrForceCollisionTrackData::fForceCollisionOperator, G4BOptrForceCollisionTrackData::fForceCollisionState, fGeomBoundary, fInitialTrackWeight, fParticleToBias, G4Track::GetAuxiliaryTrackInformation(), G4Track::GetDefinition(), G4Step::GetPreStepPoint(), G4Track::GetStep(), G4StepPoint::GetStepStatus(), G4Track::GetWeight(), G4BOptrForceCollisionTrackData::IsFreeFromBiasing(), G4Track::SetAuxiliaryTrackInformation(), G4BOptnCloning::SetCloneWeights(), and toBeCloned.

◆ ProposeOccurenceBiasingOperation()

G4VBiasingOperation * G4BOptrForceCollision::ProposeOccurenceBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
finalprivatevirtual

Implements G4VBiasingOperator.

Definition at line 129 of file G4BOptrForceCollision.cc.

130{
131 // -- does nothing if particle is not of requested type:
132 if ( track->GetDefinition() != fParticleToBias ) return 0;
133
134 // -- trying to get auxiliary track data...
135 if ( fCurrentTrackData == nullptr )
136 {
137 // ... and if the track has no aux. track data, it means its biasing is not started yet (note that cloning is to happen first):
139 if ( fCurrentTrackData == nullptr ) return nullptr;
140 }
141
142
143 // -- Send force free flight to the callingProcess:
144 // ------------------------------------------------
145 // -- The track has been cloned in the previous step, it has now to be
146 // -- forced for a free flight.
147 // -- This track will fly with 0.0 weight during its forced flight:
148 // -- it is to forbid the double counting with the force interaction track.
149 // -- Its weight is restored at the end of its free flight, this weight
150 // -- being its initial weight * the weight for the free flight travel,
151 // -- this last one being per process. The initial weight is common, and is
152 // -- arbitrary asked to the first operation to take care of it.
154 {
155 G4BOptnForceFreeFlight* operation = fFreeFlightOperations[callingProcess];
156 if ( callingProcess->GetWrappedProcess()->GetCurrentInteractionLength() < DBL_MAX/10. )
157 {
158 // -- the initial track weight will be restored only by the first DoIt free flight:
160 return operation;
161 }
162 else
163 {
164 return nullptr;
165 }
166 }
167
168
169 // -- Send force interaction operation to the callingProcess:
170 // ----------------------------------------------------------
171 // -- at this level, a copy of the track entering the volume was
172 // -- generated (borned) earlier. This copy will make the forced
173 // -- interaction in the volume.
175 {
176 // -- Remember if this calling process is the first of the physics wrapper in
177 // -- the PostStepGPIL loop (using default argument of method below):
178 G4bool isFirstPhysGPIL = callingProcess-> GetIsFirstPostStepGPILInterface();
179
180 // -- [*first process*] Initialize or update force interaction operation:
181 if ( isFirstPhysGPIL )
182 {
183 // -- first step of cloned track, initialize the forced interaction operation:
185 else
186 {
188 {
189 // -- means that some other physics process, not under control of the forced interaction operation,
190 // -- has occured, need to re-initialize the operation as distance to boundary has changed.
191 // -- [ Note the re-initialization is only possible for a Markovian law. ]
193 }
194 else
195 {
196 // -- means that some other non-physics process (biasing or not, like step limit), has occured,
197 // -- but track conserves its momentum direction, only need to reduced the maximum distance for
198 // -- forced interaction.
199 // -- [ Note the update is only possible for a Markovian law. ]
201 }
202 }
203 }
204
205 // -- [*all processes*] Sanity check : it may happen in limit cases that distance to
206 // -- out is zero, weight would be infinite in this case: abort forced interaction
207 // -- and abandon biasing.
209 {
211 return 0;
212 }
213
214 // -- [* first process*] collect cross-sections and sample force law to determine interaction length
215 // -- and winning process:
216 if ( isFirstPhysGPIL )
217 {
218 // -- collect cross-sections:
219 // -- ( Remember that the first of the G4BiasingProcessInterface triggers the update
220 // -- of these cross-sections )
221 const G4BiasingProcessSharedData* sharedData = callingProcess->GetSharedData();
222 for ( size_t i = 0 ; i < (sharedData->GetPhysicsBiasingProcessInterfaces()).size() ; i++ )
223 {
224 const G4BiasingProcessInterface* wrapperProcess = ( sharedData->GetPhysicsBiasingProcessInterfaces() )[i];
225 G4double interactionLength = wrapperProcess->GetWrappedProcess()->GetCurrentInteractionLength();
226 // -- keep only well defined cross-sections, other processes are ignored. These are not pathological cases
227 // -- but cases where a threhold effect par example (eg pair creation) may be at play. (**!**)
228 if ( interactionLength < DBL_MAX/10. )
229 fSharedForceInteractionOperation->AddCrossSection( wrapperProcess->GetWrappedProcess(), 1.0/interactionLength );
230 }
231 // -- sample the shared law (interaction length, and winning process):
233 }
234
235 // -- [*all processes*] Send operation for processes with well defined XS (see "**!**" above):
236 G4VBiasingOperation* operationToReturn = nullptr;
237 if ( callingProcess->GetWrappedProcess()->GetCurrentInteractionLength() < DBL_MAX/10. ) operationToReturn = fSharedForceInteractionOperation;
238 return operationToReturn;
239
240
241 } // -- end of "if ( fCurrentTrackData->fForceCollisionState == ForceCollisionState::toBeForced )"
242
243
244 // -- other cases here: particle appearing in the volume by some
245 // -- previous interaction : we decide to not bias these.
246 return 0;
247
248}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
void AddCrossSection(const G4VProcess *, G4double)
const G4ThreeVector & GetInitialMomentum() const
void ResetInitialTrackWeight(G4double w)
G4int GetCurrentStepNumber() const
G4ThreeVector GetMomentum() const
G4double GetCurrentInteractionLength() const
Definition: G4VProcess.hh:443
#define DBL_MIN
Definition: templates.hh:54
#define DBL_MAX
Definition: templates.hh:62

References G4BOptnForceCommonTruncatedExp::AddCrossSection(), DBL_MAX, DBL_MIN, fCurrentTrackData, fForceCollisionModelID, G4BOptrForceCollisionTrackData::fForceCollisionState, fFreeFlightOperations, fInitialTrackWeight, fParticleToBias, fSharedForceInteractionOperation, G4Track::GetAuxiliaryTrackInformation(), G4VProcess::GetCurrentInteractionLength(), G4Track::GetCurrentStepNumber(), G4Track::GetDefinition(), G4BOptnForceCommonTruncatedExp::GetInitialMomentum(), G4BOptnForceCommonTruncatedExp::GetMaximumDistance(), G4Track::GetMomentum(), G4BOptnForceCommonTruncatedExp::GetNumberOfSharing(), G4BiasingProcessSharedData::GetPhysicsBiasingProcessInterfaces(), G4BiasingProcessInterface::GetSharedData(), G4Track::GetStep(), G4BiasingProcessInterface::GetWrappedProcess(), G4BOptnForceCommonTruncatedExp::Initialize(), G4BOptrForceCollisionTrackData::Reset(), G4BOptnForceFreeFlight::ResetInitialTrackWeight(), G4BOptnForceCommonTruncatedExp::Sample(), toBeForced, toBeFreeFlight, and G4BOptnForceCommonTruncatedExp::UpdateForStep().

◆ ReportOperationApplied() [1/2]

void G4VBiasingOperator::ReportOperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation occurenceOperationApplied,
G4double  weightForOccurenceInteraction,
G4VBiasingOperation finalStateOperationApplied,
const G4VParticleChange particleChangeProduced 
)
inherited

Definition at line 138 of file G4VBiasingOperator.cc.

144{
145 fPreviousBiasingAppliedCase = biasingCase;
146 fPreviousAppliedOccurenceBiasingOperation = occurenceOperationApplied;
147 fPreviousAppliedFinalStateBiasingOperation = finalStateOperationApplied;
148 OperationApplied( callingProcess, biasingCase, occurenceOperationApplied, weightForOccurenceInteraction, finalStateOperationApplied, particleChangeProduced );
149}
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)

References G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation, G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation, G4VBiasingOperator::fPreviousBiasingAppliedCase, and G4VBiasingOperator::OperationApplied().

◆ ReportOperationApplied() [2/2]

void G4VBiasingOperator::ReportOperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation operationApplied,
const G4VParticleChange particleChangeProduced 
)
inherited

Definition at line 104 of file G4VBiasingOperator.cc.

108{
109 fPreviousBiasingAppliedCase = biasingCase;
113 switch ( biasingCase )
114 {
115 case BAC_None:
116 break;
117 case BAC_NonPhysics:
119 break;
120 case BAC_FinalState:
122 break;
123 case BAC_Occurence:
124 G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
125 "BIAS.MNG.02",
127 "Internal logic error, please report !");
128 break;
129 default:
130 G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
131 "BIAS.MNG.03",
133 "Internal logic error, please report !");
134 }
135 OperationApplied( callingProcess, biasingCase, operationApplied, particleChangeProduced );
136}
@ BAC_Occurence
@ BAC_FinalState
@ BAC_NonPhysics

References BAC_FinalState, BAC_None, BAC_NonPhysics, BAC_Occurence, G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation, G4VBiasingOperator::fPreviousAppliedNonPhysicsBiasingOperation, G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation, G4VBiasingOperator::fPreviousBiasingAppliedCase, G4Exception(), JustWarning, and G4VBiasingOperator::OperationApplied().

Referenced by G4ChannelingOptrMultiParticleChangeCrossSection::OperationApplied().

◆ StartRun()

void G4BOptrForceCollision::StartRun ( )
finalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 124 of file G4BOptrForceCollision.cc.

125{
126}

◆ StartTracking()

void G4BOptrForceCollision::StartTracking ( const G4Track track)
finalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 303 of file G4BOptrForceCollision.cc.

304{
305 fCurrentTrack = track;
306 fCurrentTrackData = nullptr;
307}

References fCurrentTrack, and fCurrentTrackData.

Field Documentation

◆ fCloningOperation

G4BOptnCloning* G4BOptrForceCollision::fCloningOperation
private

◆ fCurrentTrack

const G4Track* G4BOptrForceCollision::fCurrentTrack
private

Definition at line 88 of file G4BOptrForceCollision.hh.

Referenced by EndTracking(), and StartTracking().

◆ fCurrentTrackData

G4BOptrForceCollisionTrackData* G4BOptrForceCollision::fCurrentTrackData
private

◆ fDepthInTree

std::map< const G4LogicalVolume*, G4int > G4VBiasingOperator::fDepthInTree
privateinherited

Definition at line 333 of file G4VBiasingOperator.hh.

◆ fFinalStateBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fFinalStateBiasingOperation
privateinherited

◆ fForceCollisionModelID

G4int G4BOptrForceCollision::fForceCollisionModelID
private

◆ fFreeFlightOperations

std::map< const G4BiasingProcessInterface*, G4BOptnForceFreeFlight* > G4BOptrForceCollision::fFreeFlightOperations
private

◆ fInitialTrackWeight

G4double G4BOptrForceCollision::fInitialTrackWeight
private

◆ fLogicalToSetupMap

G4MapCache< const G4LogicalVolume *, G4VBiasingOperator * > G4VBiasingOperator::fLogicalToSetupMap
staticprivateinherited

◆ fName

const G4String G4VBiasingOperator::fName
privateinherited

Definition at line 319 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::GetName().

◆ fNonPhysicsBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fNonPhysicsBiasingOperation
privateinherited

◆ fOccurenceBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fOccurenceBiasingOperation
privateinherited

◆ fOperators

G4VectorCache< G4VBiasingOperator * > G4VBiasingOperator::fOperators
staticprivateinherited

◆ fParticleToBias

const G4ParticleDefinition* G4BOptrForceCollision::fParticleToBias
private

◆ fPreviousAppliedFinalStateBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation
privateinherited

◆ fPreviousAppliedNonPhysicsBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedNonPhysicsBiasingOperation
privateinherited

◆ fPreviousAppliedOccurenceBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation
privateinherited

◆ fPreviousBiasingAppliedCase

G4BiasingAppliedCase G4VBiasingOperator::fPreviousBiasingAppliedCase
privateinherited

◆ fPreviousProposedFinalStateBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedFinalStateBiasingOperation
privateinherited

Definition at line 342 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::ExitingBiasing().

◆ fPreviousProposedNonPhysicsBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedNonPhysicsBiasingOperation
privateinherited

Definition at line 343 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::ExitingBiasing().

◆ fPreviousProposedOccurenceBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedOccurenceBiasingOperation
privateinherited

Definition at line 341 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::ExitingBiasing().

◆ fRootVolumes

std::vector< const G4LogicalVolume* > G4VBiasingOperator::fRootVolumes
privateinherited

Definition at line 332 of file G4VBiasingOperator.hh.

◆ fSetup

G4bool G4BOptrForceCollision::fSetup
private

Definition at line 94 of file G4BOptrForceCollision.hh.

Referenced by ConfigureForWorker().

◆ fSharedForceInteractionOperation

G4BOptnForceCommonTruncatedExp* G4BOptrForceCollision::fSharedForceInteractionOperation
private

◆ fStateNotifier

G4Cache< G4BiasingOperatorStateNotifier * > G4VBiasingOperator::fStateNotifier
staticprivateinherited

Definition at line 328 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::G4VBiasingOperator().


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