Skip to content
Snippets Groups Projects
Commit 3e8d601e authored by Christoph Grüninger's avatar Christoph Grüninger
Browse files

Add UMFPack to test stokes2cni.

It is 10% faster in the test.
We have to investigate to use UMFPack more.
(reviewed by bernd)


git-svn-id: svn://svn.iws.uni-stuttgart.de/DUMUX/dumux/trunk@13068 2fb0f335-1f38-0410-981e-8018bf24f1b0
parent 9aef6f92
No related branches found
No related tags found
No related merge requests found
...@@ -26,6 +26,7 @@ ...@@ -26,6 +26,7 @@
#include <dune/istl/preconditioners.hh> #include <dune/istl/preconditioners.hh>
#include <dune/istl/solvers.hh> #include <dune/istl/solvers.hh>
#include <dune/istl/superlu.hh> #include <dune/istl/superlu.hh>
#include <dune/istl/umfpack.hh>
#include <dumux/common/parameters.hh> #include <dumux/common/parameters.hh>
#include <dumux/common/basicproperties.hh> #include <dumux/common/basicproperties.hh>
...@@ -587,7 +588,63 @@ private: ...@@ -587,7 +588,63 @@ private:
Dune::InverseOperatorResult result_; Dune::InverseOperatorResult result_;
const Problem& problem_; const Problem& problem_;
}; };
#endif #endif // HAVE_SUPERLU
#if HAVE_UMFPACK
template <class TypeTag>
class UMFPackBackend
{
typedef typename GET_PROP_TYPE(TypeTag, Problem) Problem;
public:
UMFPackBackend(const Problem& problem)
: problem_(problem)
{}
template<class Matrix, class Vector>
bool solve(const Matrix& A, Vector& x, const Vector& b)
{
Vector bTmp(b);
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
enum {blockSize = GET_PROP_VALUE(TypeTag, LinearSolverBlockSize)};
typedef typename Dune::FieldMatrix<Scalar, blockSize, blockSize> MatrixBlock;
typedef typename Dune::BCRSMatrix<MatrixBlock> ISTLMatrix;
int verbosity = GET_PARAM_FROM_GROUP(TypeTag, int, LinearSolver, Verbosity);
Dune::UMFPack<ISTLMatrix> solver(A, verbosity > 0);
solver.apply(x, bTmp, result_);
int size = x.size();
for (int i = 0; i < size; i++)
{
for (int j = 0; j < blockSize; j++)
{
if (std::isnan(x[i][j]) || std::isinf(x[i][j]))
{
result_.converged = false;
break;
}
}
}
return result_.converged;
}
const Dune::InverseOperatorResult& result() const
{
return result_;
}
private:
Dune::InverseOperatorResult result_;
const Problem& problem_;
};
#endif // HAVE_UMFPACK
} }
#endif #endif
...@@ -7,3 +7,4 @@ add_dumux_test(test_stokes2cni test_stokes2cni test_stokes2cni.cc ...@@ -7,3 +7,4 @@ add_dumux_test(test_stokes2cni test_stokes2cni test_stokes2cni.cc
-ParameterFile ${CMAKE_CURRENT_SOURCE_DIR}/test_stokes2cni.input -ParameterFile ${CMAKE_CURRENT_SOURCE_DIR}/test_stokes2cni.input
-Grid.File ${CMAKE_CURRENT_SOURCE_DIR}/grids/test_stokes2cni.dgf) -Grid.File ${CMAKE_CURRENT_SOURCE_DIR}/grids/test_stokes2cni.dgf)
add_dune_superlu_flags(test_stokes2cni) add_dune_superlu_flags(test_stokes2cni)
add_dune_umfpack_flags(test_stokes2cni)
...@@ -28,6 +28,7 @@ ...@@ -28,6 +28,7 @@
#if HAVE_PARDISO #if HAVE_PARDISO
#include <dumux/linear/pardisobackend.hh> #include <dumux/linear/pardisobackend.hh>
#endif #endif
#include <dumux/linear/seqsolverbackend.hh>
#include <dumux/material/fluidsystems/h2oairfluidsystem.hh> #include <dumux/material/fluidsystems/h2oairfluidsystem.hh>
#include <dumux/freeflow/stokesncni/stokesncnimodel.hh> #include <dumux/freeflow/stokesncni/stokesncnimodel.hh>
...@@ -58,6 +59,10 @@ SET_PROP(Stokes2cniTestProblem, FluidSystem) ...@@ -58,6 +59,10 @@ SET_PROP(Stokes2cniTestProblem, FluidSystem)
typedef Dumux::FluidSystems::H2OAir<Scalar> type; typedef Dumux::FluidSystems::H2OAir<Scalar> type;
}; };
#if HAVE_UMFPACK
SET_TYPE_PROP(Stokes2cniTestProblem, LinearSolver, UMFPackBackend<TypeTag>);
#endif
//! Scalar is set to type long double for higher accuracy //! Scalar is set to type long double for higher accuracy
//SET_TYPE_PROP(Stokes2cniTestProblem, Scalar, long double); //SET_TYPE_PROP(Stokes2cniTestProblem, Scalar, long double);
......
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