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

#include <G4BOptrForceCollision.hh>

+ Inheritance diagram for G4BOptrForceCollision:

Public Member Functions

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

Additional Inherited Members

- Static Public Member Functions inherited from G4VBiasingOperator
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators ()
 
static G4VBiasingOperatorGetBiasingOperator (const G4LogicalVolume *)
 
virtual G4VBiasingOperationProposeNonPhysicsBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
 
virtual G4VBiasingOperationProposeOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
 
virtual G4VBiasingOperationProposeFinalStateBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
 
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void ExitBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 

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.

47 : G4VBiasingOperator(name),
48 fForceCollisionModelID(-1),
49 fCurrentTrack(nullptr),
50 fCurrentTrackData(nullptr),
51 fInitialTrackWeight(-1.0),
52 fSetup(true)
53{
54 fSharedForceInteractionOperation = new G4BOptnForceCommonTruncatedExp("SharedForceInteraction");
55 fCloningOperation = new G4BOptnCloning("Cloning");
56 fParticleToBias = G4ParticleTable::GetParticleTable()->FindParticle(particleName);
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
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()

◆ G4BOptrForceCollision() [2/2]

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

Definition at line 70 of file G4BOptrForceCollision.cc.

71 : G4VBiasingOperator(name),
72 fForceCollisionModelID(-1),
73 fCurrentTrack(nullptr),
74 fCurrentTrackData(nullptr),
75 fInitialTrackWeight(-1.0),
76 fSetup(true)
77{
78 fSharedForceInteractionOperation = new G4BOptnForceCommonTruncatedExp("SharedForceInteraction");
79 fCloningOperation = new G4BOptnCloning("Cloning");
80 fParticleToBias = particle;
81}

◆ ~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;
89 delete fSharedForceInteractionOperation;
90 delete fCloningOperation;
91}

Member Function Documentation

◆ Configure()

void G4BOptrForceCollision::Configure ( )
finalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 94 of file G4BOptrForceCollision.cc.

95{
96 // -- Create ID for force collision:
97 fForceCollisionModelID = G4PhysicsModelCatalog::Register("GenBiasForceCollision");
98
99 // -- build free flight operations:
101}
virtual void ConfigureForWorker() final
static G4int Register(const G4String &)

◆ ConfigureForWorker()

void G4BOptrForceCollision::ConfigureForWorker ( )
finalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 104 of file G4BOptrForceCollision.cc.

105{
106 // -- start by remembering processes under biasing, and create needed biasing operations:
107 if ( fSetup )
108 {
109 // -- get back ID created in master thread in case of MT (or reget just created ID in sequential mode):
110 fForceCollisionModelID = G4PhysicsModelCatalog::Register("GenBiasForceCollision");
111
112 const G4ProcessManager* processManager = fParticleToBias->GetProcessManager();
113 const G4BiasingProcessSharedData* interfaceProcessSharedData = G4BiasingProcessInterface::GetSharedData( processManager );
114 if ( interfaceProcessSharedData ) // -- sharedData tested, as is can happen a user attaches an operator
115 // -- to a volume but without defining BiasingProcessInterface processes.
116 {
117 for ( size_t i = 0 ; i < (interfaceProcessSharedData->GetPhysicsBiasingProcessInterfaces()).size(); i++ )
118 {
119 const G4BiasingProcessInterface* wrapperProcess =
120 (interfaceProcessSharedData->GetPhysicsBiasingProcessInterfaces())[i];
121 G4String operationName = "FreeFlight-"+wrapperProcess->GetWrappedProcess()->GetProcessName();
122 fFreeFlightOperations[wrapperProcess] = new G4BOptnForceFreeFlight(operationName);
123 }
124 }
125 fSetup = false;
126 }
127}
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

Referenced by Configure().

◆ EndTracking()

void G4BOptrForceCollision::EndTracking ( )
finalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 316 of file G4BOptrForceCollision.cc.

