Skip to content
Snippets Groups Projects
Select Git revision
  • 5f89f48ca7ba22c22a81ac161e2e8cfedab13e02
  • master default
  • feature/umfpack-multitype
  • proximal-gradient-solver
  • feature/pn-solver
  • fix-master
  • feature/codespell
  • codespell
  • feature/cholmod-reuse-factorization
  • releases/2.8
  • test
  • fix/loopsolver-criterions
  • feature/multitype-cholmod
  • feature/two-norm
  • feature/replace-unused
  • feature/generic-transfer-operators
  • feature/cholmod-solver
  • releases/2.7-1
  • feature/P0-element-tranferoperators
  • flexible-loopsolver-max
  • releases/2.5-1
  • subversion->git
22 results

mmgstep.hh

Blame
  • Forked from agnumpde / dune-solvers
    Source project has a limited visibility.
    Code owners
    Assign users and groups as approvers for specific file changes. Learn more.
    mmgstep.hh 4.29 KiB
    // -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 4 -*-
    // vi: set et ts=8 sw=4 sts=4:
    #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/common/wrapownshare.hh>
    #include <dune/solvers/transferoperators/obstaclerestrictor.hh>
    
    /** \brief The general monotone multigrid solver
    */
    template<class MatrixType, class VectorType>
    class MonotoneMGStep : public MultigridStep<MatrixType, VectorType>
    {
    
        using Base = 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)
            : Base(mat, x, rhs),
              hasObstacleHierarchy_(0), obstacleHierarchy_(0)
        {
            oldCritical.resize(x.size(), false);
        }
    
        virtual ~MonotoneMGStep() {}
    
        virtual void setProblem(const MatrixType& mat,
                                VectorType& x,
                                const VectorType& rhs)
        {
            Base::setProblem(mat,x,rhs);
            oldCritical.resize(x.size(), false);
            oldCritical.unsetAll();
        }
    
        virtual void iterate();
    
        virtual void preprocess();
    
        virtual void nestedIteration();
    
        //! Set the hasObstacle bitfield
        [[deprecated("Setting by raw pointer is deprecated. Use l-value or r-value or a shared_ptr")]]
        void setHasObstacles(Dune::BitSetVector<dim>* hasObstacle) {
          hasObstacle_ = Dune::stackobject_to_shared_ptr(*hasObstacle);
        }
    
        //! Set the hasObstacle bitfield
        template <class BitVector, class Enable = std::enable_if_t<not std::is_pointer<BitVector>::value> >
        void setHasObstacles(BitVector&& hasObstacles) {
          hasObstacle_ = Dune::Solvers::wrap_own_share<Dune::BitSetVector<dim> >(std::forward<BitVector>(hasObstacles));
        }
    
        //! Set the obstacle field
        [[deprecated("Setting by raw pointer is deprecated. Use l-value or r-value or a shared_ptr")]]
        void setObstacles(ObstacleVectorType* obstacles) {
          obstacles_ = Dune::stackobject_to_shared_ptr(*obstacles);
        }
    
        //! Set the obstacle field
        template <class ObstacleVector, class Enable = std::enable_if_t<not std::is_pointer<ObstacleVector>::value> >
        void setObstacles(ObstacleVector&& obstacles) {
          obstacles_ = Dune::Solvers::wrap_own_share<ObstacleVectorType>(std::forward<ObstacleVector>(obstacles));
        }
    
        //! Get the obstacle field
        ObstacleVectorType& getObstacles() { return *obstacles_; }
    
        //! Set the obstacle restrictor
        template <class Restrictor>
        void setObstacleRestrictor(Restrictor&& restrictor) {
            obstacleRestrictor_ = Dune::Solvers::wrap_own_share<ObstacleRestrictor<VectorType> >(std::forward<Restrictor>(restrictor));
        }
    
    protected:
        //! Bitfield determining which fine grid nodes have an obstacle
        std::shared_ptr<Dune::BitSetVector<dim> > hasObstacle_;
    
        //! Vector containing the obstacle values of the fine grid nodes
        std::shared_ptr<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;
    
        //! The restrictor used to construct the coarse obstacles
        std::shared_ptr<ObstacleRestrictor<VectorType> > obstacleRestrictor_;
    
        //! Bitfield hierarchy containing the coarse obstacle nodes
        std::vector<std::shared_ptr<Dune::BitSetVector<dim>> > hasObstacleHierarchy_;
    
        //! For each level specify the matrix lines that need reassembly
        std::vector<Dune::BitSetVector<dim> > recompute_;
    
        //! For each level specify dofs that changed status
        std::vector<Dune::BitSetVector<dim> > oldCritical_;
    
        //! Hierarchy containing the obstacle values of the coarse obstacles
        std::vector<std::shared_ptr<ObstacleVectorType> >  obstacleHierarchy_;
    private:
        using Base::setProblem;
    
    };
    
    #include "mmgstep.cc"
    
    #endif