Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BOptnForceCommonTruncatedExp.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//
30
31#include "Randomize.hh"
33
35 : G4VBiasingOperation(name),
36 fNumberOfSharing(0),
37 fProcessToApply(nullptr),
38 fInteractionOccured(false),
39 fMaximumDistance(-1.0)
40{
41 fCommonTruncatedExpLaw = new G4ILawCommonTruncatedExp("ExpLawForOperation"+name);
42 fForceFreeFlightLaw = new G4ILawForceFreeFlight ("FFFLawForOperation"+name);
43
44 fTotalCrossSection = 0.0;
45}
46
48{
49 if ( fCommonTruncatedExpLaw ) delete fCommonTruncatedExpLaw;
50 if ( fForceFreeFlightLaw ) delete fForceFreeFlightLaw;
51}
52
55 G4ForceCondition& proposeForceCondition )
56{
57 if ( callingProcess->GetWrappedProcess() == fProcessToApply )
58 {
59 proposeForceCondition = Forced;
60 return fCommonTruncatedExpLaw;
61 }
62 else
63 {
64 proposeForceCondition = Forced;
65 return fForceFreeFlightLaw;
66 }
67}
68
69
71{
73}
74
75
77 const G4Track* track,
78 const G4Step* step,
79 G4bool& forceFinalState )
80{
81 if ( callingProcess->GetWrappedProcess() != fProcessToApply )
82 {
83 forceFinalState = true;
84 fDummyParticleChange.Initialize( *track );
85 return &fDummyParticleChange;
86 }
87 if ( fInteractionOccured )
88 {
89 forceFinalState = true;
90 fDummyParticleChange.Initialize( *track );
91 return &fDummyParticleChange;
92 }
93
94 // -- checks if process won the GPIL race:
95 G4double processGPIL = callingProcess->GetPostStepGPIL() < callingProcess->GetAlongStepGPIL() ?
96 callingProcess->GetPostStepGPIL() : callingProcess->GetAlongStepGPIL() ;
97 if ( processGPIL <= step->GetStepLength() )
98 {
99 // -- if process won, wrapped process produces the final state.
100 // -- In this case, the weight for occurrence biasing is applied
101 // -- by the callingProcess, at exit of present method. This is
102 // -- selected by "forceFinalState = false":
103 forceFinalState = false;
104 fInteractionOccured = true;
105 return callingProcess->GetWrappedProcess()->PostStepDoIt( *track, *step );
106 }
107 else
108 {
109 forceFinalState = true;
110 fDummyParticleChange.Initialize( *track );
111 return &fDummyParticleChange;
112 }
113}
114
115
117{
118 fTotalCrossSection += crossSection;
119 fCrossSections[process] = crossSection;
120 fNumberOfSharing = fCrossSections.size();
121}
122
123
125{
126 fCrossSections.clear();
127 fTotalCrossSection = 0.0;
128 fNumberOfSharing = 0;
129 fProcessToApply = 0;
130 fInteractionOccured = false;
131 fInitialMomentum = track->GetMomentum();
132
133 G4VSolid* currentSolid = track->GetVolume()->GetLogicalVolume()->GetSolid();
135 GetNavigatorForTracking()->
136 GetGlobalToLocalTransform()).TransformPoint(track->GetPosition());
138 GetNavigatorForTracking()->
139 GetGlobalToLocalTransform()).TransformAxis(track->GetMomentumDirection());
140 fMaximumDistance = currentSolid->DistanceToOut(localPosition, localDirection);
141 if ( fMaximumDistance <= DBL_MIN ) fMaximumDistance = 0.0;
142 fCommonTruncatedExpLaw->SetMaximumDistance( fMaximumDistance );
143}
144
145
147{
148 fCrossSections.clear();
149 fTotalCrossSection = 0.0;
150 fNumberOfSharing = 0;
151 fProcessToApply = 0;
152
153 fCommonTruncatedExpLaw->UpdateForStep( step->GetStepLength() );
154 fMaximumDistance = fCommonTruncatedExpLaw->GetMaximumDistance();
155}
156
157
159{
160 fCommonTruncatedExpLaw->SetForceCrossSection( fTotalCrossSection );
161 fCommonTruncatedExpLaw->Sample();
163 fCommonTruncatedExpLaw->SetSelectedProcessXSfraction(fCrossSections[fProcessToApply] / fTotalCrossSection);
164}
165
166
168{
169 G4double sigmaRand = G4UniformRand() * fTotalCrossSection;
170 G4double sigmaSelect = 0.0;
171 for ( std::map< const G4VProcess*, G4double>::const_iterator it = fCrossSections.begin();
172 it != fCrossSections.end();
173 it++)
174 {
175 sigmaSelect += (*it).second;
176 if ( sigmaRand <= sigmaSelect )
177 {
178 fProcessToApply = (*it).first;
179 break;
180 }
181 }
182}
G4ForceCondition
@ Forced
G4GPILSelection
@ NotCandidateForSelection
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4UniformRand()
Definition: Randomize.hh:52
virtual G4GPILSelection ProposeGPILSelection(const G4GPILSelection processSelection)
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *, G4ForceCondition &)
void AddCrossSection(const G4VProcess *, G4double)
virtual G4VParticleChange * ApplyFinalStateBiasing(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
G4VProcess * GetWrappedProcess() const
void SetSelectedProcessXSfraction(G4double fXS)
G4VSolid * GetSolid() const
virtual void Initialize(const G4Track &track)
Definition: G4Step.hh:62
G4double GetStepLength() const
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4ThreeVector GetMomentum() const
const G4ThreeVector & GetMomentumDirection() const
static G4TransportationManager * GetTransportationManager()
G4double UpdateForStep(G4double truePathLength)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
#define DBL_MIN
Definition: templates.hh:54