Geant4-11
source
processes
biasing
generic
include
G4BOptnChangeCrossSection.hh
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
//
27
//
28
//---------------------------------------------------------------
29
//
30
// G4BOptnChangeCrossSection
31
//
32
// Class Description:
33
// A G4VBiasingOperation to change a process cross-section.
34
//
35
//
36
//---------------------------------------------------------------
37
// Initial version Nov. 2013 M. Verderi
38
39
40
#ifndef G4BOptnChangeCrossSection_hh
41
#define G4BOptnChangeCrossSection_hh 1
42
43
#include "
G4VBiasingOperation.hh
"
44
class
G4InteractionLawPhysical
;
45
46
class
G4BOptnChangeCrossSection
:
public
G4VBiasingOperation
{
47
public
:
48
// -- Constructor :
49
G4BOptnChangeCrossSection
(
G4String
name
);
50
// -- destructor:
51
virtual
~G4BOptnChangeCrossSection
();
52
53
public
:
54
// -- Methods from G4VBiasingOperation interface:
55
// ----------------------------------------------
56
// -- Used:
57
virtual
const
G4VBiasingInteractionLaw
*
ProvideOccurenceBiasingInteractionLaw
(
const
G4BiasingProcessInterface
*,
58
G4ForceCondition
& proposeForceCondition );
59
60
// -- Unused:
61
virtual
G4VParticleChange
*
ApplyFinalStateBiasing
(
const
G4BiasingProcessInterface
*,
62
const
G4Track
*,
63
const
G4Step
*,
64
G4bool
& ) {
return
0;}
65
virtual
G4double
DistanceToApplyOperation
(
const
G4Track
*,
66
G4double
,
67
G4ForceCondition
*) {
return
DBL_MAX
;}
68
virtual
G4VParticleChange
*
GenerateBiasingFinalState
(
const
G4Track
*,
69
const
G4Step
* ) {
return
0;}
70
71
public
:
72
// -- Additional methods, specific to this class:
73
// ----------------------------------------------
74
// -- return concrete type of interaction law:
75
G4InteractionLawPhysical
*
GetBiasedExponentialLaw
() {
return
fBiasedExponentialLaw
;}
76
// -- set biased cross-section:
77
void
SetBiasedCrossSection
(
G4double
xst,
bool
updateInteractionLength =
false
);
78
G4double
GetBiasedCrossSection
()
const
;
79
// -- Sample underneath distribution:
80
void
Sample
();
81
// -- Update for a made step, without resampling:
82
void
UpdateForStep
(
G4double
stepLength );
83
// -- set/get if interaction occured.
84
// -- Interaction flag turned off in "Sample()" method.
85
G4bool
GetInteractionOccured
()
const
{
return
fInteractionOccured
; }
86
void
SetInteractionOccured
() {
fInteractionOccured
=
true
; }
87
88
89
90
private
:
91
G4InteractionLawPhysical
*
fBiasedExponentialLaw
;
92
G4bool
fInteractionOccured
;
93
};
94
95
#endif
G4ForceCondition
G4ForceCondition
Definition:
G4ForceCondition.hh:41
G4double
double G4double
Definition:
G4Types.hh:83
G4bool
bool G4bool
Definition:
G4Types.hh:86
G4VBiasingOperation.hh
G4BOptnChangeCrossSection
Definition:
G4BOptnChangeCrossSection.hh:46
G4BOptnChangeCrossSection::~G4BOptnChangeCrossSection
virtual ~G4BOptnChangeCrossSection()
Definition:
G4BOptnChangeCrossSection.cc:38
G4BOptnChangeCrossSection::GetBiasedCrossSection
G4double GetBiasedCrossSection() const
Definition:
G4BOptnChangeCrossSection.cc:54
G4BOptnChangeCrossSection::G4BOptnChangeCrossSection
G4BOptnChangeCrossSection(G4String name)
Definition:
G4BOptnChangeCrossSection.cc:31
G4BOptnChangeCrossSection::fBiasedExponentialLaw
G4InteractionLawPhysical * fBiasedExponentialLaw
Definition:
G4BOptnChangeCrossSection.hh:91
G4BOptnChangeCrossSection::ProvideOccurenceBiasingInteractionLaw
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *, G4ForceCondition &proposeForceCondition)
Definition:
G4BOptnChangeCrossSection.cc:43
G4BOptnChangeCrossSection::DistanceToApplyOperation
virtual G4double DistanceToApplyOperation(const G4Track *, G4double, G4ForceCondition *)
Definition:
G4BOptnChangeCrossSection.hh:65
G4BOptnChangeCrossSection::GenerateBiasingFinalState
virtual G4VParticleChange * GenerateBiasingFinalState(const G4Track *, const G4Step *)
Definition:
G4BOptnChangeCrossSection.hh:68
G4BOptnChangeCrossSection::Sample
void Sample()
Definition:
G4BOptnChangeCrossSection.cc:59
G4BOptnChangeCrossSection::UpdateForStep
void UpdateForStep(G4double stepLength)
Definition:
G4BOptnChangeCrossSection.cc:65
G4BOptnChangeCrossSection::GetInteractionOccured
G4bool GetInteractionOccured() const
Definition:
G4BOptnChangeCrossSection.hh:85
G4BOptnChangeCrossSection::fInteractionOccured
G4bool fInteractionOccured
Definition:
G4BOptnChangeCrossSection.hh:92
G4BOptnChangeCrossSection::SetBiasedCrossSection
void SetBiasedCrossSection(G4double xst, bool updateInteractionLength=false)
Definition:
G4BOptnChangeCrossSection.cc:48
G4BOptnChangeCrossSection::GetBiasedExponentialLaw
G4InteractionLawPhysical * GetBiasedExponentialLaw()
Definition:
G4BOptnChangeCrossSection.hh:75
G4BOptnChangeCrossSection::ApplyFinalStateBiasing
virtual G4VParticleChange * ApplyFinalStateBiasing(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
Definition:
G4BOptnChangeCrossSection.hh:61
G4BOptnChangeCrossSection::SetInteractionOccured
void SetInteractionOccured()
Definition:
G4BOptnChangeCrossSection.hh:86
G4BiasingProcessInterface
Definition:
G4BiasingProcessInterface.hh:69
G4InteractionLawPhysical
Definition:
G4InteractionLawPhysical.hh:45
G4Step
Definition:
G4Step.hh:62
G4String
Definition:
G4String.hh:62
G4Track
Definition:
G4Track.hh:67
G4VBiasingInteractionLaw
Definition:
G4VBiasingInteractionLaw.hh:61
G4VBiasingOperation
Definition:
G4VBiasingOperation.hh:76
G4VParticleChange
Definition:
G4VParticleChange.hh:71
G4InuclParticleNames::name
const char * name(G4int ptype)
Definition:
G4InuclParticleNames.hh:76
DBL_MAX
#define DBL_MAX
Definition:
templates.hh:62
Generated by
1.9.3