Skip to content
Snippets Groups Projects

Compare revisions

Changes are shown as if the source revision was being merged into the target revision. Learn more about comparing revisions.

Source

Select target project
No results found

Target

Select target project
  • podlesny/dune-tectonic
  • agnumpde/dune-tectonic
2 results
Show changes
Commits on Source (1)
...@@ -32,6 +32,7 @@ AdaptiveTimeStepper<Factory, Updaters, ErrorNorm>::AdaptiveTimeStepper( ...@@ -32,6 +32,7 @@ AdaptiveTimeStepper<Factory, Updaters, ErrorNorm>::AdaptiveTimeStepper(
parset_(parset), parset_(parset),
globalFriction_(globalFriction), globalFriction_(globalFriction),
current_(current), current_(current),
R1_(),
externalForces_(externalForces), externalForces_(externalForces),
mustRefine_(mustRefine), mustRefine_(mustRefine),
errorNorm_(errorNorm) {} errorNorm_(errorNorm) {}
...@@ -42,73 +43,68 @@ bool AdaptiveTimeStepper<Factory, Updaters, ErrorNorm>::reachedEnd() { ...@@ -42,73 +43,68 @@ bool AdaptiveTimeStepper<Factory, Updaters, ErrorNorm>::reachedEnd() {
} }
template <class Factory, class Updaters, class ErrorNorm> template <class Factory, class Updaters, class ErrorNorm>
bool AdaptiveTimeStepper<Factory, Updaters, ErrorNorm>::coarsen() { IterationRegister AdaptiveTimeStepper<Factory, Updaters, ErrorNorm>::advance() {
bool didCoarsen = false; /*
while (relativeTime_ + relativeTau_ < 1.0) { | C | We check here if making the step R1 of size tau is a
R2_ = {R1_.updaters.clone(), FixedPointIterationCounter()}; | R1 | R2 | good idea. To check if we can coarsen, we compare
step(R2_, relativeTime_ + relativeTau_, relativeTau_); |F1|F2| | the result of (R1+R2) with C, i.e. two steps of size
tau with one of size 2*tau. To check if we need to
UpdatersWithCount C{current_.clone(), FixedPointIterationCounter()}; refine, we compare the result of (F1+F2) with R1, i.e. two steps
step(C, relativeTime_, 2.0 * relativeTau_); of size tau/2 with one of size tau. The method makes multiple
coarsening/refining attempts, with coarsening coming first. */
if (R1_.updaters == Updaters())
R1_ = step(current_, relativeTime_, relativeTau_);
if (!mustRefine_(C.updaters, R2_.updaters)) { bool didCoarsen = false;
R2_ = {}; iterationRegister_.reset();
R1_ = C; UpdatersWithCount R2;
relativeTau_ *= 2.0; UpdatersWithCount C;
didCoarsen = true; while (relativeTime_ + relativeTau_ <= 1.0) {
} else { R2 = step(R1_.updaters, relativeTime_ + relativeTau_, relativeTau_);
C = step(current_, relativeTime_, 2 * relativeTau_);
if (mustRefine_(R2.updaters, C.updaters))
break; break;
}
}
return didCoarsen;
}
template <class Factory, class Updaters, class ErrorNorm> didCoarsen = true;
void AdaptiveTimeStepper<Factory, Updaters, ErrorNorm>::refine() { relativeTau_ *= 2;
while (true) { R1_ = C;
UpdatersWithCount F1{current_.clone(), FixedPointIterationCounter()}; }
step(F1, relativeTime_, relativeTau_ / 2.0); UpdatersWithCount F1;
UpdatersWithCount F2{F1.updaters.clone(), FixedPointIterationCounter()}; UpdatersWithCount F2;
step(F2, relativeTime_ + relativeTau_ / 2.0, relativeTau_ / 2.0); if (!didCoarsen) {
while (true) {
F1 = step(current_, relativeTime_, relativeTau_ / 2.0);
F2 = step(F1.updaters, relativeTime_ + relativeTau_ / 2.0,
relativeTau_ / 2.0);
if (!mustRefine_(F2.updaters, R1_.updaters))
break;
if (!mustRefine_(R1_.updaters, F2.updaters)) {
break;
} else {
R1_ = F1;
R2_ = F2;
relativeTau_ /= 2.0; relativeTau_ /= 2.0;
R1_ = F1;
R2 = F2;
} }
} }
}
template <class Factory, class Updaters, class ErrorNorm>
IterationRegister AdaptiveTimeStepper<Factory, Updaters, ErrorNorm>::advance() {
if (R2_.updaters == Updaters()) {
R1_ = {current_.clone(), FixedPointIterationCounter()};
step(R1_, relativeTime_, relativeTau_); // counting again upon restart
} else {
R1_ = R2_;
}
iterationRegister_.reset();
if (!coarsen())
refine();
current_ = R1_.updaters;
iterationRegister_.registerFinalCount(R1_.count); iterationRegister_.registerFinalCount(R1_.count);
relativeTime_ += relativeTau_; relativeTime_ += relativeTau_;
current_ = R1_.updaters;
R1_ = R2;
return iterationRegister_; return iterationRegister_;
} }
template <class Factory, class Updaters, class ErrorNorm> template <class Factory, class Updaters, class ErrorNorm>
void AdaptiveTimeStepper<Factory, Updaters, ErrorNorm>::step( typename AdaptiveTimeStepper<Factory, Updaters, ErrorNorm>::UpdatersWithCount
UpdatersWithCount &updatersAndCount, double rTime, double rTau) { AdaptiveTimeStepper<Factory, Updaters, ErrorNorm>::step(
updatersAndCount.count = Updaters const &oldUpdaters, double rTime, double rTau) {
UpdatersWithCount newUpdatersAndCount = {oldUpdaters.clone(), {}};
newUpdatersAndCount.count =
MyCoupledTimeStepper(finalTime_, factory_, parset_, globalFriction_, MyCoupledTimeStepper(finalTime_, factory_, parset_, globalFriction_,
updatersAndCount.updaters, errorNorm_, newUpdatersAndCount.updaters, errorNorm_,
externalForces_).step(rTime, rTau); externalForces_)
iterationRegister_.registerCount(updatersAndCount.count); .step(rTime, rTau);
iterationRegister_.registerCount(newUpdatersAndCount.count);
return newUpdatersAndCount;
} }
#include "adaptivetimestepper_tmpl.cc" #include "adaptivetimestepper_tmpl.cc"
...@@ -38,15 +38,14 @@ class AdaptiveTimeStepper { ...@@ -38,15 +38,14 @@ class AdaptiveTimeStepper {
std::function<bool(Updaters &, Updaters &)> mustRefine); std::function<bool(Updaters &, Updaters &)> mustRefine);
bool reachedEnd(); bool reachedEnd();
bool coarsen();
void refine();
IterationRegister advance(); IterationRegister advance();
double relativeTime_; double relativeTime_;
double relativeTau_; double relativeTau_;
private: private:
void step(UpdatersWithCount &updatersWithCount, double rTime, double rTau); UpdatersWithCount step(Updaters const &oldUpdaters, double rTime,
double rTau);
double finalTime_; double finalTime_;
Factory &factory_; Factory &factory_;
...@@ -54,7 +53,6 @@ class AdaptiveTimeStepper { ...@@ -54,7 +53,6 @@ class AdaptiveTimeStepper {
std::shared_ptr<Nonlinearity> globalFriction_; std::shared_ptr<Nonlinearity> globalFriction_;
Updaters &current_; Updaters &current_;
UpdatersWithCount R1_; UpdatersWithCount R1_;
UpdatersWithCount R2_;
std::function<void(double, Vector &)> externalForces_; std::function<void(double, Vector &)> externalForces_;
std::function<bool(Updaters &, Updaters &)> mustRefine_; std::function<bool(Updaters &, Updaters &)> mustRefine_;
ErrorNorm const &errorNorm_; ErrorNorm const &errorNorm_;
......