Geant4
11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BOptnChangeCrossSection.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
// G4BOptnChangeCrossSection
27
// --------------------------------------------------------------------
28
29
#include "
G4BOptnChangeCrossSection.hh
"
30
#include "
G4InteractionLawPhysical.hh
"
31
32
G4BOptnChangeCrossSection::G4BOptnChangeCrossSection
(
const
G4String
& name)
33
:
G4VBiasingOperation
( name )
34
{
35
fBiasedExponentialLaw =
new
G4InteractionLawPhysical
(
"LawForOperation"
+name);
36
}
37
38
G4BOptnChangeCrossSection::~G4BOptnChangeCrossSection
()
39
{
40
delete
fBiasedExponentialLaw;
41
}
42
43
const
G4VBiasingInteractionLaw
*
G4BOptnChangeCrossSection::ProvideOccurenceBiasingInteractionLaw
(
const
G4BiasingProcessInterface
*,
G4ForceCondition
& )
44
{
45
return
fBiasedExponentialLaw;
46
}
47
48
void
G4BOptnChangeCrossSection::SetBiasedCrossSection
(
G4double
xst,
G4bool
updateInteractionLength )
49
{
50
fBiasedExponentialLaw->SetPhysicalCrossSection( xst );
51
if
( updateInteractionLength )
UpdateForStep
( 0.0 );
52
}
53
54
G4double
G4BOptnChangeCrossSection::GetBiasedCrossSection
()
const
55
{
56
return
fBiasedExponentialLaw->GetPhysicalCrossSection();
57
}
58
59
void
G4BOptnChangeCrossSection::Sample
()
60
{
61
fInteractionOccured =
false
;
62
fBiasedExponentialLaw->Sample();
63
}
64
65
void
G4BOptnChangeCrossSection::UpdateForStep
(
G4double
truePathLength )
66
{
67
fBiasedExponentialLaw->UpdateForStep( truePathLength );
68
}
G4BOptnChangeCrossSection.hh
G4ForceCondition
G4ForceCondition
Definition
G4ForceCondition.hh:41
G4InteractionLawPhysical.hh
G4double
double G4double
Definition
G4Types.hh:83
G4bool
bool G4bool
Definition
G4Types.hh:86
G4BOptnChangeCrossSection::~G4BOptnChangeCrossSection
virtual ~G4BOptnChangeCrossSection()
Definition
G4BOptnChangeCrossSection.cc:38
G4BOptnChangeCrossSection::GetBiasedCrossSection
G4double GetBiasedCrossSection() const
Definition
G4BOptnChangeCrossSection.cc:54
G4BOptnChangeCrossSection::G4BOptnChangeCrossSection
G4BOptnChangeCrossSection(const G4String &name)
Definition
G4BOptnChangeCrossSection.cc:32
G4BOptnChangeCrossSection::ProvideOccurenceBiasingInteractionLaw
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *, G4ForceCondition &proposeForceCondition)
Definition
G4BOptnChangeCrossSection.cc:43
G4BOptnChangeCrossSection::Sample
void Sample()
Definition
G4BOptnChangeCrossSection.cc:59
G4BOptnChangeCrossSection::UpdateForStep
void UpdateForStep(G4double stepLength)
Definition
G4BOptnChangeCrossSection.cc:65
G4BOptnChangeCrossSection::SetBiasedCrossSection
void SetBiasedCrossSection(G4double xst, G4bool updateInteractionLength=false)
Definition
G4BOptnChangeCrossSection.cc:48
G4BiasingProcessInterface
Definition
G4BiasingProcessInterface.hh:64
G4InteractionLawPhysical
Definition
G4InteractionLawPhysical.hh:41
G4String
Definition
G4String.hh:62
G4VBiasingInteractionLaw
Definition
G4VBiasingInteractionLaw.hh:54
G4VBiasingOperation::G4VBiasingOperation
G4VBiasingOperation(const G4String &name)
Definition
G4VBiasingOperation.cc:32
geant4-v11.3.0
source
processes
biasing
generic
src
G4BOptnChangeCrossSection.cc
Generated by
1.13.2