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

[Cleanup] Reinvent MySolver interface

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