Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:30

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 // G4InterpolationDriver
0027 //
0028 // Class description:
0029 //
0030 // Driver class which uses Runge-Kutta stepper with interpolation property
0031 // to integrate track with error control
0032 
0033 // Created: D.Sorokin, 2018
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* /*track*/ = nullptr) override;
0067     G4bool DoesReIntegrate() const override { return false; }
0068       // Interpolation driver does not recalculate when AccurateAdvance
0069       // is called -- reintegration would require other calls
0070 
0071     G4bool AccurateAdvance(G4FieldTrack& track,
0072                            G4double hstep,
0073                            G4double eps,  // Requested y_err/hstep
0074                            G4double hinitial = 0) override;
0075       // Integrates ODE from current s (s=s0) to s=s0+h with accuracy eps.
0076       // On output track is replaced by value at end of interval.
0077       // The concept is similar to the odeint routine from NRC p.721-722.
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       // This takes one Step that is of size htry, or as large
0105       // as possible while satisfying the accuracy criterion of:
0106       //     yerr < eps * |y_end-y_start|
0107       // return hdid
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       // Memory of last good step size for integration
0145 
0146     G4double fMinimumStep;
0147       // Minimum Step allowed in a Step (in units of length)   // Parameter
0148 
0149     G4double fChordStepEstimate = DBL_MAX;
0150     const G4double fFractionNextEstimate = 0.98;  // Constant
0151     const G4double fSmallestCurveFraction = 0.01;  // Constant
0152 
0153     G4int fVerboseLevel;  // Parameter
0154 
0155     field_utils::State fdydx;
0156     G4bool fFirstStep = true;
0157 
0158     G4int fMaxTrials = 100;  // Constant
0159     G4int fTotalStepsForTrack = 0;
0160 
0161     // statistics
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