Geant4  10.01.p02
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 #include "G4BiasingTrackData.hh"
31 
32 
36 
37 
39  : fName(name)
40 {
41  fOperators.Push_back(this);
42 
44 
45 }
46 
48 {
49 }
50 
52 {
54  it = fLogicalToSetupMap.Find(logical);
55  if ( it == fLogicalToSetupMap.End() ) fLogicalToSetupMap[logical] = this;
56  else if ( (*it).second != this )
57  {
59  ed << "Biasing operator `" << GetName()
60  << "' can not be attached to Logical volume `"
61  << logical->GetName() << "' which is already used by an other operator !" << G4endl;
62  G4Exception("G4VBiasingOperator::AttachTo(...)",
63  "BIAS.MNG.01",
65  ed);
66  }
67 }
68 
69 
71 {
73  it = fLogicalToSetupMap.Find(logical);
74  if ( it == fLogicalToSetupMap.End() ) return 0;
75  else return (*it).second;
76 }
77 
79 {
82 }
83 
85 {
88 }
89 
91 {
94 }
95 
97 {
99  if ( biasingData != 0 ) return biasingData->GetBirthOperation();
100  else return 0;
101 }
102 
103 
105  G4BiasingAppliedCase biasingCase,
106  G4VBiasingOperation* operationApplied,
107  const G4VParticleChange* particleChangeProduced )
108 {
109  fPreviousBiasingAppliedCase = biasingCase;
113  switch ( biasingCase )
114  {
115  case BAC_None:
116  break;
117  case BAC_NonPhysics:
118  fPreviousAppliedNonPhysicsBiasingOperation = operationApplied ;
119  break;
120  case BAC_DenyInteraction:
122  break;
123  case BAC_FinalState:
125  break;
126  case BAC_Occurence:
127  G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
128  "BIAS.MNG.02",
129  JustWarning,
130  "Internal logic error, please report !");
131  break;
132  default:
133  G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
134  "BIAS.MNG.03",
135  JustWarning,
136  "Internal logic error, please report !");
137  }
138  OperationApplied( callingProcess, biasingCase, operationApplied, particleChangeProduced );
139 }
140 
142  G4BiasingAppliedCase biasingCase,
143  G4VBiasingOperation* occurenceOperationApplied,
144  G4double weightForOccurenceInteraction,
145  G4VBiasingOperation* finalStateOperationApplied,
146  const G4VParticleChange* particleChangeProduced )
147 {
148  fPreviousBiasingAppliedCase = biasingCase;
149  fPreviousAppliedOccurenceBiasingOperation = occurenceOperationApplied;
150  fPreviousAppliedFinalStateBiasingOperation = finalStateOperationApplied;
151  OperationApplied( callingProcess, biasingCase, occurenceOperationApplied, weightForOccurenceInteraction, finalStateOperationApplied, particleChangeProduced );
152 }
153 
154 
155 void G4VBiasingOperator::ExitingBiasing( const G4Track* track, const G4BiasingProcessInterface* callingProcess )
156 {
157  ExitBiasing( track, callingProcess );
158 
159  // -- reset all data members:
170 }
171 
172 
174  const G4VBiasingOperation* operationApplied,
175  const G4VParticleChange* particleChangeProduced)
176 {
177  for (G4int i2nd = 0; i2nd < particleChangeProduced->GetNumberOfSecondaries (); i2nd++)
178  new G4BiasingTrackData( particleChangeProduced->GetSecondary (i2nd),
179  operationApplied,
180  this,
181  callingProcess);
182 }
183 
185 {
187  if ( biasingData != 0 ) delete biasingData;
188 }
189 
190 
191 // -- dummy empty implementations to allow letting arguments visible in the .hh
192 // -- but avoiding annoying warning messages about unused variables
193 // -- methods to inform operator that its biasing control is over:
195 {}
198 {}
202 {}
203 
204 
205 // ----------------------------------------------------------------------------
206 // -- state machine to get biasing operators messaged at the beginning of runs:
207 // ----------------------------------------------------------------------------
208 
211 {
213 }
214 
216 {}
217 
219 {
220  if ( ( fPreviousState == G4State_Idle ) && ( requestedState == G4State_GeomClosed ) )
221  {
222  for ( size_t i = 0 ; i < G4VBiasingOperator::fOperators.Size() ; i++ ) G4VBiasingOperator::fOperators[i]->StartRun();
223  }
224 
225  fPreviousState = requestedState;
226 
227  return true;
228 }
const G4VBiasingOperation * fPreviousProposedOccurenceBiasingOperation
G4VBiasingOperation * GetProposedFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
G4int GetNumberOfSecondaries() const
G4Track * GetSecondary(G4int anIndex) const
static G4BiasingTrackDataStore * GetInstance()
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4String GetName() const
iterator End()
Definition: G4Cache.hh:411
G4String fName
Definition: G4AttUtils.hh:55
size_type Size()
Definition: G4Cache.hh:159
G4String name
Definition: TRTMaterials.hh:40
value_type & Get() const
Definition: G4Cache.hh:282
const G4VBiasingOperation * GetBirthOperation() const
G4BiasingAppliedCase fPreviousBiasingAppliedCase
G4VBiasingOperation * GetProposedOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
G4bool Notify(G4ApplicationState requestedState)
int G4int
Definition: G4Types.hh:78
G4BiasingTrackData * GetBiasingTrackData(const G4Track *track)
virtual void ExitBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
G4BiasingAppliedCase
const G4VBiasingOperation * fPreviousAppliedOccurenceBiasingOperation
static G4VectorCache< G4VBiasingOperator * > fOperators
static G4MapCache< const G4LogicalVolume *, G4VBiasingOperator * > fLogicalToSetupMap
bool G4bool
Definition: G4Types.hh:79
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
const G4VBiasingOperation * fPreviousAppliedNonPhysicsBiasingOperation
const G4VBiasingOperation * fPreviousProposedFinalStateBiasingOperation
void AttachTo(const G4LogicalVolume *)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4VBiasingOperation * fPreviousAppliedFinalStateBiasingOperation
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
static G4VBiasingOperator * GetBiasingOperator(const G4LogicalVolume *)
void ForgetTrack(const G4Track *track)
void Push_back(const value_type &val)
Definition: G4Cache.hh:333
iterator Find(const key_type &k)
Definition: G4Cache.hh:417
G4VBiasingOperator(G4String name)
void RememberSecondaries(const G4BiasingProcessInterface *callingProcess, const G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
G4VBiasingOperation * fFinalStateBiasingOperation
#define G4endl
Definition: G4ios.hh:61
static G4Cache< G4BiasingOperatorStateNotifier * > fStateNotifier
G4VBiasingOperation * fOccurenceBiasingOperation
double G4double
Definition: G4Types.hh:76
void ReportOperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
const G4VBiasingOperation * fPreviousProposedNonPhysicsBiasingOperation
G4VBiasingOperation * GetProposedNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
void ExitingBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
const G4VBiasingOperation * GetBirthOperation(const G4Track *)
void Put(const value_type &val) const
Definition: G4Cache.hh:286
G4VBiasingOperation * fNonPhysicsBiasingOperation
G4ApplicationState
const G4String GetName() const