-
Jonathan Youett authoredJonathan Youett authored
Code owners
Assign users and groups as approvers for specific file changes. Learn more.
mmgstep.hh 2.26 KiB
#ifndef DUNE_MONOTONE_MULTIGRID_STEP_HH
#define DUNE_MONOTONE_MULTIGRID_STEP_HH
#include <dune/solvers/iterationsteps/multigridstep.hh>
#include <dune/solvers/iterationsteps/projectedblockgsstep.hh>
#include <dune/solvers/transferoperators/multigridtransfer.hh>
#include <dune/solvers/common/boxconstraint.hh>
#include <dune/solvers/transferoperators/obstaclerestrictor.hh>
/** \brief The general monotone multigrid solver
*/
template<class MatrixType, class VectorType>
class MonotoneMGStep : public MultigridStep<MatrixType, VectorType>
{
static const int dim = VectorType::block_type::dimension;
typedef typename VectorType::field_type field_type;
typedef std::vector<BoxConstraint<field_type,dim> > ObstacleVectorType;
public:
MonotoneMGStep() : hasObstacleHierarchy_(0), obstacleHierarchy_(0)
{}
MonotoneMGStep(const MatrixType& mat,
VectorType& x,
VectorType& rhs)
: MultigridStep<MatrixType, VectorType>(mat, x, rhs),
hasObstacleHierarchy_(0), obstacleHierarchy_(0)
{
oldCritical.resize(x.size(), false);
}
virtual ~MonotoneMGStep() {}
virtual void setProblem(const MatrixType& mat,
VectorType& x,
VectorType& rhs)
{
MultigridStep<MatrixType, VectorType>::setProblem(mat,x,rhs);
oldCritical.resize(x.size(), false);
}
virtual void iterate();
virtual void preprocess();
virtual void nestedIteration();
ObstacleRestrictor<VectorType>* obstacleRestrictor_;
//! Bitfield determining which fine grid nodes have an obstacle
Dune::BitSetVector<dim>* hasObstacle_;
//! Vector containing the obstacle values of the fine grid nodes
ObstacleVectorType * obstacles_;
// Needed to track changes in the set of critical bits, and allows
// to check which dofs where critical after the last iteration
Dune::BitSetVector<dim> oldCritical;
protected:
//! Bitfield hierarchy containing the coarse obstacle nodes
std::vector<Dune::BitSetVector<dim>* > hasObstacleHierarchy_;
//! Hierarchy containing the obstacle values of the coarse obstacles
std::vector<ObstacleVectorType*> obstacleHierarchy_;
};
#include "mmgstep.cc"
#endif