Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VBiasingOperator.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// G4VBiasingOperator
27// --------------------------------------------------------------------
28
29#include "G4VBiasingOperator.hh"
31#include "G4VParticleChange.hh"
32
33G4MapCache< const G4LogicalVolume*, G4VBiasingOperator* > G4VBiasingOperator::fLogicalToSetupMap;
34G4VectorCache< G4VBiasingOperator* > G4VBiasingOperator::fOperators;
35G4Cache< G4BiasingOperatorStateNotifier* > G4VBiasingOperator::fStateNotifier(nullptr);
36
38 : fName( name )
39{
40 fOperators.Push_back(this);
41
42 if ( fStateNotifier.Get() == nullptr )
43 fStateNotifier.Put( new G4BiasingOperatorStateNotifier() );
44}
45
47{
48 auto it = fLogicalToSetupMap.Find(logical);
49 if ( it == fLogicalToSetupMap.End() )
50 {
51 fLogicalToSetupMap[logical] = this;
52 }
53 else if ( (*it).second != this )
54 {
56 ed << "Biasing operator `" << GetName()
57 << "' can not be attached to Logical volume `"
58 << logical->GetName() << "' which is already used by another operator !"
59 << G4endl;
60 G4Exception("G4VBiasingOperator::AttachTo(...)",
61 "BIAS.MNG.01", JustWarning, ed);
62 }
63}
64
65const std::vector < G4VBiasingOperator* >&
67{
68 return fOperators.Get();
69}
70
72{
73 auto it = fLogicalToSetupMap.Find(logical);
74 if ( it == fLogicalToSetupMap.End() ) { return nullptr; }
75 else { return (*it).second; }
76}
77
79{
80 fOccurenceBiasingOperation = ProposeOccurenceBiasingOperation(track, callingProcess);
81 return fOccurenceBiasingOperation;
82}
83
85{
86 fFinalStateBiasingOperation = ProposeFinalStateBiasingOperation(track, callingProcess);
87 return fFinalStateBiasingOperation;
88}
89
91{
92 fNonPhysicsBiasingOperation = ProposeNonPhysicsBiasingOperation(track, callingProcess);
93 return fNonPhysicsBiasingOperation;
94}
95
98 G4BiasingAppliedCase biasingCase,
99 G4VBiasingOperation* operationApplied,
100 const G4VParticleChange* particleChangeProduced )
101{
102 fPreviousBiasingAppliedCase = biasingCase;
103 fPreviousAppliedOccurenceBiasingOperation = nullptr;
104 fPreviousAppliedFinalStateBiasingOperation = nullptr;
105 fPreviousAppliedNonPhysicsBiasingOperation = nullptr;
106 switch ( biasingCase )
107 {
108 case BAC_None:
109 break;
110 case BAC_NonPhysics:
111 fPreviousAppliedNonPhysicsBiasingOperation = operationApplied ;
112 break;
113 case BAC_FinalState:
114 fPreviousAppliedFinalStateBiasingOperation = operationApplied;
115 break;
116 case BAC_Occurence:
117 G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
118 "BIAS.MNG.02", JustWarning,
119 "Internal logic error, please report !");
120 break;
121 default:
122 G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
123 "BIAS.MNG.03", JustWarning,
124 "Internal logic error, please report !");
125 }
126 OperationApplied( callingProcess, biasingCase, operationApplied, particleChangeProduced );
127}
128
131 G4BiasingAppliedCase biasingCase,
132 G4VBiasingOperation* occurenceOperationApplied,
133 G4double weightForOccurenceInteraction,
134 G4VBiasingOperation* finalStateOperationApplied,
135 const G4VParticleChange* particleChangeProduced )
136{
137 fPreviousBiasingAppliedCase = biasingCase;
138 fPreviousAppliedOccurenceBiasingOperation = occurenceOperationApplied;
139 fPreviousAppliedFinalStateBiasingOperation = finalStateOperationApplied;
140 OperationApplied( callingProcess, biasingCase, occurenceOperationApplied, weightForOccurenceInteraction, finalStateOperationApplied, particleChangeProduced );
141}
142
144ExitingBiasing( const G4Track* track, const G4BiasingProcessInterface* callingProcess )
145{
146 ExitBiasing( track, callingProcess );
147
148 // -- reset all data members:
149 fOccurenceBiasingOperation = nullptr ;
150 fFinalStateBiasingOperation = nullptr ;
151 fNonPhysicsBiasingOperation = nullptr ;
152 fPreviousProposedOccurenceBiasingOperation = nullptr ;
153 fPreviousProposedFinalStateBiasingOperation = nullptr ;
154 fPreviousProposedNonPhysicsBiasingOperation = nullptr ;
155 fPreviousAppliedOccurenceBiasingOperation = nullptr ;
156 fPreviousAppliedFinalStateBiasingOperation = nullptr ;
157 fPreviousAppliedNonPhysicsBiasingOperation = nullptr ;
158 fPreviousBiasingAppliedCase = BAC_None ;
159}
160
161// -- dummy empty implementations to allow letting arguments visible in the .hh
162// -- but avoiding annoying warning messages about unused variables
163// -- methods to inform operator that its biasing control is over:
178
179// ----------------------------------------------------------------------------
180// -- state machine to get biasing operators messaged at the beginning of runs:
181// ----------------------------------------------------------------------------
182
188
190Notify( G4ApplicationState requestedState )
191{
192 if ( ( fPreviousState == G4State_Idle )
193 && ( requestedState == G4State_GeomClosed ) )
194 {
195 for ( auto i = 0; i < (G4int)G4VBiasingOperator::fOperators.Size(); ++i )
196 {
197 G4VBiasingOperator::fOperators[i]->StartRun();
198 }
199 }
200 fPreviousState = requestedState;
201
202 return true;
203}
G4ApplicationState
@ G4State_Idle
@ G4State_GeomClosed
@ G4State_PreInit
G4BiasingAppliedCase
@ BAC_Occurence
@ BAC_FinalState
@ BAC_NonPhysics
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4bool Notify(G4ApplicationState requestedState)
const G4String & GetName() const
virtual void ExitBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
const G4String & GetName() const
friend class G4BiasingOperatorStateNotifier
G4VBiasingOperation * GetProposedOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
void AttachTo(const G4LogicalVolume *)
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators()
G4VBiasingOperation * GetProposedNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
static G4VBiasingOperator * GetBiasingOperator(const G4LogicalVolume *)
void ExitingBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
G4VBiasingOperator(const G4String &name)
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
void ReportOperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
G4VBiasingOperation * GetProposedFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
G4VStateDependent(G4bool bottom=false)