diff --git a/src/mysolver.cc b/src/mysolver.cc index 3c2d6a804442d8ccf793bc06eda0bb9c0e0d0286..576e312e37249150685f97cddb256388b2413786 100644 --- a/src/mysolver.cc +++ b/src/mysolver.cc @@ -11,8 +11,8 @@ #include "mysolver.hh" -template <int dim, class MatrixType, class VectorType, class GridType> -MySolver<dim, MatrixType, VectorType, GridType>::MySolver( +template <int dim, class BlockProblemType, class GridType> +MySolver<dim, BlockProblemType, GridType>::MySolver( Dune::ParameterTree const &parset, int refinements, double solver_tolerance, GridType const &grid, Dune::BitSetVector<dim> const &ignoreNodes) : baseEnergyNorm(linearBaseSolverStep), @@ -42,17 +42,17 @@ MySolver<dim, MatrixType, VectorType, GridType>::MySolver( multigridStep->ignoreNodes_ = &ignoreNodes; } -template <int dim, class MatrixType, class VectorType, class GridType> -MySolver<dim, MatrixType, VectorType, GridType>::~MySolver() { +template <int dim, class BlockProblemType, class GridType> +MySolver<dim, BlockProblemType, GridType>::~MySolver() { for (auto &x : transferOperators) delete x; delete multigridStep; } -template <int dim, class MatrixType, class VectorType, class GridType> -typename MySolver<dim, MatrixType, VectorType, GridType>::SolverType * -MySolver<dim, MatrixType, VectorType, GridType>::getSolver() { +template <int dim, class BlockProblemType, class GridType> +typename MySolver<dim, BlockProblemType, GridType>::SolverType * +MySolver<dim, BlockProblemType, GridType>::getSolver() { return multigridStep; } diff --git a/src/mysolver.hh b/src/mysolver.hh index 74eaa62116b9486c62d251de3875cf82b696d086..51f36bae611fef2c1a78bfc51d38a21672f898d9 100644 --- a/src/mysolver.hh +++ b/src/mysolver.hh @@ -11,17 +11,16 @@ #include <dune/tnnmg/iterationsteps/genericnonlineargs.hh> #include <dune/solvers/transferoperators/compressedmultigridtransfer.hh> #include <dune/tnnmg/iterationsteps/tnnmgstep.hh> -#include <dune/tnnmg/problem-classes/convexproblem.hh> -#include "dune/tectonic/myblockproblem.hh" -#include "dune/tectonic/globalnonlinearity.hh" - -template <int dim, class MatrixType, class VectorType, class GridType> +template <int dim, class BlockProblemTypeTEMPLATE, class GridType> class MySolver { +public: + using BlockProblemType = BlockProblemTypeTEMPLATE; + using ConvexProblemType = typename BlockProblemType::ConvexProblemType; + using VectorType = typename BlockProblemType::VectorType; + using MatrixType = typename BlockProblemType::MatrixType; + private: - using ConvexProblemType = ConvexProblem< - Dune::GlobalNonlinearity<MatrixType, VectorType>, MatrixType>; - using BlockProblemType = MyBlockProblem<ConvexProblemType>; using NonlinearSmootherType = GenericNonlinearGS<BlockProblemType>; using SolverType = TruncatedNonsmoothNewtonMultigrid<BlockProblemType, NonlinearSmootherType>; diff --git a/src/mysolver_tmpl.cc b/src/mysolver_tmpl.cc index 5e9809b6a625b9593af2227e9be4c5fb7ccae2bc..a4e7953fdfb4ca84f1dae84895c3a15ec95fcabb 100644 --- a/src/mysolver_tmpl.cc +++ b/src/mysolver_tmpl.cc @@ -13,6 +13,11 @@ #include <dune/istl/bcrsmatrix.hh> #include <dune/istl/bvector.hh> +#include <dune/tnnmg/problem-classes/convexproblem.hh> + +#include <dune/tectonic/globalnonlinearity.hh> +#include <dune/tectonic/myblockproblem.hh> + using SmallVector = Dune::FieldVector<double, DIM>; using SmallMatrix = Dune::FieldMatrix<double, DIM, DIM>; using VectorType = Dune::BlockVector<SmallVector>; @@ -20,4 +25,7 @@ using MatrixType = Dune::BCRSMatrix<SmallMatrix>; using GridType = Dune::ALUGrid<DIM, DIM, Dune::simplex, Dune::nonconforming>; -template class MySolver<DIM, MatrixType, VectorType, GridType>; +template class MySolver< + DIM, MyBlockProblem<ConvexProblem< + Dune::GlobalNonlinearity<MatrixType, VectorType>, MatrixType>>, + GridType>; diff --git a/src/one-body-sample.cc b/src/one-body-sample.cc index 4e2759fb50debb2f6ea0cf5d421a1399f0d49e0f..46311dafb84c505eaf7455843e2d8886e9219e75 100644 --- a/src/one-body-sample.cc +++ b/src/one-body-sample.cc @@ -72,6 +72,7 @@ #include <dune/tnnmg/problem-classes/convexproblem.hh> #include <dune/tectonic/myblockproblem.hh> +#include <dune/tectonic/globalnonlinearity.hh> #include "assemblers.hh" #include "mysolver.hh" @@ -168,6 +169,7 @@ int main(int argc, char *argv[]) { using SingletonMatrixType = Dune::BCRSMatrix<SmallSingletonMatrix>; using VectorType = Dune::BlockVector<SmallVector>; using SingletonVectorType = Dune::BlockVector<Dune::FieldVector<double, 1>>; + using NonlinearityType = Dune::GlobalNonlinearity<MatrixType, VectorType>; auto const E = parset.get<double>("body.E"); auto const nu = parset.get<double>("body.nu"); @@ -390,10 +392,14 @@ int main(int argc, char *argv[]) { // Set up TNNMG solver auto const solverTolerance = parset.get<double>("solver.tolerance"); - MySolver<dims, MatrixType, VectorType, GridType> solverHost( - parset.sub("solver.tnnmg"), refinements, solverTolerance, *grid, - velocityDirichletNodes); - auto multigridStep = solverHost.getSolver(); + using FactoryType = MySolver< + dims, + MyBlockProblem<ConvexProblem< + Dune::GlobalNonlinearity<MatrixType, VectorType>, MatrixType>>, + GridType>; + FactoryType factory(parset.sub("solver.tnnmg"), refinements, + solverTolerance, *grid, velocityDirichletNodes); + auto multigridStep = factory.getSolver(); Solver::VerbosityMode const verbosity = parset.get<bool>("verbose") ? Solver::FULL : Solver::QUIET; @@ -460,14 +466,12 @@ int main(int argc, char *argv[]) { SingletonVectorType const &_alpha) { myGlobalNonlinearity->updateState(_alpha); - using ConvexProblemType = ConvexProblem< - Dune::GlobalNonlinearity<MatrixType, VectorType>, MatrixType>; // FIXME: Do we really need to pass u here? - ConvexProblemType const myConvexProblem(1.0, problem_AM, - *myGlobalNonlinearity, - problem_rhs, _problem_iterate); - MyBlockProblem<ConvexProblemType> velocityProblem(parset, - myConvexProblem); + typename FactoryType::ConvexProblemType const myConvexProblem( + 1.0, problem_AM, *myGlobalNonlinearity, problem_rhs, + _problem_iterate); + typename FactoryType::BlockProblemType velocityProblem(parset, + myConvexProblem); multigridStep->setProblem(_problem_iterate, velocityProblem); velocityProblemSolver.preprocess();