Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
45
47public:
48 // -- Constructor :
50 // -- destructor:
52
53public:
54 // -- Methods from G4VBiasingOperation interface:
55 // ----------------------------------------------
56 // -- Used:
58 G4ForceCondition& proposeForceCondition );
59
60 // -- Unused:
62 const G4Track*,
63 const G4Step*,
64 G4bool& ) {return 0;}
67 G4ForceCondition*) {return DBL_MAX;}
69 const G4Step* ) {return 0;}
70
71public:
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);
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
90private:
91 G4InteractionLawPhysical* fBiasedExponentialLaw;
92 G4bool fInteractionOccured;
93};
94
95#endif
G4ForceCondition
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *, G4ForceCondition &proposeForceCondition)
virtual G4double DistanceToApplyOperation(const G4Track *, G4double, G4ForceCondition *)
virtual G4VParticleChange * GenerateBiasingFinalState(const G4Track *, const G4Step *)
void UpdateForStep(G4double stepLength)
void SetBiasedCrossSection(G4double xst, bool updateInteractionLength=false)
G4InteractionLawPhysical * GetBiasedExponentialLaw()
virtual G4VParticleChange * ApplyFinalStateBiasing(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
Definition: G4Step.hh:62
#define DBL_MAX
Definition: templates.hh:62