[831] | 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 | // |
---|
[850] | 27 | // $Id: G4ChordFinder.hh,v 1.19 2008/07/15 14:02:06 japost Exp $ |
---|
| 28 | // GEANT4 tag $Name: HEAD $ |
---|
[831] | 29 | // |
---|
| 30 | // |
---|
| 31 | // class G4ChordFinder |
---|
| 32 | // |
---|
| 33 | // Class description: |
---|
| 34 | // |
---|
| 35 | // A class that provides RK integration of motion ODE (as does g4magtr) |
---|
| 36 | // and also has a method that returns an Approximate point on the curve |
---|
| 37 | // near to a (chord) point. |
---|
| 38 | |
---|
| 39 | // History: |
---|
| 40 | // - 25.02.97 John Apostolakis, design and implementation |
---|
| 41 | // - 05.03.97 V. Grichine , makeup to G4 'standard' |
---|
| 42 | // ------------------------------------------------------------------- |
---|
| 43 | |
---|
| 44 | #ifndef G4CHORDFINDER_HH |
---|
| 45 | #define G4CHORDFINDER_HH |
---|
| 46 | |
---|
| 47 | #include "G4MagIntegratorDriver.hh" |
---|
| 48 | #include "G4FieldTrack.hh" |
---|
| 49 | |
---|
| 50 | class G4MagneticField; |
---|
| 51 | |
---|
| 52 | class G4ChordFinder |
---|
| 53 | { |
---|
| 54 | public: // with description |
---|
| 55 | |
---|
| 56 | G4ChordFinder( G4MagInt_Driver* pIntegrationDriver ); |
---|
| 57 | |
---|
| 58 | G4ChordFinder( G4MagneticField* itsMagField, |
---|
| 59 | G4double stepMinimum = 1.0e-2 * mm, |
---|
| 60 | G4MagIntegratorStepper* pItsStepper = 0 ); |
---|
| 61 | // A constructor that creates defaults for all "children" classes. |
---|
| 62 | |
---|
| 63 | virtual ~G4ChordFinder(); |
---|
| 64 | |
---|
| 65 | |
---|
| 66 | G4double AdvanceChordLimited( G4FieldTrack& yCurrent, |
---|
| 67 | G4double stepInitial, |
---|
| 68 | G4double epsStep_Relative, |
---|
| 69 | const G4ThreeVector latestSafetyOrigin, |
---|
| 70 | G4double lasestSafetyRadius); |
---|
| 71 | // Uses ODE solver's driver to find the endpoint that satisfies |
---|
| 72 | // the chord criterion: that d_chord < delta_chord |
---|
| 73 | // -> Returns Length of Step taken. |
---|
[850] | 74 | |
---|
| 75 | G4FieldTrack ApproxCurvePointS(const G4FieldTrack& curveAPointVelocity, |
---|
| 76 | const G4FieldTrack& curveBPointVelocity, |
---|
| 77 | const G4ThreeVector& currentEPoint, |
---|
| 78 | const G4ThreeVector& currentFPoint, |
---|
| 79 | const G4ThreeVector& PointG, |
---|
| 80 | G4bool first,G4double epsStep); |
---|
| 81 | |
---|
[831] | 82 | G4FieldTrack ApproxCurvePointV(const G4FieldTrack& curveAPointVelocity, |
---|
| 83 | const G4FieldTrack& curveBPointVelocity, |
---|
| 84 | const G4ThreeVector& currentEPoint, |
---|
| 85 | G4double epsStep); |
---|
| 86 | |
---|
[850] | 87 | inline G4double InvParabolic( const G4double xa, const G4double ya, |
---|
| 88 | const G4double xb, const G4double yb, |
---|
| 89 | const G4double xc, const G4double yc ); |
---|
| 90 | |
---|
[831] | 91 | inline G4double GetDeltaChord() const; |
---|
| 92 | inline void SetDeltaChord(G4double newval); |
---|
| 93 | |
---|
| 94 | inline void SetChargeMomentumMass(G4double pCharge, // in e+ units |
---|
| 95 | G4double pMomentum, |
---|
| 96 | G4double pMass ); |
---|
| 97 | // Function to inform integration driver of charge, speed. |
---|
| 98 | |
---|
| 99 | inline void SetIntegrationDriver(G4MagInt_Driver* IntegrationDriver); |
---|
| 100 | inline G4MagInt_Driver* GetIntegrationDriver(); |
---|
| 101 | // Access and set Driver. |
---|
| 102 | |
---|
| 103 | inline void ResetStepEstimate(); |
---|
| 104 | // Clear internal state (last step estimate) |
---|
| 105 | |
---|
| 106 | inline G4int GetNoCalls(); |
---|
| 107 | inline G4int GetNoTrials(); // Total number of trials |
---|
| 108 | inline G4int GetNoMaxTrials(); // Maximum # of trials for one call |
---|
| 109 | // Get statistics about number of calls & trials in FindNextChord |
---|
| 110 | |
---|
| 111 | virtual void PrintStatistics(); |
---|
| 112 | // A report with the above -- and possibly other stats |
---|
| 113 | inline G4int SetVerbose( G4int newvalue=1); |
---|
| 114 | // Set verbosity and return old value |
---|
| 115 | |
---|
| 116 | protected: // ......................................................... |
---|
| 117 | |
---|
| 118 | inline void AccumulateStatistics( G4int noTrials ); |
---|
| 119 | // Accumulate the basic statistics |
---|
| 120 | // - other specialised ones must be kept by derived classes |
---|
| 121 | |
---|
| 122 | inline G4bool AcceptableMissDist(G4double dChordStep) const; |
---|
| 123 | |
---|
| 124 | G4double NewStep( G4double stepTrialOld, |
---|
| 125 | G4double dChordStep, // Current dchord estimate |
---|
| 126 | G4double& stepEstimate_Unconstrained ) ; |
---|
| 127 | |
---|
[850] | 128 | virtual G4double FindNextChord( const G4FieldTrack& yStart, |
---|
[831] | 129 | G4double stepMax, |
---|
| 130 | G4FieldTrack& yEnd, |
---|
| 131 | G4double& dyErr, // Error of endpoint |
---|
| 132 | G4double epsStep, |
---|
| 133 | G4double* pNextStepForAccuracy, // = 0, |
---|
| 134 | const G4ThreeVector latestSafetyOrigin, |
---|
| 135 | G4double latestSafetyRadius |
---|
| 136 | ); |
---|
| 137 | |
---|
| 138 | void PrintDchordTrial(G4int noTrials, |
---|
| 139 | G4double stepTrial, |
---|
| 140 | G4double oldStepTrial, |
---|
| 141 | G4double dChordStep); |
---|
| 142 | public: // no description |
---|
| 143 | void TestChordPrint( G4int noTrials, |
---|
| 144 | G4int lastStepTrial, |
---|
| 145 | G4double dChordStep, |
---|
| 146 | G4double nextStepTrial ); |
---|
| 147 | // Printing for monitoring ... |
---|
| 148 | |
---|
| 149 | inline G4double GetFirstFraction(); // Originally 0.999 |
---|
| 150 | inline G4double GetFractionLast(); // Originally 1.000 |
---|
| 151 | inline G4double GetFractionNextEstimate(); // Originally 0.980 |
---|
| 152 | inline G4double GetMultipleRadius(); // No original value |
---|
| 153 | // Parameters for adapting performance ... use with great care |
---|
| 154 | |
---|
| 155 | public: // with description |
---|
| 156 | void SetFractions_Last_Next( G4double fractLast= 0.90, |
---|
| 157 | G4double fractNext= 0.95 ); |
---|
| 158 | // Parameters for performance ... change with great care |
---|
| 159 | |
---|
| 160 | inline void SetFirstFraction(G4double fractFirst); |
---|
| 161 | // Parameter for performance ... change with great care |
---|
| 162 | |
---|
| 163 | protected: |
---|
| 164 | inline G4double GetLastStepEstimateUnc(); |
---|
| 165 | inline void SetLastStepEstimateUnc( G4double stepEst ); |
---|
| 166 | |
---|
| 167 | private: // ............................................................ |
---|
| 168 | |
---|
| 169 | G4ChordFinder(const G4ChordFinder&); |
---|
| 170 | G4ChordFinder& operator=(const G4ChordFinder&); |
---|
| 171 | // Private copy constructor and assignment operator. |
---|
| 172 | |
---|
| 173 | private: // ............................................................ |
---|
| 174 | // G4int nOK, nBAD; |
---|
| 175 | G4MagInt_Driver* fIntgrDriver; |
---|
| 176 | |
---|
| 177 | const G4double fDefaultDeltaChord; // SET in G4ChordFinder.cc = 0.25 mm |
---|
| 178 | |
---|
| 179 | G4double fDeltaChord; // Maximum miss distance |
---|
| 180 | |
---|
| 181 | G4double fLastStepEstimate_Unconstrained; // State information for efficiency |
---|
| 182 | // Variables used in construction/destruction |
---|
| 183 | G4bool fAllocatedStepper; |
---|
| 184 | G4EquationOfMotion* fEquation; |
---|
| 185 | G4MagIntegratorStepper* fDriversStepper; |
---|
| 186 | |
---|
| 187 | // Parameters |
---|
| 188 | G4double fFirstFraction, fFractionLast, fFractionNextEstimate; |
---|
| 189 | G4double fMultipleRadius; |
---|
| 190 | |
---|
| 191 | // For Statistics |
---|
| 192 | // -- G4int fNoTrials, fNoCalls; |
---|
| 193 | G4int fTotalNoTrials_FNC, fNoCalls_FNC, fmaxTrials_FNC; // fnoTimesMaxTrFNC; |
---|
| 194 | G4int fStatsVerbose; // if > 0, print Statistics in destructor |
---|
| 195 | }; |
---|
| 196 | |
---|
| 197 | // Inline function implementation: |
---|
| 198 | |
---|
| 199 | #include "G4ChordFinder.icc" |
---|
| 200 | |
---|
| 201 | #endif // G4CHORDFINDER_HH |
---|