File indexing completed on 2025-02-23 09:20:45
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029 #include "GB02BOptrMultiParticleForceCollision.hh"
0030
0031 #include "G4BOptrForceCollision.hh"
0032 #include "G4BiasingProcessInterface.hh"
0033 #include "G4ParticleDefinition.hh"
0034 #include "G4ParticleTable.hh"
0035
0036
0037
0038 GB02BOptrMultiParticleForceCollision::GB02BOptrMultiParticleForceCollision()
0039 : G4VBiasingOperator("TestManyForceCollision")
0040 {}
0041
0042
0043
0044 void GB02BOptrMultiParticleForceCollision::AddParticle(G4String particleName)
0045 {
0046 const G4ParticleDefinition* particle =
0047 G4ParticleTable::GetParticleTable()->FindParticle(particleName);
0048
0049 if (particle == 0) {
0050 G4ExceptionDescription ed;
0051 ed << "Particle `" << particleName << "' not found !" << G4endl;
0052 G4Exception("GB02BOptrMultiParticleForceCollision::AddParticle(...)", "exGB02.01", JustWarning,
0053 ed);
0054 return;
0055 }
0056
0057 G4BOptrForceCollision* optr =
0058 new G4BOptrForceCollision(particleName, "ForceCollisionFor" + particleName);
0059 fParticlesToBias.push_back(particle);
0060 fBOptrForParticle[particle] = optr;
0061 }
0062
0063
0064
0065 G4VBiasingOperation* GB02BOptrMultiParticleForceCollision::ProposeOccurenceBiasingOperation(
0066 const G4Track* track, const G4BiasingProcessInterface* callingProcess)
0067 {
0068 if (fCurrentOperator)
0069 return fCurrentOperator->GetProposedOccurenceBiasingOperation(track, callingProcess);
0070 else
0071 return 0;
0072 }
0073
0074
0075
0076 G4VBiasingOperation* GB02BOptrMultiParticleForceCollision::ProposeNonPhysicsBiasingOperation(
0077 const G4Track* track, const G4BiasingProcessInterface* callingProcess)
0078 {
0079 if (fCurrentOperator)
0080 return fCurrentOperator->GetProposedNonPhysicsBiasingOperation(track, callingProcess);
0081 else
0082 return 0;
0083 }
0084
0085
0086
0087 G4VBiasingOperation* GB02BOptrMultiParticleForceCollision::ProposeFinalStateBiasingOperation(
0088 const G4Track* track, const G4BiasingProcessInterface* callingProcess)
0089 {
0090 if (fCurrentOperator)
0091 return fCurrentOperator->GetProposedFinalStateBiasingOperation(track, callingProcess);
0092 else
0093 return 0;
0094 }
0095
0096
0097
0098 void GB02BOptrMultiParticleForceCollision::StartTracking(const G4Track* track)
0099 {
0100 const G4ParticleDefinition* definition = track->GetParticleDefinition();
0101 std::map<const G4ParticleDefinition*, G4BOptrForceCollision*>::iterator it =
0102 fBOptrForParticle.find(definition);
0103 fCurrentOperator = 0;
0104 if (it != fBOptrForParticle.end()) fCurrentOperator = (*it).second;
0105 }
0106
0107
0108
0109 void GB02BOptrMultiParticleForceCollision::OperationApplied(
0110 const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase,
0111 G4VBiasingOperation* operationApplied, const G4VParticleChange* particleChangeProduced)
0112 {
0113 if (fCurrentOperator)
0114 fCurrentOperator->ReportOperationApplied(callingProcess, biasingCase, operationApplied,
0115 particleChangeProduced);
0116 }
0117
0118
0119
0120 void GB02BOptrMultiParticleForceCollision::OperationApplied(
0121 const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase,
0122 G4VBiasingOperation* occurenceOperationApplied, G4double weightForOccurenceInteraction,
0123 G4VBiasingOperation* finalStateOperationApplied, const G4VParticleChange* particleChangeProduced)
0124 {
0125 if (fCurrentOperator)
0126 fCurrentOperator->ReportOperationApplied(callingProcess, biasingCase, occurenceOperationApplied,
0127 weightForOccurenceInteraction,
0128 finalStateOperationApplied, particleChangeProduced);
0129 }
0130
0131
0132
0133 void GB02BOptrMultiParticleForceCollision::ExitBiasing(
0134 const G4Track* track, const G4BiasingProcessInterface* callingProcess)
0135 {
0136 if (fCurrentOperator) fCurrentOperator->ExitingBiasing(track, callingProcess);
0137 }
0138
0139