source: trunk/source/geometry/magneticfield/include/G4ChordFinder.hh@ 834

Last change on this file since 834 was 831, checked in by garnier, 17 years ago

import all except CVS

File size: 8.1 KB
RevLine 
[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//
27// $Id: G4ChordFinder.hh,v 1.17 2006/06/29 18:21:02 gunter Exp $
28// GEANT4 tag $Name: $
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
50class G4MagneticField;
51
52class 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.
74
75 G4FieldTrack ApproxCurvePointV(const G4FieldTrack& curveAPointVelocity,
76 const G4FieldTrack& curveBPointVelocity,
77 const G4ThreeVector& currentEPoint,
78 G4double epsStep);
79
80 inline G4double GetDeltaChord() const;
81 inline void SetDeltaChord(G4double newval);
82
83 inline void SetChargeMomentumMass(G4double pCharge, // in e+ units
84 G4double pMomentum,
85 G4double pMass );
86 // Function to inform integration driver of charge, speed.
87
88 inline void SetIntegrationDriver(G4MagInt_Driver* IntegrationDriver);
89 inline G4MagInt_Driver* GetIntegrationDriver();
90 // Access and set Driver.
91
92 inline void ResetStepEstimate();
93 // Clear internal state (last step estimate)
94
95 inline G4int GetNoCalls();
96 inline G4int GetNoTrials(); // Total number of trials
97 inline G4int GetNoMaxTrials(); // Maximum # of trials for one call
98 // Get statistics about number of calls & trials in FindNextChord
99
100 virtual void PrintStatistics();
101 // A report with the above -- and possibly other stats
102 inline G4int SetVerbose( G4int newvalue=1);
103 // Set verbosity and return old value
104
105 protected: // .........................................................
106
107 inline void AccumulateStatistics( G4int noTrials );
108 // Accumulate the basic statistics
109 // - other specialised ones must be kept by derived classes
110
111 inline G4bool AcceptableMissDist(G4double dChordStep) const;
112
113 G4double NewStep( G4double stepTrialOld,
114 G4double dChordStep, // Current dchord estimate
115 G4double& stepEstimate_Unconstrained ) ;
116
117 virtual G4double FindNextChord( const G4FieldTrack yStart,
118 G4double stepMax,
119 G4FieldTrack& yEnd,
120 G4double& dyErr, // Error of endpoint
121 G4double epsStep,
122 G4double* pNextStepForAccuracy, // = 0,
123 const G4ThreeVector latestSafetyOrigin,
124 G4double latestSafetyRadius
125 );
126
127 void PrintDchordTrial(G4int noTrials,
128 G4double stepTrial,
129 G4double oldStepTrial,
130 G4double dChordStep);
131 public: // no description
132 void TestChordPrint( G4int noTrials,
133 G4int lastStepTrial,
134 G4double dChordStep,
135 G4double nextStepTrial );
136 // Printing for monitoring ...
137
138 inline G4double GetFirstFraction(); // Originally 0.999
139 inline G4double GetFractionLast(); // Originally 1.000
140 inline G4double GetFractionNextEstimate(); // Originally 0.980
141 inline G4double GetMultipleRadius(); // No original value
142 // Parameters for adapting performance ... use with great care
143
144 public: // with description
145 void SetFractions_Last_Next( G4double fractLast= 0.90,
146 G4double fractNext= 0.95 );
147 // Parameters for performance ... change with great care
148
149 inline void SetFirstFraction(G4double fractFirst);
150 // Parameter for performance ... change with great care
151
152 protected:
153 inline G4double GetLastStepEstimateUnc();
154 inline void SetLastStepEstimateUnc( G4double stepEst );
155
156 private: // ............................................................
157
158 G4ChordFinder(const G4ChordFinder&);
159 G4ChordFinder& operator=(const G4ChordFinder&);
160 // Private copy constructor and assignment operator.
161
162 private: // ............................................................
163 // G4int nOK, nBAD;
164 G4MagInt_Driver* fIntgrDriver;
165
166 const G4double fDefaultDeltaChord; // SET in G4ChordFinder.cc = 0.25 mm
167
168 G4double fDeltaChord; // Maximum miss distance
169
170 G4double fLastStepEstimate_Unconstrained; // State information for efficiency
171 // Variables used in construction/destruction
172 G4bool fAllocatedStepper;
173 G4EquationOfMotion* fEquation;
174 G4MagIntegratorStepper* fDriversStepper;
175
176 // Parameters
177 G4double fFirstFraction, fFractionLast, fFractionNextEstimate;
178 G4double fMultipleRadius;
179
180 // For Statistics
181 // -- G4int fNoTrials, fNoCalls;
182 G4int fTotalNoTrials_FNC, fNoCalls_FNC, fmaxTrials_FNC; // fnoTimesMaxTrFNC;
183 G4int fStatsVerbose; // if > 0, print Statistics in destructor
184};
185
186// Inline function implementation:
187
188#include "G4ChordFinder.icc"
189
190#endif // G4CHORDFINDER_HH
Note: See TracBrowser for help on using the repository browser.