Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 #include "G4VBiasingOperator.hh"
27 #include "G4VBiasingOperation.hh"
28 #include "G4VParticleChange.hh"
29 
30 
31 G4MapCache< const G4LogicalVolume*, G4VBiasingOperator* > G4VBiasingOperator::fLogicalToSetupMap;
32 G4VectorCache< G4VBiasingOperator* > G4VBiasingOperator::fOperators;
33 G4Cache< G4BiasingOperatorStateNotifier* > G4VBiasingOperator::fStateNotifier(0);
34 
35 
37  : fName ( name ),
38  fOccurenceBiasingOperation ( nullptr ),
39  fFinalStateBiasingOperation ( nullptr ),
40  fNonPhysicsBiasingOperation ( nullptr ),
41  fPreviousProposedOccurenceBiasingOperation ( nullptr ),
42  fPreviousProposedFinalStateBiasingOperation( nullptr ),
43  fPreviousProposedNonPhysicsBiasingOperation( nullptr ),
44  fPreviousAppliedOccurenceBiasingOperation ( nullptr ),
45  fPreviousAppliedFinalStateBiasingOperation ( nullptr ),
46  fPreviousAppliedNonPhysicsBiasingOperation ( nullptr ),
47  fPreviousBiasingAppliedCase ( BAC_None )
48 {
49  fOperators.Push_back(this);
50 
51  if ( fStateNotifier.Get() == 0 ) fStateNotifier.Put( new G4BiasingOperatorStateNotifier() );
52 
53 }
54 
56 {
57 }
58 
60 {
62  it = fLogicalToSetupMap.Find(logical);
63  if ( it == fLogicalToSetupMap.End() ) fLogicalToSetupMap[logical] = this;
64  else if ( (*it).second != this )
65  {
67  ed << "Biasing operator `" << GetName()
68  << "' can not be attached to Logical volume `"
69  << logical->GetName() << "' which is already used by an other operator !" << G4endl;
70  G4Exception("G4VBiasingOperator::AttachTo(...)",
71  "BIAS.MNG.01",
73  ed);
74  }
75 }
76 
77 
79 {
81  it = fLogicalToSetupMap.Find(logical);
82  if ( it == fLogicalToSetupMap.End() ) return nullptr;
83  else return (*it).second;
84 }
85 
87 {
88  fOccurenceBiasingOperation = ProposeOccurenceBiasingOperation(track, callingProcess);
89  return fOccurenceBiasingOperation;
90 }
91 
93 {
94  fFinalStateBiasingOperation = ProposeFinalStateBiasingOperation(track, callingProcess);
95  return fFinalStateBiasingOperation;
96 }
97 
99 {
100  fNonPhysicsBiasingOperation = ProposeNonPhysicsBiasingOperation(track, callingProcess);
101  return fNonPhysicsBiasingOperation;
102 }
103 
105  G4BiasingAppliedCase biasingCase,
106  G4VBiasingOperation* operationApplied,
107  const G4VParticleChange* particleChangeProduced )
108 {
109  fPreviousBiasingAppliedCase = biasingCase;
110  fPreviousAppliedOccurenceBiasingOperation = nullptr;
111  fPreviousAppliedFinalStateBiasingOperation = nullptr;
112  fPreviousAppliedNonPhysicsBiasingOperation = nullptr;
113  switch ( biasingCase )
114  {
115  case BAC_None:
116  break;
117  case BAC_NonPhysics:
118  fPreviousAppliedNonPhysicsBiasingOperation = operationApplied ;
119  break;
120  case BAC_FinalState:
121  fPreviousAppliedFinalStateBiasingOperation = operationApplied;
122  break;
123  case BAC_Occurence:
124  G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
125  "BIAS.MNG.02",
126  JustWarning,
127  "Internal logic error, please report !");
128  break;
129  default:
130  G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
131  "BIAS.MNG.03",
132  JustWarning,
133  "Internal logic error, please report !");
134  }
135  OperationApplied( callingProcess, biasingCase, operationApplied, particleChangeProduced );
136 }
137 
139  G4BiasingAppliedCase biasingCase,
140  G4VBiasingOperation* occurenceOperationApplied,
141  G4double weightForOccurenceInteraction,
142  G4VBiasingOperation* finalStateOperationApplied,
143  const G4VParticleChange* particleChangeProduced )
144 {
145  fPreviousBiasingAppliedCase = biasingCase;
146  fPreviousAppliedOccurenceBiasingOperation = occurenceOperationApplied;
147  fPreviousAppliedFinalStateBiasingOperation = finalStateOperationApplied;
148  OperationApplied( callingProcess, biasingCase, occurenceOperationApplied, weightForOccurenceInteraction, finalStateOperationApplied, particleChangeProduced );
149 }
150 
151 
152 void G4VBiasingOperator::ExitingBiasing( const G4Track* track, const G4BiasingProcessInterface* callingProcess )
153 {
154  ExitBiasing( track, callingProcess );
155 
156  // -- reset all data members:
157  fOccurenceBiasingOperation = nullptr ;
158  fFinalStateBiasingOperation = nullptr ;
159  fNonPhysicsBiasingOperation = nullptr ;
160  fPreviousProposedOccurenceBiasingOperation = nullptr ;
161  fPreviousProposedFinalStateBiasingOperation = nullptr ;
162  fPreviousProposedNonPhysicsBiasingOperation = nullptr ;
163  fPreviousAppliedOccurenceBiasingOperation = nullptr ;
164  fPreviousAppliedFinalStateBiasingOperation = nullptr ;
165  fPreviousAppliedNonPhysicsBiasingOperation = nullptr ;
166  fPreviousBiasingAppliedCase = BAC_None ;
167 }
168 
169 
170 // -- dummy empty implementations to allow letting arguments visible in the .hh
171 // -- but avoiding annoying warning messages about unused variables
172 // -- methods to inform operator that its biasing control is over:
174 {}
177 {
178 }
182 {
183 }
184 
185 
186 // ----------------------------------------------------------------------------
187 // -- state machine to get biasing operators messaged at the beginning of runs:
188 // ----------------------------------------------------------------------------
189 
192 {
193  fPreviousState = G4State_PreInit;
194 }
195 
197 {}
198 
200 {
201  if ( ( fPreviousState == G4State_Idle ) && ( requestedState == G4State_GeomClosed ) )
202  {
203  for ( size_t i = 0 ; i < G4VBiasingOperator::fOperators.Size() ; i++ ) G4VBiasingOperator::fOperators[i]->StartRun();
204  }
205 
206  fPreviousState = requestedState;
207 
208  return true;
209 }
const XML_Char * name
Definition: expat.h:151
G4VBiasingOperation * GetProposedFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
iterator End()
Definition: G4Cache.hh:411
G4String fName
Definition: G4AttUtils.hh:55
size_type Size()
Definition: G4Cache.hh:159
value_type & Get() const
Definition: G4Cache.hh:282
G4VBiasingOperation * GetProposedOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
G4bool Notify(G4ApplicationState requestedState)
virtual void ExitBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
G4BiasingAppliedCase
bool G4bool
Definition: G4Types.hh:79
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
void AttachTo(const G4LogicalVolume *)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
static G4VBiasingOperator * GetBiasingOperator(const G4LogicalVolume *)
void Push_back(const value_type &val)
Definition: G4Cache.hh:333
iterator Find(const key_type &k)
Definition: G4Cache.hh:417
G4VBiasingOperator(G4String name)
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76
void ReportOperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
G4VBiasingOperation * GetProposedNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
void ExitingBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
void Put(const value_type &val) const
Definition: G4Cache.hh:286
G4ApplicationState
const G4String GetName() const