Skip to content
Snippets Groups Projects
Commit a19626f5 authored by Elias Pipping's avatar Elias Pipping
Browse files

[Cleanup] Use explicit template instantiation

parent 1a3a36cf
No related branches found
No related tags found
No related merge requests found
...@@ -3,6 +3,7 @@ bin_PROGRAMS = sand-wedge-2D sand-wedge-3D ...@@ -3,6 +3,7 @@ bin_PROGRAMS = sand-wedge-2D sand-wedge-3D
common_sources = \ common_sources = \
assemblers.cc \ assemblers.cc \
boundary_writer.cc \ boundary_writer.cc \
coupledtimestepper.cc \
enumparser.cc \ enumparser.cc \
fixedpointiterator.cc \ fixedpointiterator.cc \
friction_writer.cc \ friction_writer.cc \
......
#ifdef HAVE_CONFIG_H
#include "config.h"
#endif
#include <dune/solvers/norms/energynorm.hh>
#include "coupledtimestepper.hh"
#include "fixedpointiterator.hh"
template <class Factory, class StateUpdater, class VelocityUpdater>
CoupledTimeStepper<Factory, StateUpdater, VelocityUpdater>::CoupledTimeStepper(
double finalTime, Factory &factory, Dune::ParameterTree const &parset,
std::shared_ptr<Nonlinearity> globalFriction,
std::shared_ptr<StateUpdater> stateUpdater,
std::shared_ptr<VelocityUpdater> velocityUpdater,
std::function<void(double, Vector &)> externalForces)
: finalTime_(finalTime),
factory_(factory),
parset_(parset),
globalFriction_(globalFriction),
stateUpdater_(stateUpdater),
velocityUpdater_(velocityUpdater),
externalForces_(externalForces) {}
template <class Factory, class StateUpdater, class VelocityUpdater>
int CoupledTimeStepper<Factory, StateUpdater, VelocityUpdater>::step(
double relativeTime, double relativeTau) {
stateUpdater_->nextTimeStep();
velocityUpdater_->nextTimeStep();
auto const newRelativeTime = relativeTime + relativeTau;
Vector ell;
externalForces_(newRelativeTime, ell);
Matrix velocityMatrix;
Vector velocityRHS;
Vector velocityIterate;
auto const tau = relativeTau * finalTime_;
stateUpdater_->setup(tau);
velocityUpdater_->setup(ell, tau, newRelativeTime, velocityRHS,
velocityIterate, velocityMatrix);
EnergyNorm<Matrix, Vector> const velocityMatrixNorm(velocityMatrix);
FixedPointIterator<Factory, StateUpdater, VelocityUpdater> fixedPointIterator(
factory_, parset_, globalFriction_);
auto const iterations =
fixedPointIterator.run(stateUpdater_, velocityUpdater_, velocityMatrix,
velocityMatrixNorm, velocityRHS, velocityIterate);
return iterations;
}
#include "coupledtimestepper_tmpl.cc"
#ifndef SRC_COUPLEDTIMESTEPPER_HH
#define SRC_COUPLEDTIMESTEPPER_HH
#include <functional>
#include <memory>
#include <dune/common/parametertree.hh>
template <class Factory, class StateUpdater, class VelocityUpdater>
class CoupledTimeStepper {
using Vector = typename Factory::Vector;
using Matrix = typename Factory::Matrix;
using ConvexProblem = typename Factory::ConvexProblem;
using Nonlinearity = typename ConvexProblem::NonlinearityType;
public:
CoupledTimeStepper(double finalTime, Factory &factory,
Dune::ParameterTree const &parset,
std::shared_ptr<Nonlinearity> globalFriction,
std::shared_ptr<StateUpdater> stateUpdater,
std::shared_ptr<VelocityUpdater> velocityUpdater,
std::function<void(double, Vector &)> externalForces);
int step(double relativeTime, double relativeTau);
private:
double finalTime_;
Factory &factory_;
Dune::ParameterTree const &parset_;
std::shared_ptr<Nonlinearity> globalFriction_;
std::shared_ptr<StateUpdater> stateUpdater_;
std::shared_ptr<VelocityUpdater> velocityUpdater_;
std::function<void(double, Vector &)> externalForces_;
};
#endif
#ifndef DIM
#error DIM unset
#endif
#include <dune/common/function.hh>
#include <dune/tnnmg/problem-classes/convexproblem.hh>
#include <dune/tectonic/globalfriction.hh>
#include <dune/tectonic/myblockproblem.hh>
#include "explicitgrid.hh"
#include "explicitvectors.hh"
#include "solverfactory.hh"
#include "state/stateupdater.hh"
#include "timestepping.hh"
using Function = Dune::VirtualFunction<double, double>;
using Factory = SolverFactory<
DIM, MyBlockProblem<ConvexProblem<GlobalFriction<Matrix, Vector>, Matrix>>,
Grid>;
using MyStateUpdater = StateUpdater<ScalarVector, Vector>;
using MyVelocityUpdater = TimeSteppingScheme<Vector, Matrix, Function, DIM>;
template class CoupledTimeStepper<Factory, MyStateUpdater, MyVelocityUpdater>;
...@@ -68,6 +68,7 @@ ...@@ -68,6 +68,7 @@
#include "assemblers.hh" #include "assemblers.hh"
#include "tobool.hh" #include "tobool.hh"
#include "coupledtimestepper.hh"
#include "enumparser.hh" #include "enumparser.hh"
#include "enums.hh" #include "enums.hh"
#include "fixedpointiterator.hh" #include "fixedpointiterator.hh"
...@@ -93,64 +94,6 @@ void initPython() { ...@@ -93,64 +94,6 @@ void initPython() {
Python::run("sys.path.append('" datadir "')"); Python::run("sys.path.append('" datadir "')");
} }
template <class Factory, class StateUpdater, class VelocityUpdater>
class CoupledTimeStepper {
using Vector = typename Factory::Vector;
using Matrix = typename Factory::Matrix;
using ConvexProblem = typename Factory::ConvexProblem;
using Nonlinearity = typename ConvexProblem::NonlinearityType;
public:
CoupledTimeStepper(double finalTime, Factory &factory,
Dune::ParameterTree const &parset,
std::shared_ptr<Nonlinearity> globalFriction,
std::shared_ptr<StateUpdater> stateUpdater,
std::shared_ptr<VelocityUpdater> velocityUpdater,
std::function<void(double, Vector &)> externalForces)
: finalTime_(finalTime),
factory_(factory),
parset_(parset),
globalFriction_(globalFriction),
stateUpdater_(stateUpdater),
velocityUpdater_(velocityUpdater),
externalForces_(externalForces) {}
int step(double relativeTime, double relativeTau) {
stateUpdater_->nextTimeStep();
velocityUpdater_->nextTimeStep();
auto const newRelativeTime = relativeTime + relativeTau;
Vector ell;
externalForces_(newRelativeTime, ell);
Matrix velocityMatrix;
Vector velocityRHS;
Vector velocityIterate;
auto const tau = relativeTau * finalTime_;
stateUpdater_->setup(tau);
velocityUpdater_->setup(ell, tau, newRelativeTime, velocityRHS,
velocityIterate, velocityMatrix);
EnergyNorm<Matrix, Vector> const velocityMatrixNorm(velocityMatrix);
FixedPointIterator<Factory, StateUpdater, VelocityUpdater>
fixedPointIterator(factory_, parset_, globalFriction_);
auto const iterations = fixedPointIterator.run(
stateUpdater_, velocityUpdater_, velocityMatrix, velocityMatrixNorm,
velocityRHS, velocityIterate);
return iterations;
}
private:
double finalTime_;
Factory &factory_;
Dune::ParameterTree const &parset_;
std::shared_ptr<Nonlinearity> globalFriction_;
std::shared_ptr<StateUpdater> stateUpdater_;
std::shared_ptr<VelocityUpdater> velocityUpdater_;
std::function<void(double, Vector &)> externalForces_;
};
int main(int argc, char *argv[]) { int main(int argc, char *argv[]) {
try { try {
Dune::ParameterTree parset; Dune::ParameterTree parset;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment