Newer
Older
#ifndef RIEMANNIAN_TRUST_REGION_SOLVER_HH
#define RIEMANNIAN_TRUST_REGION_SOLVER_HH
#include <dune/istl/bcrsmatrix.hh>
#include <dune/istl/bvector.hh>
#include <dune/solvers/common/boxconstraint.hh>
#include <dune/solvers/norms/h1seminorm.hh>
#include <dune/solvers/solvers/loopsolver.hh>
#include <dune/fufem/functionspacebases/p1nodalbasis.hh>
#include <dune/fufem/functionspacebases/p2nodalbasis.hh>
#include <dune/fufem/functionspacebases/p3nodalbasis.hh>

Oliver Sander
committed
#include "geodesicfeassembler.hh"
#include <dune/grid/utility/globalindexset.hh>

Oliver Sander
committed
#include <dune/gfe/parallel/globalp2mapper.hh>

Oliver Sander
committed
/** \brief Riemannian trust-region solver for geodesic finite-element problems */
template <class GridType, class TargetSpace>
: public IterativeSolver<std::vector<TargetSpace>,
Dune::BitSetVector<TargetSpace::TangentVector::dimension> >
const static int blocksize = TargetSpace::TangentVector::dimension;

Oliver Sander
committed
const static int gridDim = GridType::dimension;
// Centralize the field type here
typedef double field_type;
typedef Dune::BCRSMatrix<Dune::FieldMatrix<field_type, blocksize, blocksize> > MatrixType;
typedef Dune::BlockVector<Dune::FieldVector<field_type, blocksize> > CorrectionType;

Oliver Sander
committed
typedef std::vector<TargetSpace> SolutionType;

Oliver Sander
committed
#ifdef THIRD_ORDER
#error RiemannianTrustRegionSolver cannot currently be used for third-order spaces
#elif defined SECOND_ORDER
typedef Dune::GlobalP2Mapper<typename GridType::LeafGridView> GUIndex;
#else
typedef GlobalUniqueIndex<typename GridType::LeafGridView, gridDim> GUIndex;
#ifdef THIRD_ORDER
typedef P3NodalBasis<typename GridType::LeafGridView,double> BasisType;
#elif defined SECOND_ORDER
typedef P2NodalBasis<typename GridType::LeafGridView,double> BasisType;
#else
typedef P1NodalBasis<typename GridType::LeafGridView,double> BasisType;
#endif
RiemannianTrustRegionSolver()

Oliver Sander
committed
: IterativeSolver<std::vector<TargetSpace>, Dune::BitSetVector<blocksize> >(0,100,NumProc::FULL),
hessianMatrix_(std::auto_ptr<MatrixType>(NULL)), h1SemiNorm_(NULL)

Oliver Sander
committed
/** \brief Set up the solver using a monotone multigrid method as the inner solver */
const GeodesicFEAssembler<BasisType, TargetSpace>* assembler,
const Dune::BitSetVector<blocksize>& dirichletNodes,
int maxTrustRegionSteps,
double initialTrustRegionRadius,
int multigridIterations,
double mgTolerance,
int nu1,
int nu2,
int baseIterations,
double baseTolerance,
bool instrumented);
void setIgnoreNodes(const Dune::BitSetVector<blocksize>& ignoreNodes)
{
ignoreNodes_ = &ignoreNodes;
Dune::shared_ptr<LoopSolver<CorrectionType> > loopSolver = std::dynamic_pointer_cast<LoopSolver<CorrectionType> >(innerSolver_);
assert(loopSolver);
loopSolver->iterationStep_->ignoreNodes_ = ignoreNodes_;
}
void setInitialSolution(const SolutionType& x) DUNE_DEPRECATED {
x_ = x;
}
void setInitialIterate(const SolutionType& x) {
x_ = x;
}
SolutionType getSol() const {return x_;}
protected:
std::unique_ptr<GUIndex> guIndex_;

Oliver Sander
committed
/** \brief The grid */
const GridType* grid_;
/** \brief The solution vector */
SolutionType x_;
/** \brief The initial trust-region radius in the maximum-norm */
double initialTrustRegionRadius_;
/** \brief Maximum number of trust-region steps */
int maxTrustRegionSteps_;
/** \brief Maximum number of multigrid iterations */

Oliver Sander
committed
int innerIterations_;
/** \brief Error tolerance of the multigrid QP solver */

Oliver Sander
committed
double innerTolerance_;
std::auto_ptr<MatrixType> hessianMatrix_;
/** \brief The assembler for the material law */
const GeodesicFEAssembler<BasisType, TargetSpace>* assembler_;

Oliver Sander
committed
/** \brief The solver for the quadratic inner problems */
std::shared_ptr<Solver> innerSolver_;
/** \brief Contains 'true' everywhere -- the trust-region is bounded */
Dune::BitSetVector<blocksize> hasObstacle_;
/** \brief The Dirichlet nodes */
const Dune::BitSetVector<blocksize>* ignoreNodes_;
/** \brief The norm used to measure multigrid convergence */
H1SemiNorm<CorrectionType>* h1SemiNorm_;
/** \brief An L2-norm, really. The H1SemiNorm class is badly named */
std::shared_ptr<H1SemiNorm<CorrectionType> > l2Norm_;
/** \brief If set to true we log convergence speed and other stuff */
bool instrumented_;

Oliver Sander
committed
#include "riemanniantrsolver.cc"