Skip to content
Snippets Groups Projects
Commit d0a4d569 authored by Jonathan Youett's avatar Jonathan Youett
Browse files

Remove scalar version, it is not needed anymore and its presence leads to an...

Remove scalar version, it is not needed anymore and its presence leads to an compiler error if the grid dimension is equal to one.
parent c196b700
No related branches found
No related tags found
No related merge requests found
#include <dune/solvers/transferoperators/compressedmultigridtransfer.hh>
template <class DiscFuncType>
void MandelObstacleRestrictor<DiscFuncType>::
restrict(const std::vector<BoxConstraint<typename DiscFuncType::field_type,blocksize> >&f,
std::vector<BoxConstraint<typename DiscFuncType::field_type,blocksize> > &t,
const Dune::BitSetVector<1>& fHasObstacle,
const Dune::BitSetVector<1>& tHasObstacle,
const MultigridTransfer<DiscFuncType>& transfer,
const Dune::BitSetVector<blocksize>& critical)
{
if (!(dynamic_cast<const CompressedMultigridTransfer<DiscFuncType>*>(&transfer)))
DUNE_THROW(Dune::NotImplemented, "Implementation assumes transfer to be of type CompressedMultigridTransfer!");
const Dune::BCRSMatrix<Dune::FieldMatrix<field_type, 1, 1> >& transferMat =
dynamic_cast<const CompressedMultigridTransfer<DiscFuncType>*>(&transfer)->getMatrix();
t.resize(transferMat.M());
for (size_t i=0; i<t.size(); i++)
t[i].clear();
// The Mandel restriction computes the maximum $m$ of all nodal base coefficients
// in the support of a given coarse grid base function \f$ \phi_c \f$. If the
// coefficient of \f$ \phi_c \f$ is smaller than \f$ m \f$ , it is set to \f$ m \f$.
// We use the multigrid transfer matrix to find all relevant fine grid dofs
typedef typename Dune::BCRSMatrix<Dune::FieldMatrix<field_type,1,1> >::row_type RowType;
typedef typename RowType::ConstIterator ColumnIterator;
// Loop over all coarse grid dofs
for (size_t i=0; i<transferMat.N(); i++) {
if (fHasObstacle[i].none())
continue;
const RowType& row = transferMat[i];
ColumnIterator cIt = row.begin();
ColumnIterator cEndIt = row.end();
// Loop over all fine grid base functions in the support of tObsIt
for(; cIt!=cEndIt; ++cIt) {
if (tHasObstacle[cIt.index()].none())
continue;
// Sort out spurious entries due to numerical dirt
if ((*cIt)[0][0] < 1e-3)
continue;
for (int j=0; j<blocksize; j++)
if (!critical[i][j]) {
t[cIt.index()].lower(j) = std::max(t[cIt.index()].lower(j), f[i].lower(j));
t[cIt.index()].upper(j) = std::min(t[cIt.index()].upper(j), f[i].upper(j));
}
}
}
}
template <class DiscFuncType>
void MandelObstacleRestrictor<DiscFuncType>::
restrict(const std::vector<BoxConstraint<typename DiscFuncType::field_type,blocksize> >&f,
......
......@@ -15,13 +15,6 @@ class MandelObstacleRestrictor : public ObstacleRestrictor<DiscFuncType>
public:
virtual void restrict(const std::vector<BoxConstraint<field_type,blocksize> >& f,
std::vector<BoxConstraint<field_type,blocksize> >& t,
const Dune::BitSetVector<1>& fHasObstacle,
const Dune::BitSetVector<1>& tHasObstacle,
const MultigridTransfer<DiscFuncType>& transfer,
const Dune::BitSetVector<blocksize>& critical);
virtual void restrict(const std::vector<BoxConstraint<field_type,blocksize> >& f,
std::vector<BoxConstraint<field_type,blocksize> >& t,
const Dune::BitSetVector<blocksize>& fHasObstacle,
......
......@@ -18,14 +18,7 @@ class ObstacleRestrictor
typedef typename DiscFuncType::field_type field_type;
public:
virtual void restrict(const std::vector<BoxConstraint<field_type,blocksize> >& f,
std::vector<BoxConstraint<field_type,blocksize> >& t,
const Dune::BitSetVector<1>& fHasObstacle,
const Dune::BitSetVector<1>& tHasObstacle,
const MultigridTransfer<DiscFuncType>& transfer,
const Dune::BitSetVector<blocksize>& critical) = 0;
virtual void restrict(const std::vector<BoxConstraint<field_type,blocksize> >& f,
std::vector<BoxConstraint<field_type,blocksize> >& t,
const Dune::BitSetVector<blocksize>& fHasObstacle,
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment