File indexing completed on 2025-01-18 09:57:59
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049 #ifndef G4BOptnForceCommonTruncatedExp_hh
0050 #define G4BOptnForceCommonTruncatedExp_hh 1
0051
0052 #include "G4VBiasingOperation.hh"
0053 #include "G4ThreeVector.hh"
0054 #include "G4ParticleChangeForNothing.hh"
0055 class G4ILawCommonTruncatedExp;
0056 class G4ILawForceFreeFlight;
0057 #include <map>
0058
0059 class G4BOptnForceCommonTruncatedExp : public G4VBiasingOperation {
0060 public:
0061
0062 G4BOptnForceCommonTruncatedExp(G4String name);
0063
0064 virtual ~G4BOptnForceCommonTruncatedExp();
0065
0066 public:
0067
0068
0069
0070 virtual const G4VBiasingInteractionLaw* ProvideOccurenceBiasingInteractionLaw( const G4BiasingProcessInterface*, G4ForceCondition& );
0071 virtual G4double ProposeAlongStepLimit( const G4BiasingProcessInterface* ) { return DBL_MAX; }
0072 virtual G4GPILSelection ProposeGPILSelection( const G4GPILSelection processSelection );
0073 virtual G4VParticleChange* ApplyFinalStateBiasing( const G4BiasingProcessInterface*,
0074 const G4Track*,
0075 const G4Step*,
0076 G4bool& );
0077
0078 virtual G4double DistanceToApplyOperation( const G4Track*,
0079 G4double,
0080 G4ForceCondition*) {return DBL_MAX;}
0081 virtual G4VParticleChange* GenerateBiasingFinalState( const G4Track*,
0082 const G4Step* ) {return 0;}
0083
0084 public:
0085
0086
0087
0088 G4ILawCommonTruncatedExp* GetCommonTruncatedExpLaw()
0089 {
0090 return fCommonTruncatedExpLaw;
0091 }
0092 G4ILawForceFreeFlight* GetForceFreeFlightLaw()
0093 {
0094 return fForceFreeFlightLaw;
0095 }
0096
0097 void Initialize( const G4Track* );
0098 void UpdateForStep( const G4Step* );
0099 void Sample();
0100 const G4ThreeVector& GetInitialMomentum() const { return fInitialMomentum; }
0101 G4double GetMaximumDistance() const { return fMaximumDistance; }
0102 void ChooseProcessToApply();
0103 const G4VProcess* GetProcessToApply() const { return fProcessToApply; }
0104 void AddCrossSection( const G4VProcess*, G4double );
0105 size_t GetNumberOfSharing() const { return fNumberOfSharing;}
0106 void SetInteractionOccured( G4bool b ) { fInteractionOccured = b; }
0107 G4bool GetInteractionOccured() const { return fInteractionOccured; }
0108
0109 private:
0110 G4ILawCommonTruncatedExp* fCommonTruncatedExpLaw;
0111 G4ILawForceFreeFlight* fForceFreeFlightLaw;
0112 G4double fTotalCrossSection;
0113 std::map < const G4VProcess*, G4double > fCrossSections;
0114 size_t fNumberOfSharing;
0115 const G4VProcess* fProcessToApply;
0116 G4bool fInteractionOccured;
0117 G4ThreeVector fInitialMomentum;
0118 G4double fMaximumDistance;
0119 G4ParticleChangeForNothing fDummyParticleChange;
0120 };
0121
0122 #endif