Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4BOptrForceCollision Class Reference

#include <G4BOptrForceCollision.hh>

Inheritance diagram for G4BOptrForceCollision:
G4VBiasingOperator

Public Member Functions

 G4BOptrForceCollision (G4String particleToForce, G4String name="ForceCollision")
 
 G4BOptrForceCollision (const G4ParticleDefinition *particleToForce, G4String name="ForceCollision")
 
 ~G4BOptrForceCollision ()
 
virtual void StartTracking (const G4Track *track)
 
virtual void ExitBiasing (const G4Track *, const G4BiasingProcessInterface *)
 
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
- Public Member Functions inherited from G4VBiasingOperator
 G4VBiasingOperator (G4String name)
 
virtual ~G4VBiasingOperator ()
 
virtual void StartRun ()
 
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 ()
 
const G4VBiasingOperationGetBirthOperation (const G4Track *)
 

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
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 
void RememberSecondaries (const G4BiasingProcessInterface *callingProcess, const G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
void ForgetTrack (const G4Track *track)
 

Detailed Description

Definition at line 58 of file G4BOptrForceCollision.hh.

Constructor & Destructor Documentation

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

Definition at line 43 of file G4BOptrForceCollision.cc.

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

44  : G4VBiasingOperator(name),
45  fFirstProcess(0), fLastProcess(0),
46  fSetup(true)
47 {
48  fSharedForceInteractionOperation = new G4BOptnForceCommonTruncatedExp("SharedForceInteraction");
49  fCloningOperation = new G4BOptnCloning("Cloning");
50  fParticle = G4ParticleTable::GetParticleTable()->FindParticle(particleName);
51 
52  if ( fParticle == 0 )
53  {
55  ed << " Particle `" << particleName << "' not found !" << G4endl;
56  G4Exception(" G4BOptrForceCollision::G4BOptrForceCollision(...)",
57  "BIAS.GEN.07",
59  ed);
60  }
61 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
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
G4BOptrForceCollision::G4BOptrForceCollision ( const G4ParticleDefinition particleToForce,
G4String  name = "ForceCollision" 
)

Definition at line 63 of file G4BOptrForceCollision.cc.

64  : G4VBiasingOperator(name),
65  fFirstProcess(0), fLastProcess(0),
66  fSetup(true)
67 {
68  fSharedForceInteractionOperation = new G4BOptnForceCommonTruncatedExp("SharedForceInteraction");
69  fCloningOperation = new G4BOptnCloning("Cloning");
70  fParticle = particle;
71 }
G4VBiasingOperator(G4String name)
G4BOptrForceCollision::~G4BOptrForceCollision ( )

Definition at line 73 of file G4BOptrForceCollision.cc.

74 {
75  for ( std::map< const G4BiasingProcessInterface*, G4BOptnForceFreeFlight* >::iterator it = fFreeFlightOperations.begin() ;
76  it != fFreeFlightOperations.end() ;
77  it++ ) delete (*it).second;
78  delete fSharedForceInteractionOperation;
79  delete fCloningOperation;
80 }

Member Function Documentation

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

Reimplemented from G4VBiasingOperator.

Definition at line 234 of file G4BOptrForceCollision.cc.

References G4VBiasingOperator::ForgetTrack().

235 {
236  fPreviousOperationApplied = 0;
237  // -- clean-up track for what happens with this operator:
238  ForgetTrack ( track );
239 }
void ForgetTrack(const G4Track *track)
void G4BOptrForceCollision::OperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation operationApplied,
const G4VParticleChange particleChangeProduced 
)
virtual

Reimplemented from G4VBiasingOperator.

Definition at line 242 of file G4BOptrForceCollision.cc.

References G4VBiasingOperator::RememberSecondaries().

244 {
245  fPreviousOperationApplied = operationApplied;
246  if ( operationApplied == fCloningOperation )
247  RememberSecondaries( callingProcess, operationApplied, particleChangeProduced );
248 }
void RememberSecondaries(const G4BiasingProcessInterface *callingProcess, const G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
void G4BOptrForceCollision::StartTracking ( const G4Track track)
virtual

Reimplemented from G4VBiasingOperator.

Definition at line 228 of file G4BOptrForceCollision.cc.

229 {
230  fPreviousOperationApplied = 0;
231 }

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