File indexing completed on 2025-01-18 09:58:30
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 #ifndef G4INTERPOLATION_DRIVER_HH
0036 #define G4INTERPOLATION_DRIVER_HH
0037
0038 #include "G4FieldUtils.hh"
0039 #include "G4RKIntegrationDriver.hh"
0040 #include "globals.hh"
0041
0042 #include <memory>
0043 #include <vector>
0044
0045 template <class T, G4bool StepperCachesDchord = true>
0046 class G4InterpolationDriver : public G4RKIntegrationDriver<T>
0047 {
0048 public:
0049
0050 G4InterpolationDriver(G4double hminimum,
0051 T* stepper,
0052 G4int numberOfComponents = 6,
0053 G4int statisticsVerbosity = 0);
0054
0055 ~G4InterpolationDriver() override;
0056
0057 G4InterpolationDriver(const G4InterpolationDriver&) = delete;
0058 const G4InterpolationDriver& operator=(const G4InterpolationDriver&) = delete;
0059
0060 G4double AdvanceChordLimited(G4FieldTrack& track,
0061 G4double hstep,
0062 G4double eps,
0063 G4double chordDistance) override;
0064
0065 void OnStartTracking() override;
0066 void OnComputeStep(const G4FieldTrack* = nullptr) override;
0067 G4bool DoesReIntegrate() const override { return false; }
0068
0069
0070
0071 G4bool AccurateAdvance(G4FieldTrack& track,
0072 G4double hstep,
0073 G4double eps,
0074 G4double hinitial = 0) override;
0075
0076
0077
0078
0079 void SetVerboseLevel(G4int level) override;
0080 G4int GetVerboseLevel() const override;
0081
0082 void StreamInfo(std::ostream& os) const override;
0083
0084 protected:
0085
0086 struct InterpStepper
0087 {
0088 std::unique_ptr<T> stepper;
0089 G4double begin;
0090 G4double end;
0091 G4double inverseLength;
0092 };
0093
0094 using StepperIterator = typename std::vector<InterpStepper>::iterator;
0095 using ConstStepperIterator = typename std::vector<InterpStepper>::const_iterator;
0096
0097 virtual G4double OneGoodStep(StepperIterator it,
0098 field_utils::State& y,
0099 field_utils::State& dydx,
0100 G4double& hstep,
0101 G4double eps,
0102 G4double curveLength,
0103 G4FieldTrack* track = nullptr);
0104
0105
0106
0107
0108
0109 void Interpolate(G4double curveLength, field_utils::State& y) const;
0110
0111 void InterpolateImpl(G4double curveLength,
0112 ConstStepperIterator it,
0113 field_utils::State& y) const;
0114
0115 G4double DistChord(const field_utils::State& yBegin,
0116 G4double curveLengthBegin,
0117 const field_utils::State& yEnd,
0118 G4double curveLengthEnd) const;
0119
0120 G4double FindNextChord(const field_utils::State& yBegin,
0121 G4double curveLengthBegin,
0122 field_utils::State& yEnd,
0123 G4double curveLengthEnd,
0124 G4double dChord,
0125 G4double maxChordDistance);
0126
0127 G4double CalcChordStep(G4double stepTrialOld,
0128 G4double dChordStep,
0129 G4double fDeltaChord);
0130
0131 void PrintState() const;
0132
0133 void CheckState() const;
0134
0135 void AccumulateStatistics(G4int noTrials);
0136
0137 protected:
0138
0139 std::vector<InterpStepper> fSteppers;
0140 StepperIterator fLastStepper;
0141 G4bool fKeepLastStepper = false;
0142
0143 G4double fhnext = DBL_MAX;
0144
0145
0146 G4double fMinimumStep;
0147
0148
0149 G4double fChordStepEstimate = DBL_MAX;
0150 const G4double fFractionNextEstimate = 0.98;
0151 const G4double fSmallestCurveFraction = 0.01;
0152
0153 G4int fVerboseLevel;
0154
0155 field_utils::State fdydx;
0156 G4bool fFirstStep = true;
0157
0158 G4int fMaxTrials = 100;
0159 G4int fTotalStepsForTrack = 0;
0160
0161
0162 G4int fTotalNoTrials = 0;
0163 G4int fNoCalls = 0;
0164 G4int fmaxTrials = 0;
0165
0166 using Base = G4RKIntegrationDriver<T>;
0167 };
0168
0169 #include "G4InterpolationDriver.icc"
0170
0171 #endif