Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BOptnForceCommonTruncatedExp.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// G4BOptnForceCommonTruncatedExp
31//
32// Class Description:
33// A G4VBiasingOperation physics-based biasing operation. It
34// handles several processes together, biasing on the total
35// cross-section of these processes (instead of biasing them
36// individually).
37// The biasing interaction law is a truncated exponential
38// one, driven by the total cross-section and which extends
39// in the range [0,L].
40// Process are registered with the AddCrossSection method.
41// As cross-sections are all known at the end of the
42// PostStepGPIL loop, the step limitation is made at the
43// AlongStepGPIL level.
44//
45//---------------------------------------------------------------
46// Initial version Nov. 2013 M. Verderi
47
48
49#ifndef G4BOptnForceCommonTruncatedExp_hh
50#define G4BOptnForceCommonTruncatedExp_hh 1
51
53#include "G4ThreeVector.hh"
57#include <map>
58
60public:
61 // -- Constructor :
63 // -- destructor:
65
66public:
67 // -- Methods from G4VBiasingOperation interface:
68 // -------------------------------------------
69 // -- Used:
72 virtual G4GPILSelection ProposeGPILSelection( const G4GPILSelection processSelection );
74 const G4Track*,
75 const G4Step*,
76 G4bool& );
77 // -- Unused:
80 G4ForceCondition*) {return DBL_MAX;}
82 const G4Step* ) {return 0;}
83
84public:
85 // -- Additional methods, specific to this class:
86 // ----------------------------------------------
87 // -- return concrete type of interaction laws:
89 {
90 return fCommonTruncatedExpLaw;
91 }
93 {
94 return fForceFreeFlightLaw;
95 }
96
97 void Initialize( const G4Track* );
98 void UpdateForStep( const G4Step* );
99 void Sample();
100 const G4ThreeVector& GetInitialMomentum() const { return fInitialMomentum; }
101 G4double GetMaximumDistance() const { return fMaximumDistance; }
103 const G4VProcess* GetProcessToApply() const { return fProcessToApply; }
104 void AddCrossSection( const G4VProcess*, G4double );
105 size_t GetNumberOfSharing() const { return fNumberOfSharing;}
106 void SetInteractionOccured( G4bool b ) { fInteractionOccured = b; }
107 G4bool GetInteractionOccured() const { return fInteractionOccured; }
108
109private:
110 G4ILawCommonTruncatedExp* fCommonTruncatedExpLaw;
111 G4ILawForceFreeFlight* fForceFreeFlightLaw;
112 G4double fTotalCrossSection;
113 std::map < const G4VProcess*, G4double > fCrossSections;
114 size_t fNumberOfSharing;
115 const G4VProcess* fProcessToApply;
116 G4bool fInteractionOccured;
117 G4ThreeVector fInitialMomentum;
118 G4double fMaximumDistance;
119 G4ParticleChangeForNothing fDummyParticleChange;
120};
121
122#endif
G4ForceCondition
G4GPILSelection
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
virtual G4double ProposeAlongStepLimit(const G4BiasingProcessInterface *)
G4ILawCommonTruncatedExp * GetCommonTruncatedExpLaw()
virtual G4GPILSelection ProposeGPILSelection(const G4GPILSelection processSelection)
virtual G4VParticleChange * GenerateBiasingFinalState(const G4Track *, const G4Step *)
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *, G4ForceCondition &)
void AddCrossSection(const G4VProcess *, G4double)
const G4ThreeVector & GetInitialMomentum() const
virtual G4double DistanceToApplyOperation(const G4Track *, G4double, G4ForceCondition *)
virtual G4VParticleChange * ApplyFinalStateBiasing(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
Definition: G4Step.hh:62
#define DBL_MAX
Definition: templates.hh:62