File indexing completed on 2025-01-18 09:59:19
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 #ifndef G4VContinuousProcess_hh
0038 #define G4VContinuousProcess_hh 1
0039
0040 #include "globals.hh"
0041 #include "G4ios.hh"
0042
0043 #include "G4VProcess.hh"
0044
0045 class G4VContinuousProcess : public G4VProcess
0046 {
0047 public:
0048
0049 G4VContinuousProcess(const G4String& aName,
0050 G4ProcessType aType = fNotDefined );
0051 G4VContinuousProcess(G4VContinuousProcess &);
0052
0053 virtual ~G4VContinuousProcess();
0054
0055 G4VContinuousProcess& operator=(const G4VContinuousProcess&) = delete;
0056
0057 virtual G4double AlongStepGetPhysicalInteractionLength(
0058 const G4Track& track,
0059 G4double previousStepSize,
0060 G4double currentMinimumStep,
0061 G4double& proposedSafety,
0062 G4GPILSelection* selection
0063 );
0064
0065 virtual G4VParticleChange* AlongStepDoIt(
0066 const G4Track& ,
0067 const G4Step&
0068 );
0069
0070
0071
0072 virtual G4double PostStepGetPhysicalInteractionLength(
0073 const G4Track&,
0074 G4double,
0075 G4ForceCondition*
0076 ) { return -1.0; }
0077
0078 virtual G4double AtRestGetPhysicalInteractionLength(
0079 const G4Track& ,
0080 G4ForceCondition*
0081 ) { return -1.0; }
0082
0083
0084
0085 virtual G4VParticleChange* AtRestDoIt(
0086 const G4Track& ,
0087 const G4Step&
0088 ) { return 0; }
0089
0090 virtual G4VParticleChange* PostStepDoIt(
0091 const G4Track& ,
0092 const G4Step&
0093 ) { return 0; }
0094
0095 protected:
0096
0097 virtual G4double GetContinuousStepLimit( const G4Track& aTrack,
0098 G4double previousStepSize,
0099 G4double currentMinimumStep,
0100 G4double& currentSafety ) = 0;
0101
0102
0103
0104 inline void SetGPILSelection(G4GPILSelection selection)
0105 { valueGPILSelection = selection; }
0106
0107 inline G4GPILSelection GetGPILSelection() const
0108 { return valueGPILSelection; }
0109
0110 private:
0111
0112 G4VContinuousProcess();
0113
0114
0115 G4GPILSelection valueGPILSelection = CandidateForSelection;
0116
0117
0118 };
0119
0120 #endif