317{
318 // -- check for consistency, operator should have cleaned the track:
319 if ( fCurrentTrackData != nullptr )
320 {
321 if ( !fCurrentTrackData->IsFreeFromBiasing() )
322 {
323 if ( (fCurrentTrack->GetTrackStatus() == fStopAndKill) || (fCurrentTrack->GetTrackStatus() == fKillTrackAndSecondaries) )
324 {
326 ed << "Current track deleted while under biasing by " << GetName() << ". Will result in inconsistencies.";
327 G4Exception(" G4BOptrForceCollision::EndTracking()",
328 "BIAS.GEN.18",
330 ed);
331 }
332 }
333 }
334}
@ fKillTrackAndSecondaries
@ fStopAndKill
G4TrackStatus GetTrackStatus() const
const G4String GetName() const

◆ ExitBiasing()

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

Reimplemented from G4VBiasingOperator.

Definition at line 75 of file G4BOptrForceCollision.hh.

75{};

◆ 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 407 of file G4BOptrForceCollision.cc.

410{
411
412 if ( fCurrentTrackData->fForceCollisionState == ForceCollisionState::toBeForced )
413 {
414 if ( finalStateOperationApplied != fSharedForceInteractionOperation )
415 {
417 ed << " Internal inconsistency : please submit bug report. " << G4endl;
418 G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
419 "BIAS.GEN.20.5",
421 ed);
422 }
423 if ( fSharedForceInteractionOperation->GetInteractionOccured() ) fCurrentTrackData->Reset(); // -- off biasing for this track
424 }
425 else
426 {
428 ed << " Internal inconsistency : please submit bug report. " << G4endl;
429 G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
430 "BIAS.GEN.20.6",
432 ed);
433 }
434
435}

◆ OperationApplied() [2/2]

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

Reimplemented from G4VBiasingOperator.

Definition at line 337 of file G4BOptrForceCollision.cc.

341{
342
343 if ( fCurrentTrackData == nullptr )
344 {
345 if ( BAC != BAC_None )
346 {
348 ed << " Internal inconsistency : please submit bug report. " << G4endl;
349 G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
350 "BIAS.GEN.20.1",
352 ed);
353 }
354 return;
355 }
356
357 if ( fCurrentTrackData->fForceCollisionState == ForceCollisionState::toBeCloned )
358 {
359 fCurrentTrackData->fForceCollisionState = ForceCollisionState::toBeFreeFlight;
360 auto cloneData = new G4BOptrForceCollisionTrackData( this );
361 cloneData->fForceCollisionState = ForceCollisionState::toBeForced;
362 fCloningOperation->GetCloneTrack()->SetAuxiliaryTrackInformation(fForceCollisionModelID, cloneData);
363 }
364 else if ( fCurrentTrackData->fForceCollisionState == ForceCollisionState::toBeFreeFlight )
365 {
366 if ( fFreeFlightOperations[callingProcess]->OperationComplete() ) fCurrentTrackData->Reset(); // -- off biasing for this track
367 }
368 else if ( fCurrentTrackData->fForceCollisionState == ForceCollisionState::toBeForced )
369 {
370 if ( operationApplied != fSharedForceInteractionOperation )
371 {
373 ed << " Internal inconsistency : please submit bug report. " << G4endl;
374 G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
375 "BIAS.GEN.20.2",
377 ed);
378 }
379 if ( fSharedForceInteractionOperation->GetInteractionOccured() )
380 {
381 if ( operationApplied != fSharedForceInteractionOperation )
382 {
384 ed << " Internal inconsistency : please submit bug report. " << G4endl;
385 G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
386 "BIAS.GEN.20.3",
388 ed);
389 }
390 }
391 }
392 else
393 {
394 if ( fCurrentTrackData->fForceCollisionState != ForceCollisionState::free )
395 {
397 ed << " Internal inconsistency : please submit bug report. " << G4endl;
398 G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
399 "BIAS.GEN.20.4",
401 ed);
402 }
403 }
404}
G4Track * GetCloneTrack() const
void SetAuxiliaryTrackInformation(G4int idx, G4VAuxiliaryTrackInformation *info) const
Definition: G4Track.cc:198

◆ StartRun()

void G4BOptrForceCollision::StartRun ( )
finalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 130 of file G4BOptrForceCollision.cc.

131{
132}

◆ StartTracking()

void G4BOptrForceCollision::StartTracking ( const G4Track track)
finalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 309 of file G4BOptrForceCollision.cc.

310{
311 fCurrentTrack = track;
312 fCurrentTrackData = nullptr;
313}

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