Geant4 11.2.2
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 ()
 
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 *)
 
- Protected Member Functions inherited from G4VBiasingOperator

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(G4PhysicsModelCatalog::GetModelID("model_GenBiasForceCollision")),
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)
std::ostringstream G4ExceptionDescription
#define G4endl
Definition G4ios.hh:67
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4int GetModelID(const G4int modelIndex)
G4VBiasingOperator(G4String name)

◆ G4BOptrForceCollision() [2/2]

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

Definition at line 70 of file G4BOptrForceCollision.cc.

71 : G4VBiasingOperator(name),
72 fForceCollisionModelID(G4PhysicsModelCatalog::GetModelID("model_GenBiasForceCollision")),
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 // -- build free flight operations:
98}
virtual void ConfigureForWorker() final

◆ 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
const std::vector< const G4BiasingProcessInterface * > & GetPhysicsBiasingProcessInterfaces() const
G4ProcessManager * GetProcessManager() const
const G4String & GetProcessName() const

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 {
315 if ( !fCurrentTrackData->IsFreeFromBiasing() )
316 {
317 if ( (fCurrentTrack->GetTrackStatus() == fStopAndKill) || (fCurrentTrack->GetTrackStatus() == fKillTrackAndSecondaries) )
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
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 401 of file G4BOptrForceCollision.cc.

404{
405
406 if ( fCurrentTrackData->fForceCollisionState == ForceCollisionState::toBeForced )
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}

◆ 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
351 if ( fCurrentTrackData->fForceCollisionState == ForceCollisionState::toBeCloned )
352 {
353 fCurrentTrackData->fForceCollisionState = ForceCollisionState::toBeFreeFlight;
354 auto cloneData = new G4BOptrForceCollisionTrackData( this );
355 cloneData->fForceCollisionState = ForceCollisionState::toBeForced;
356 fCloningOperation->GetCloneTrack()->SetAuxiliaryTrackInformation(fForceCollisionModelID, cloneData);
357 }
358 else if ( fCurrentTrackData->fForceCollisionState == ForceCollisionState::toBeFreeFlight )
359 {
360 if ( fFreeFlightOperations[callingProcess]->OperationComplete() ) fCurrentTrackData->Reset(); // -- off biasing for this track
361 }
362 else if ( fCurrentTrackData->fForceCollisionState == ForceCollisionState::toBeForced )
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 }
373 if ( fSharedForceInteractionOperation->GetInteractionOccured() )
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 {
388 if ( fCurrentTrackData->fForceCollisionState != ForceCollisionState::free )
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:209

◆ 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}

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