Loading...
Searching...
No Matches
GB02BOptrMultiParticleForceCollision.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26/// \file GB02/src/GB02BOptrMultiParticleForceCollision.cc
27/// \brief Implementation of the GB02BOptrMultiParticleForceCollision class
28//
30#include "G4BiasingProcessInterface.hh"
31
32#include "G4BOptrForceCollision.hh"
33#include "G4ParticleDefinition.hh"
34#include "G4ParticleTable.hh"
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
37
42
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44
46{
47 const G4ParticleDefinition* particle =
48 G4ParticleTable::GetParticleTable()->FindParticle( particleName );
49
50 if ( particle == 0 )
51 {
52 G4ExceptionDescription ed;
53 ed << "Particle `" << particleName << "' not found !" << G4endl;
54 G4Exception("GB02BOptrMultiParticleForceCollision::AddParticle(...)",
55 "exGB02.01",
56 JustWarning,
57 ed);
58 return;
59 }
60
61 G4BOptrForceCollision* optr = new G4BOptrForceCollision(particleName,
62 "ForceCollisionFor"+particleName);
63 fParticlesToBias.push_back( particle );
64 fBOptrForParticle[ particle ] = optr;
65
66}
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69
73 const G4BiasingProcessInterface* callingProcess)
74{
76 GetProposedOccurenceBiasingOperation(track, callingProcess);
77 else return 0;
78}
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
81
85 const G4BiasingProcessInterface* callingProcess)
86{
88 GetProposedNonPhysicsBiasingOperation(track, callingProcess);
89 else return 0;
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93
97 const G4BiasingProcessInterface* callingProcess)
98{
100 GetProposedFinalStateBiasingOperation(track, callingProcess);
101 else return 0;
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105
107{
108 const G4ParticleDefinition* definition = track->GetParticleDefinition();
109 std::map < const G4ParticleDefinition*, G4BOptrForceCollision* > :: iterator
110 it = fBOptrForParticle.find( definition );
112 if ( it != fBOptrForParticle.end() ) fCurrentOperator = (*it).second;
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116
117void
119OperationApplied( const G4BiasingProcessInterface* callingProcess,
120 G4BiasingAppliedCase biasingCase,
121 G4VBiasingOperation* operationApplied,
122 const G4VParticleChange* particleChangeProduced )
123{
124 if ( fCurrentOperator ) fCurrentOperator->ReportOperationApplied( callingProcess,
125 biasingCase,
126 operationApplied,
127 particleChangeProduced );
128
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
132
133void
135OperationApplied( const G4BiasingProcessInterface* callingProcess,
136 G4BiasingAppliedCase biasingCase,
137 G4VBiasingOperation* occurenceOperationApplied,
138 G4double weightForOccurenceInteraction,
139 G4VBiasingOperation* finalStateOperationApplied,
140 const G4VParticleChange* particleChangeProduced )
141{
142 if ( fCurrentOperator ) fCurrentOperator->ReportOperationApplied( callingProcess,
143 biasingCase,
144 occurenceOperationApplied,
145 weightForOccurenceInteraction,
146 finalStateOperationApplied,
147 particleChangeProduced );
148}
149
150//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
151
152void
154ExitBiasing( const G4Track* track,
155 const G4BiasingProcessInterface* callingProcess )
156{
157 if ( fCurrentOperator ) fCurrentOperator->ExitingBiasing( track, callingProcess );
158}
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Definition of the GB02BOptrMultiParticleForceCollision class.
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess) final
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess) final
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess) final
virtual void StartTracking(const G4Track *track) final
void ExitBiasing(const G4Track *, const G4BiasingProcessInterface *) final
std::map< const G4ParticleDefinition *, G4BOptrForceCollision * > fBOptrForParticle
void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced) final
std::vector< const G4ParticleDefinition * > fParticlesToBias

Applications | User Support | Publications | Collaboration