File indexing completed on 2025-01-18 09:58:31
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 #include "G4Exception.hh"
0032 #include "G4FieldUtils.hh"
0033 #include "G4LineSection.hh"
0034
0035 #include <CLHEP/Units/SystemOfUnits.h>
0036
0037 #include <algorithm>
0038
0039 template <class T, G4bool StepperCachesDchord>
0040 G4InterpolationDriver<T, StepperCachesDchord>::G4InterpolationDriver(
0041 G4double hminimum, T* pStepper, G4int numComponents, G4int statisticsVerbose)
0042 : G4RKIntegrationDriver<T>(pStepper), fMinimumStep(hminimum), fVerboseLevel(statisticsVerbose)
0043 {
0044 if (numComponents != Base::GetStepper()->GetNumberOfVariables()) {
0045 std::ostringstream message;
0046 message << "Driver's number of integrated components " << numComponents
0047 << " != Stepper's number of components " << pStepper->GetNumberOfVariables();
0048 G4Exception("G4InterpolationDriver", "GeomField0002", FatalException, message);
0049 }
0050
0051 for (G4int i = 0; i < Base::GetMaxNoSteps(); ++i) {
0052 fSteppers.push_back(
0053 {std::unique_ptr<T>(
0054 new T(pStepper->GetSpecificEquation(),
0055 pStepper->GetNumberOfVariables())),
0056 DBL_MAX, -DBL_MAX, 0.0});
0057 }
0058
0059 fLastStepper = fSteppers.end();
0060 }
0061
0062 template <class T, G4bool StepperCachesDchord>
0063 G4InterpolationDriver<T, StepperCachesDchord>::~G4InterpolationDriver()
0064 {
0065 #ifdef G4VERBOSE
0066 if (fVerboseLevel > 0) {
0067 G4cout << "G4ChordFinder statistics report: \n"
0068 << " No trials: " << fTotalNoTrials << " No Calls: " << fNoCalls
0069 << " Max-trial: " << fmaxTrials << G4endl;
0070 }
0071 #endif
0072 }
0073
0074 template <class T, G4bool StepperCachesDchord>
0075 void G4InterpolationDriver<T, StepperCachesDchord>::OnStartTracking()
0076 {
0077 fChordStepEstimate = DBL_MAX;
0078 fhnext = DBL_MAX;
0079 fTotalStepsForTrack = 0;
0080 }
0081
0082 template <class T, G4bool StepperCachesDchord>
0083 void G4InterpolationDriver<T, StepperCachesDchord>::OnComputeStep(const G4FieldTrack* )
0084 {
0085 fKeepLastStepper = false;
0086 fFirstStep = true;
0087 fLastStepper = fSteppers.end();
0088 }
0089
0090 template <class T, G4bool StepperCachesDchord>
0091 void G4InterpolationDriver<T, StepperCachesDchord>::SetVerboseLevel(G4int level)
0092 {
0093 fVerboseLevel = level;
0094 }
0095
0096 template <class T, G4bool StepperCachesDchord>
0097 G4int G4InterpolationDriver<T, StepperCachesDchord>::GetVerboseLevel() const
0098 {
0099 return fVerboseLevel;
0100 }
0101
0102 template <class T, G4bool StepperCachesDchord>
0103 void G4InterpolationDriver<T, StepperCachesDchord>::Interpolate(
0104 G4double curveLength, field_utils::State& y) const
0105 {
0106 if (fLastStepper == fSteppers.end()) {
0107 std::ostringstream message;
0108 message << "LOGICK ERROR: fLastStepper == end";
0109 G4Exception("G4InterpolationDriver::Interpolate()", "GeomField1001", FatalException, message);
0110 return;
0111 }
0112
0113 ConstStepperIterator end = fLastStepper + 1;
0114
0115 auto it = std::lower_bound(fSteppers.cbegin(), end, curveLength,
0116 [](const InterpStepper& stepper, G4double value) { return stepper.end < value; });
0117 if (it == end) {
0118 if (curveLength - fLastStepper->end > CLHEP::perMillion) {
0119 std::ostringstream message;
0120 message << "curveLength = " << curveLength << " > " << fLastStepper->end;
0121 G4Exception("G4InterpolationDriver::Interpolate()", "GeomField1001", JustWarning, message);
0122 }
0123
0124 return fLastStepper->stepper->Interpolate(1, y);
0125 }
0126
0127 if (curveLength < it->begin) {
0128 if (it->begin - curveLength > CLHEP::perMillion) {
0129 std::ostringstream message;
0130 message << "curveLength = " << curveLength << " < " << it->begin;
0131 G4Exception("G4InterpolationDriver::Interpolate()", "GeomField1001", JustWarning, message);
0132 }
0133
0134 return it->stepper->Interpolate(0, y);
0135 }
0136
0137 return InterpolateImpl(curveLength, it, y);
0138 }
0139
0140 template <class T, G4bool StepperCachesDchord>
0141 void G4InterpolationDriver<T, StepperCachesDchord>::InterpolateImpl(
0142 G4double curveLength, ConstStepperIterator it, field_utils::State& y) const
0143 {
0144 const G4double tau = (curveLength - it->begin) * it->inverseLength;
0145 return it->stepper->Interpolate(field_utils::clamp(tau, 0., 1.), y);
0146 }
0147
0148 template <class T, G4bool StepperCachesDchord>
0149 G4double G4InterpolationDriver<T, StepperCachesDchord>::DistChord(const field_utils::State& yBegin,
0150 G4double curveLengthBegin, const field_utils::State& yEnd, G4double curveLengthEnd) const
0151 {
0152 if (StepperCachesDchord) {
0153
0154
0155 if (curveLengthBegin == fLastStepper->begin && curveLengthEnd == fLastStepper->end) {
0156 return fLastStepper->stepper
0157 ->DistChord();
0158 }
0159 }
0160
0161 const G4double curveLengthMid = 0.5 * (curveLengthBegin + curveLengthEnd);
0162 field_utils::State yMid;
0163
0164 Interpolate(curveLengthMid, yMid);
0165
0166 return G4LineSection::Distline(field_utils::makeVector(yMid, field_utils::Value3D::Position),
0167 field_utils::makeVector(yBegin, field_utils::Value3D::Position),
0168 field_utils::makeVector(yEnd, field_utils::Value3D::Position));
0169 }
0170
0171 template <class T, G4bool StepperCachesDchord>
0172 G4double G4InterpolationDriver<T, StepperCachesDchord>::AdvanceChordLimited(
0173 G4FieldTrack& track, G4double hstep, G4double epsStep, G4double chordDistance)
0174 {
0175 ++fTotalStepsForTrack;
0176
0177 const G4double curveLengthBegin = track.GetCurveLength();
0178 const G4double hend = std::min(hstep, fChordStepEstimate);
0179 G4double hdid = 0.0;
0180 auto it = fSteppers.begin();
0181 G4double dChordStep = 0.0;
0182
0183 field_utils::State yBegin, y;
0184 track.DumpToArray(yBegin);
0185 track.DumpToArray(y);
0186
0187 if (fFirstStep) {
0188 Base::GetEquationOfMotion()->RightHandSide(y, fdydx);
0189 fFirstStep = false;
0190 }
0191
0192 if (fKeepLastStepper) {
0193 std::swap(*fSteppers.begin(), *fLastStepper);
0194 it = fSteppers.begin();
0195 fLastStepper = it;
0196 hdid = it->end - curveLengthBegin;
0197 if (hdid > hend) {
0198 hdid = hend;
0199 InterpolateImpl(curveLengthBegin + hdid, it, y);
0200 }
0201 else {
0202 field_utils::copy(y, it->stepper->GetYOut());
0203 }
0204
0205 dChordStep = DistChord(yBegin, curveLengthBegin, y, curveLengthBegin + hdid);
0206
0207 ++it;
0208 }
0209
0210
0211 G4double h = fhnext;
0212 for (; hdid < hend && dChordStep < chordDistance && it != fSteppers.end(); ++it) {
0213 h = std::min(h, hstep - hdid);
0214
0215
0216 hdid += OneGoodStep(it, y, fdydx, h, epsStep, curveLengthBegin + hdid, &track);
0217
0218
0219 fLastStepper = it;
0220
0221
0222 dChordStep =
0223 std::max(dChordStep, DistChord(yBegin, curveLengthBegin, y, curveLengthBegin + hdid));
0224 }
0225
0226
0227
0228
0229
0230
0231
0232 if (h > fMinimumStep) {
0233 fhnext = h;
0234 }
0235
0236
0237
0238
0239
0240 hdid =
0241 FindNextChord(yBegin, curveLengthBegin, y, curveLengthBegin + hdid, dChordStep, chordDistance);
0242
0243 const G4double curveLengthEnd = curveLengthBegin + hdid;
0244 fKeepLastStepper = fLastStepper->end - curveLengthEnd > CLHEP::perMillion;
0245 track.LoadFromArray(y, fLastStepper->stepper->GetNumberOfVariables());
0246 track.SetCurveLength(curveLengthBegin + hdid);
0247
0248 return hdid;
0249 }
0250
0251 template <class T, G4bool StepperCachesDchord>
0252 G4double G4InterpolationDriver<T, StepperCachesDchord>::FindNextChord(
0253 const field_utils::State& yBegin, G4double curveLengthBegin, field_utils::State& yEnd,
0254 G4double curveLengthEnd, G4double dChord, G4double chordDistance)
0255 {
0256 G4double hstep = curveLengthEnd - curveLengthBegin;
0257 G4double curveLength = curveLengthEnd;
0258
0259 G4int i = 1;
0260 for (; i < fMaxTrials && dChord > chordDistance && curveLength > fLastStepper->begin; ++i) {
0261
0262 hstep = CalcChordStep(hstep, dChord, chordDistance);
0263
0264
0265 hstep = std::max(hstep, fLastStepper->begin - curveLengthBegin);
0266 curveLength = curveLengthBegin + hstep;
0267
0268
0269 InterpolateImpl(curveLength, fLastStepper, yEnd);
0270
0271
0272 dChord = DistChord(yBegin, curveLengthBegin, yEnd, curveLength);
0273 }
0274
0275
0276
0277 if (dChord > 0.0) {
0278 fChordStepEstimate = hstep * std::sqrt(chordDistance / dChord);
0279 }
0280
0281 if (i == fMaxTrials) {
0282 G4Exception(
0283 "G4InterpolationDriver::FindNextChord()", "GeomField1001", JustWarning, "cannot converge");
0284 }
0285
0286 AccumulateStatistics(i);
0287
0288 return hstep;
0289 }
0290
0291
0292
0293
0294
0295 template <class T, G4bool StepperCachesDchord>
0296 G4double G4InterpolationDriver<T, StepperCachesDchord>::CalcChordStep(
0297 G4double stepTrialOld, G4double dChordStep, G4double chordDistance)
0298 {
0299 const G4double chordStepEstimate = stepTrialOld * std::sqrt(chordDistance / dChordStep);
0300 G4double stepTrial = fFractionNextEstimate * chordStepEstimate;
0301
0302 if (stepTrial <= 0.001 * stepTrialOld) {
0303 if (dChordStep > 1000.0 * chordDistance) {
0304 stepTrial = stepTrialOld * 0.03;
0305 }
0306 else {
0307 if (dChordStep > 100. * chordDistance) {
0308 stepTrial = stepTrialOld * 0.1;
0309 }
0310 else
0311 {
0312 stepTrial = stepTrialOld * 0.5;
0313 }
0314 }
0315 }
0316 else if (stepTrial > 1000.0 * stepTrialOld) {
0317 stepTrial = 1000.0 * stepTrialOld;
0318 }
0319
0320 if (stepTrial == 0.0) {
0321 stepTrial = 0.000001;
0322 }
0323
0324
0325
0326
0327
0328
0329
0330
0331
0332
0333 return stepTrial;
0334 }
0335
0336 template <class T, G4bool StepperCachesDchord>
0337 G4bool G4InterpolationDriver<T, StepperCachesDchord>::AccurateAdvance(
0338 G4FieldTrack& track, G4double hstep, G4double , G4double
0339 )
0340 {
0341 if (hstep == 0.0) {
0342 std::ostringstream message;
0343 message << "Proposed step is zero; hstep = " << hstep << " !";
0344 G4Exception("G4InterpolationDriver::AccurateAdvance()", "GeomField1001", JustWarning, message);
0345 return true;
0346 }
0347
0348 if (hstep < 0) {
0349 std::ostringstream message;
0350 message << "Invalid run condition." << G4endl << "Proposed step is negative; hstep = " << hstep
0351 << "." << G4endl << "Requested step cannot be negative! Aborting event.";
0352 G4Exception(
0353 "G4InterpolationDriver::AccurateAdvance()", "GeomField0003", EventMustBeAborted, message);
0354 return false;
0355 }
0356
0357 const G4double curveLength = track.GetCurveLength();
0358 const G4double curveLengthEnd = curveLength + hstep;
0359
0360 field_utils::State y;
0361 Interpolate(curveLengthEnd, y);
0362
0363 track.LoadFromArray(y, Base::GetStepper()->GetNumberOfVariables());
0364 track.SetCurveLength(curveLengthEnd);
0365
0366 return true;
0367 }
0368
0369
0370
0371
0372
0373
0374
0375
0376
0377
0378
0379
0380
0381
0382 template <class T, G4bool StepperCachesDchord>
0383 G4double G4InterpolationDriver<T, StepperCachesDchord>::OneGoodStep(StepperIterator it,
0384 field_utils::State& y, field_utils::State& dydx, G4double& hstep, G4double epsStep,
0385 G4double curveLength, G4FieldTrack* )
0386
0387 {
0388 G4double error2 = DBL_MAX;
0389 field_utils::State yerr, ytemp, dydxtemp;
0390 G4double h = hstep;
0391
0392 G4int i = 0;
0393 for (; i < fMaxTrials; ++i) {
0394 it->stepper->Stepper(y, dydx, h, ytemp, yerr, dydxtemp);
0395 error2 = field_utils::relativeError2(y, yerr, h, epsStep);
0396
0397 if (error2 <= 1.0) {
0398 hstep = std::max(Base::GrowStepSize2(h, error2), fMinimumStep);
0399 break;
0400 }
0401
0402
0403 if (h <= fMinimumStep) {
0404 hstep = fMinimumStep;
0405 break;
0406 }
0407
0408 h = std::max(Base::ShrinkStepSize2(h, error2), fMinimumStep);
0409 }
0410
0411 if (i == fMaxTrials) {
0412 G4Exception(
0413 "G4InterpolationDriver::OneGoodStep()", "GeomField1001", JustWarning, "cannot converge");
0414 hstep = std::max(Base::ShrinkStepSize2(h, error2), fMinimumStep);
0415 }
0416
0417
0418 it->begin = curveLength;
0419 it->end = curveLength + h;
0420 it->inverseLength = 1. / h;
0421
0422
0423 it->stepper->SetupInterpolation();
0424
0425 field_utils::copy(dydx, dydxtemp);
0426 field_utils::copy(y, ytemp);
0427
0428 return h;
0429 }
0430
0431 template <class T, G4bool StepperCachesDchord>
0432 void G4InterpolationDriver<T, StepperCachesDchord>::PrintState() const
0433 {
0434 using namespace field_utils;
0435 State prevEnd, currBegin;
0436 auto prev = fSteppers.begin();
0437
0438 G4cout << "====== curr state ========" << G4endl;
0439 for (auto i = fSteppers.begin(); i <= fLastStepper; ++i) {
0440 i->stepper->Interpolate(0, currBegin);
0441
0442 G4cout << "cl_begin: " << i->begin << " "
0443 << "cl_end: " << i->end << " ";
0444
0445 if (prev != i) {
0446 prev->stepper->Interpolate(1, prevEnd);
0447 auto prevPos = makeVector(prevEnd, Value3D::Position);
0448 auto currPos = makeVector(currBegin, Value3D::Position);
0449 G4cout << "diff_begin: " << (prevPos - currPos).mag();
0450 }
0451
0452 G4cout << G4endl;
0453 prev = i;
0454 }
0455
0456 const G4double clBegin = fSteppers.begin()->begin;
0457 const G4double clEnd = fLastStepper->end;
0458 const G4double hstep = (clEnd - clBegin) / 10.;
0459 State yBegin, yCurr;
0460 Interpolate(0, yBegin);
0461 for (G4double cl = clBegin; cl <= clEnd + 1e-12; cl += hstep) {
0462 Interpolate(cl, yCurr);
0463 auto d = DistChord(yBegin, clBegin, yCurr, cl);
0464 G4cout << "cl: " << cl << " chord_distance: " << d << G4endl;
0465 }
0466
0467 G4cout << "==========================" << G4endl;
0468 }
0469
0470 template <class T, G4bool StepperCachesDchord>
0471 void G4InterpolationDriver<T, StepperCachesDchord>::CheckState() const
0472 {
0473 G4int smallSteps = 0;
0474 for (auto i = fSteppers.begin(); i <= fLastStepper; ++i) {
0475 G4double stepLength = i->end - i->begin;
0476 if (stepLength < fMinimumStep) {
0477 ++smallSteps;
0478 }
0479 }
0480
0481 if (smallSteps > 1) {
0482 std::ostringstream message;
0483 message << "====== curr state ========\n";
0484 for (auto i = fSteppers.begin(); i <= fLastStepper; ++i) {
0485 message << "cl_begin: " << i->begin << " "
0486 << "cl_end: " << i->end << "\n";
0487 }
0488
0489 G4Exception("G4InterpolationDriver::CheckState()", "GeomField0003", FatalException, message);
0490 }
0491 }
0492
0493 template <class T, G4bool StepperCachesDchord>
0494 void G4InterpolationDriver<T, StepperCachesDchord>::AccumulateStatistics(G4int noTrials)
0495 {
0496 fTotalNoTrials += noTrials;
0497 ++fNoCalls;
0498
0499 if (noTrials > fmaxTrials) {
0500 fmaxTrials = noTrials;
0501 }
0502 }
0503
0504 template <class T, G4bool StepperCachesDchord>
0505 void G4InterpolationDriver<T, StepperCachesDchord>::StreamInfo(std::ostream& os) const
0506 {
0507 os << "State of G4InterpolationDriver: " << std::endl;
0508 os << "--Base state (G4RKIntegrationDriver): " << std::endl;
0509 Base::StreamInfo(os);
0510 os << " fMinimumStep = " << fMinimumStep << std::endl;
0511
0512
0513
0514
0515
0516
0517
0518 os << " Max num of Trials = " << fMaxTrials << std::endl;
0519 os << " Fract Next Estimate = " << fFractionNextEstimate << std::endl;
0520 os << " Smallest Curve Fract= " << fSmallestCurveFraction << std::endl;
0521
0522 os << " VerboseLevel = " << fVerboseLevel << std::endl;
0523 os << " KeepLastStepper = " << fKeepLastStepper << std::endl;
0524 }