Geant4 11.3.0
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//
26// G4BOptnForceCommonTruncatedExp
27// --------------------------------------------------------------------
32
33#include "Randomize.hh"
35
39{
40 fCommonTruncatedExpLaw = new G4ILawCommonTruncatedExp("ExpLawForOperation"+name);
41 fForceFreeFlightLaw = new G4ILawForceFreeFlight ("FFFLawForOperation"+name);
42}
43
45{
46 delete fCommonTruncatedExpLaw;
47 delete fForceFreeFlightLaw;
48}
49
52 G4ForceCondition& proposeForceCondition )
53{
54 if ( callingProcess->GetWrappedProcess() == fProcessToApply )
55 {
56 proposeForceCondition = Forced;
57 return fCommonTruncatedExpLaw;
58 }
59 else
60 {
61 proposeForceCondition = Forced;
62 return fForceFreeFlightLaw;
63 }
64}
65
71
74 const G4Track* track,
75 const G4Step* step,
76 G4bool& forceFinalState )
77{
78 if ( callingProcess->GetWrappedProcess() != fProcessToApply )
79 {
80 forceFinalState = true;
81 fDummyParticleChange.Initialize( *track );
82 return &fDummyParticleChange;
83 }
84 if ( fInteractionOccured )
85 {
86 forceFinalState = true;
87 fDummyParticleChange.Initialize( *track );
88 return &fDummyParticleChange;
89 }
90
91 // -- checks if process won the GPIL race:
92 G4double processGPIL = callingProcess->GetPostStepGPIL()
93 < callingProcess->GetAlongStepGPIL()
94 ? callingProcess->GetPostStepGPIL()
95 : callingProcess->GetAlongStepGPIL();
96 if ( processGPIL <= step->GetStepLength() )
97 {
98 // -- if process won, wrapped process produces the final state.
99 // -- In this case, the weight for occurrence biasing is applied
100 // -- by the callingProcess, at exit of present method. This is
101 // -- selected by "forceFinalState = false":
102 forceFinalState = false;
103 fInteractionOccured = true;
104 return callingProcess->GetWrappedProcess()->PostStepDoIt( *track, *step );
105 }
106 else
107 {
108 forceFinalState = true;
109 fDummyParticleChange.Initialize( *track );
110 return &fDummyParticleChange;
111 }
112}
113
115AddCrossSection( const G4VProcess* process, G4double crossSection )
116{
117 fTotalCrossSection += crossSection;
118 fCrossSections[process] = crossSection;
119 fNumberOfSharing = fCrossSections.size();
120}
121
123{
124 fCrossSections.clear();
125 fTotalCrossSection = 0.0;
126 fNumberOfSharing = 0;
127 fProcessToApply = 0;
128 fInteractionOccured = false;
129 fInitialMomentum = track->GetMomentum();
130
131 G4VSolid* currentSolid = track->GetVolume()->GetLogicalVolume()->GetSolid();
133 GetNavigatorForTracking()->
134 GetGlobalToLocalTransform()).TransformPoint(track->GetPosition());
136 GetNavigatorForTracking()->
137 GetGlobalToLocalTransform()).TransformAxis(track->GetMomentumDirection());
138 fMaximumDistance = currentSolid->DistanceToOut(localPosition, localDirection);
139 if ( fMaximumDistance <= DBL_MIN ) { fMaximumDistance = 0.0; }
140 fCommonTruncatedExpLaw->SetMaximumDistance( fMaximumDistance );
141}
142
144{
145 fCrossSections.clear();
146 fTotalCrossSection = 0.0;
147 fNumberOfSharing = 0;
148 fProcessToApply = 0;
149
150 fCommonTruncatedExpLaw->UpdateForStep( step->GetStepLength() );
151 fMaximumDistance = fCommonTruncatedExpLaw->GetMaximumDistance();
152}
153
155{
156 fCommonTruncatedExpLaw->SetForceCrossSection( fTotalCrossSection );
157 fCommonTruncatedExpLaw->Sample();
159 fCommonTruncatedExpLaw->SetSelectedProcessXSfraction(fCrossSections[fProcessToApply] / fTotalCrossSection);
160}
161
163{
164 G4double sigmaRand = G4UniformRand() * fTotalCrossSection;
165 G4double sigmaSelect = 0.0;
166 for ( auto it = fCrossSections.cbegin(); it != fCrossSections.cend(); ++it)
167 {
168 sigmaSelect += (*it).second;
169 if ( sigmaRand <= sigmaSelect )
170 {
171 fProcessToApply = (*it).first;
172 break;
173 }
174 }
175}
G4ForceCondition
@ Forced
G4GPILSelection
@ NotCandidateForSelection
CLHEP::Hep3Vector G4ThreeVector
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 &)
G4VSolid * GetSolid() const
G4double GetStepLength() const
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4ThreeVector GetMomentum() const
const G4ThreeVector & GetMomentumDirection() const
static G4TransportationManager * GetTransportationManager()
G4VBiasingOperation(const G4String &name)
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