Skip to content
Snippets Groups Projects
Commit 17be7684 authored by Markus Blatt's avatar Markus Blatt
Browse files

Merge branch 'feature/fix-fastamg-buildhierarchy' into 'master'

[bugfix] Fixes access to non-existing parallel info in FastAMG::recalculateHierachy

recalculateHierarchy() assumed that there is information avaiable
that is only there in parallel runs. This failed miserably.

This commit fixes this be implementing dummy versions for
redistributing matrices in  a sequential run and does not
query any parallel stuff any more

This closes #11.

See merge request !34
parents 969fa8a4 3fa6173f
No related branches found
No related tags found
1 merge request!34[bugfix] Fixes access to non-existing parallel info in FastAMG::recalculateHierachy
......@@ -7,6 +7,7 @@
#include <dune/common/parallel/indexset.hh>
#include <dune/common/unused.hh>
#include <dune/istl/owneroverlapcopy.hh>
#include <dune/istl/paamg/pinfo.hh>
/**
* @file
* @brief Functionality for redistributing a sparse matrix.
......@@ -844,5 +845,22 @@ namespace Dune
redistributeMatrixEntries(origMatrix, newMatrix, origComm, newComm, ri);
}
#endif
template<typename M>
void redistributeMatrixEntries(M& origMatrix, M& newMatrix,
Dune::Amg::SequentialInformation& origComm,
Dune::Amg::SequentialInformation& newComm,
RedistributeInformation<Dune::Amg::SequentialInformation>& ri)
{
DUNE_THROW(InvalidStateException, "Trying to redistribute in sequential program!");
}
template<typename M>
void redistributeMatrix(M& origMatrix, M& newMatrix,
Dune::Amg::SequentialInformation& origComm,
Dune::Amg::SequentialInformation& newComm,
RedistributeInformation<Dune::Amg::SequentialInformation>& ri)
{
DUNE_THROW(InvalidStateException, "Trying to redistribute in sequential program!");
}
}
#endif
......@@ -1192,7 +1192,7 @@ namespace Dune
typename RedistributeInfoList::iterator riIter = redistributes_.begin();
Iterator level = matrices_.finest(), coarsest=matrices_.coarsest();
if(level.isRedistributed()) {
info->buildGlobalLookup(info->indexSet().size());
info->buildGlobalLookup(level->getmat().N());
redistributeMatrixEntries(const_cast<Matrix&>(level->getmat()),
const_cast<Matrix&>(level.getRedistributed().getmat()),
*info,info.getRedistributed(), *riIter);
......@@ -1206,7 +1206,7 @@ namespace Dune
++riIter;
productBuilder.calculate(fine, *(*amap), const_cast<Matrix&>(level->getmat()), *info, copyFlags);
if(level.isRedistributed()) {
info->buildGlobalLookup(info->indexSet().size());
info->buildGlobalLookup(level->getmat().N());
redistributeMatrixEntries(const_cast<Matrix&>(level->getmat()),
const_cast<Matrix&>(level.getRedistributed().getmat()), *info,
info.getRedistributed(), *riIter);
......
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