Geant4  10.02.p03
G4BOptrForceCollision Class Reference

#include <G4BOptrForceCollision.hh>

Inheritance diagram for G4BOptrForceCollision:
Collaboration 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 ()
 

Private Member Functions

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

Private Attributes

G4int fForceCollisionModelID
 
const G4Track * fCurrentTrack
 
G4BOptrForceCollisionTrackDatafCurrentTrackData
 
std::map< const G4BiasingProcessInterface *, G4BOptnForceFreeFlight *> fFreeFlightOperations
 
G4BOptnForceCommonTruncatedExpfSharedForceInteractionOperation
 
G4BOptnCloningfCloningOperation
 
G4double fInitialTrackWeight
 
G4bool fSetup
 
const G4ParticleDefinitionfParticleToBias
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VBiasingOperator
static const std::vector< G4VBiasingOperator *> & GetBiasingOperators ()
 
static G4VBiasingOperatorGetBiasingOperator (const G4LogicalVolume *)
 

Detailed Description

Definition at line 59 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),
49  fCurrentTrack(nullptr),
50  fCurrentTrackData(nullptr),
51  fInitialTrackWeight(-1.0),
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 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4BOptrForceCollisionTrackData * fCurrentTrackData
G4BOptnForceCommonTruncatedExp * fSharedForceInteractionOperation
G4BOptnCloning * fCloningOperation
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4ParticleTable * GetParticleTable()
G4VBiasingOperator(G4String name)
#define G4endl
Definition: G4ios.hh:61
const G4ParticleDefinition * fParticleToBias
Here is the call graph for this function:

◆ G4BOptrForceCollision() [2/2]

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

Definition at line 70 of file G4BOptrForceCollision.cc.

71  : G4VBiasingOperator(name),
73  fCurrentTrack(nullptr),
74  fCurrentTrackData(nullptr),
75  fInitialTrackWeight(-1.0),
76  fSetup(true)
77 {
79  fCloningOperation = new G4BOptnCloning("Cloning");
80  fParticleToBias = particle;
81 }
G4BOptrForceCollisionTrackData * fCurrentTrackData
G4BOptnForceCommonTruncatedExp * fSharedForceInteractionOperation
G4BOptnCloning * fCloningOperation
G4VBiasingOperator(G4String name)
const G4ParticleDefinition * fParticleToBias

◆ ~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 }
G4BOptnForceCommonTruncatedExp * fSharedForceInteractionOperation
G4BOptnCloning * fCloningOperation
std::map< const G4BiasingProcessInterface *, G4BOptnForceFreeFlight *> fFreeFlightOperations

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:
98 
99  // -- build free flight operations:
101 }
static G4int Register(const G4String &)
virtual void ConfigureForWorker() final
Here is the call graph for this function:

◆ 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 std::vector< const G4BiasingProcessInterface *> & GetPhysicsBiasingProcessInterfaces() const
G4ProcessManager * GetProcessManager() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
const G4BiasingProcessSharedData * GetSharedData() const
std::map< const G4BiasingProcessInterface *, G4BOptnForceFreeFlight *> fFreeFlightOperations
static G4int Register(const G4String &)
G4VProcess * GetWrappedProcess() const
const G4ParticleDefinition * fParticleToBias
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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  {
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",
329  JustWarning,
330  ed);
331  }
332  }
333  }
334 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4BOptrForceCollisionTrackData * fCurrentTrackData
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4String GetName() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ExitBiasing()

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

Reimplemented from G4VBiasingOperator.

Definition at line 76 of file G4BOptrForceCollision.hh.

76 {};
Here is the call graph for this function:

◆ OperationApplied() [1/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",
351  JustWarning,
352  ed);
353  }
354  return;
355  }
356 
358  {
360  auto cloneData = new G4BOptrForceCollisionTrackData( this );
361  cloneData->fForceCollisionState = ForceCollisionState::toBeForced;
362  fCloningOperation->GetCloneTrack()->SetAuxiliaryTrackInformation(fForceCollisionModelID, cloneData);
363  }
365  {
366  if ( fFreeFlightOperations[callingProcess]->OperationComplete() ) fCurrentTrackData->Reset(); // -- off biasing for this track
367  }
369  {
370  if ( operationApplied != fSharedForceInteractionOperation )
371  {
373  ed << " Internal inconsistency : please submit bug report. " << G4endl;
374  G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
375  "BIAS.GEN.20.2",
376  JustWarning,
377  ed);
378  }
380  {
381  if ( operationApplied != fSharedForceInteractionOperation )
382  {
384  ed << " Internal inconsistency : please submit bug report. " << G4endl;
385  G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
386  "BIAS.GEN.20.3",
387  JustWarning,
388  ed);
389  }
390  }
391  }
392  else
393  {
395  {
397  ed << " Internal inconsistency : please submit bug report. " << G4endl;
398  G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
399  "BIAS.GEN.20.4",
400  JustWarning,
401  ed);
402  }
403  }
404 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4BOptrForceCollisionTrackData * fCurrentTrackData
G4BOptnForceCommonTruncatedExp * fSharedForceInteractionOperation
G4BOptnCloning * fCloningOperation
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::map< const G4BiasingProcessInterface *, G4BOptnForceFreeFlight *> fFreeFlightOperations
#define G4endl
Definition: G4ios.hh:61
G4Track * GetCloneTrack() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ OperationApplied() [2/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 
413  {
414  if ( finalStateOperationApplied != fSharedForceInteractionOperation )
415  {
417  ed << " Internal inconsistency : please submit bug report. " << G4endl;
418  G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
419  "BIAS.GEN.20.5",
420  JustWarning,
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",
431  JustWarning,
432  ed);
433  }
434 
435 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4BOptrForceCollisionTrackData * fCurrentTrackData
G4BOptnForceCommonTruncatedExp * fSharedForceInteractionOperation
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

◆ ProposeFinalStateBiasingOperation()

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

Implements G4VBiasingOperator.

Definition at line 300 of file G4BOptrForceCollision.cc.

301 {
302  // -- Propose at final state generation the same operation which was proposed at GPIL level,
303  // -- (which is either the force free flight or the force interaction operation).
304  // -- count on the process interface to collect this operation:
305  return callingProcess->GetCurrentOccurenceBiasingOperation();
306 }
G4VBiasingOperation * GetCurrentOccurenceBiasingOperation() const
Here is the call graph for this function:

◆ ProposeNonPhysicsBiasingOperation()

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

Implements G4VBiasingOperator.

Definition at line 257 of file G4BOptrForceCollision.cc.

259 {
260  if ( track->GetDefinition() != fParticleToBias ) return nullptr;
261 
262  // -- Apply biasing scheme only to tracks entering the volume.
263  // -- A "cloning" is done:
264  // -- - the primary will be forced flight under a zero weight up to volume exit,
265  // -- where the weight will be restored with proper weight for free flight
266  // -- - the clone will be forced to interact in the volume.
267  if ( track->GetStep()->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ) // -- ยงยงยง extent to case of a track shoot on the boundary ?
268  {
269  // -- check that track is free of undergoing biasing scheme ( no biasing data, or no active active )
270  // -- Get possible track data:
272  if ( fCurrentTrackData != nullptr )
273  {
275  {
276  // -- takes "ownership" (some track data created before, left free, reuse it):
278  }
279  else
280  {
281  // ยงยงยง Would something be really wrong in this case ? Could this be that a process made a zero step ?
282  }
283  }
284  else
285  {
287  track->SetAuxiliaryTrackInformation(fForceCollisionModelID, fCurrentTrackData);
288  }
290  fInitialTrackWeight = track->GetWeight();
292  return fCloningOperation;
293  }
294 
295  // --
296  return nullptr;
297 }
const G4BOptrForceCollision * fForceCollisionOperator
G4BOptrForceCollisionTrackData * fCurrentTrackData
G4BOptnCloning * fCloningOperation
void SetCloneWeights(G4double clone1Weight, G4double clone2Weight)
const G4ParticleDefinition * fParticleToBias
Here is the call graph for this function:

◆ ProposeOccurenceBiasingOperation()

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

Implements G4VBiasingOperator.

Definition at line 135 of file G4BOptrForceCollision.cc.

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

◆ 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 {
312  fCurrentTrackData = nullptr;
313 }
G4BOptrForceCollisionTrackData * fCurrentTrackData

Member Data Documentation

◆ fCloningOperation

G4BOptnCloning* G4BOptrForceCollision::fCloningOperation
private

Definition at line 93 of file G4BOptrForceCollision.hh.

◆ fCurrentTrack

const G4Track* G4BOptrForceCollision::fCurrentTrack
private

Definition at line 89 of file G4BOptrForceCollision.hh.

◆ fCurrentTrackData

G4BOptrForceCollisionTrackData* G4BOptrForceCollision::fCurrentTrackData
private

Definition at line 90 of file G4BOptrForceCollision.hh.

◆ fForceCollisionModelID

G4int G4BOptrForceCollision::fForceCollisionModelID
private

Definition at line 88 of file G4BOptrForceCollision.hh.

◆ fFreeFlightOperations

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

Definition at line 91 of file G4BOptrForceCollision.hh.

◆ fInitialTrackWeight

G4double G4BOptrForceCollision::fInitialTrackWeight
private

Definition at line 94 of file G4BOptrForceCollision.hh.

◆ fParticleToBias

const G4ParticleDefinition* G4BOptrForceCollision::fParticleToBias
private

Definition at line 96 of file G4BOptrForceCollision.hh.

◆ fSetup

G4bool G4BOptrForceCollision::fSetup
private

Definition at line 95 of file G4BOptrForceCollision.hh.

◆ fSharedForceInteractionOperation

G4BOptnForceCommonTruncatedExp* G4BOptrForceCollision::fSharedForceInteractionOperation
private

Definition at line 92 of file G4BOptrForceCollision.hh.


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