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

Added a template meta program, that checks whether a solver is a

direct solver.
Added a method to AMG to check whether it uses a direct solver on the
coarse level.

[[Imported from SVN: r1547]]
parent f26d64b6
No related branches found
No related tags found
No related merge requests found
......@@ -12,6 +12,7 @@
#include <dune/istl/solvers.hh>
#include <dune/istl/scalarproducts.hh>
#include <dune/istl/superlu.hh>
#include <dune/istl/solvertype.hh>
#include <dune/common/typetraits.hh>
#include <dune/common/exceptions.hh>
......@@ -204,6 +205,12 @@ namespace Dune
matrices_->recalculateGalerkin(NegateSet<typename PI::OwnerSet>());
}
/**
* @brief Check whether the coarse solver used is a direct solver.
* @return True if the coarse level solver is a direct solver.
*/
bool usesDirectCoarseLevelSolver() const;
private:
/** @brief Multigrid cycle on a level. */
void mgc();
......@@ -653,6 +660,12 @@ namespace Dune
}
template<class M, class X, class S, class PI, class A>
bool AMG<M,X,S,PI,A>::usesDirectCoarseLevelSolver() const
{
return IsDirectSolver< CoarseSolver>::value;
}
template<class M, class X, class S, class PI, class A>
void AMG<M,X,S,PI,A>::mgc(){
if(matrix == matrices_->matrices().coarsest() && levels()==maxlevels()) {
......
......@@ -4,7 +4,7 @@
#define DUNE_PARDISO_HH
#include <dune/istl/preconditioners.hh>
#include <dune/istl/solvertype.hh>
/* Change this, if your Fortran compiler does not append underscores. */
/* e.g. the AIX compiler: #define F77_FUNC(func) func */
......@@ -209,11 +209,12 @@ namespace Dune {
int num_procs_; //!< number of processors.
};
}
template<class M, class X, class Y>
struct IsDirectSolver<SeqPardiso<M,X,Y> >
{
enum { value=true};
};
} // end namespace Dune
#endif
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_ISTL_SOLVERTYPE_HH
#define DUNE_ISTL_SOLVERTYPE_HH
/**
* @file
* @brief Templates characterizing the type of a solver.
*/
namespace Dune
{
template<typename Solver>
struct IsDirectSolver
{
enum
{
/**
* @brief Whether this is a direct solver.
*
* If Solver is a direct solver, this is true.
*/
value =false
};
};
} // end namespace Dune
#endif
......@@ -26,6 +26,7 @@
#include <dune/common/fmatrix.hh>
#include <dune/common/fvector.hh>
#include <dune/common/stdstreams.hh>
#include <dune/istl/solvertype.hh>
namespace Dune
{
......@@ -447,6 +448,12 @@ namespace Dune
StatFree(&stat);
}
/** @} */
template<typename T, typename A, int n, int m>
struct IsDirectSolver<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
{
enum { value=true};
};
}
#endif // HAVE_SUPERLU
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment