5#ifndef DUNE_ISTL_SUPERLU_HH
6#define DUNE_ISTL_SUPERLU_HH
18#include <dune/common/fmatrix.hh>
19#include <dune/common/fvector.hh>
20#include <dune/common/stdstreams.hh>
37 template<
class M,
class T,
class TM,
class TD,
class TA>
38 class SeqOverlappingSchwarz;
40 template<
class T,
bool tag>
41 struct SeqOverlappingSchwarzAssemblerHelper;
59#if __has_include("slu_sdefs.h")
63 static void create(SuperMatrix *
mat,
int n,
int m,
float *dat,
int n1,
64 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
66 sCreate_Dense_Matrix(
mat, n, m, dat, n1, stype, dtype, mtype);
70 static void destroy(SuperMatrix*)
75 struct SuperLUSolveChooser<float>
77 static void solve(superlu_options_t *options, SuperMatrix *
mat,
int *perm_c,
int *perm_r,
int *etree,
78 char *equed,
float *R,
float *C, SuperMatrix *L, SuperMatrix *U,
79 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
80 float *rpg,
float *rcond,
float *ferr,
float *berr,
81 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
84 sgssvx(options,
mat, perm_c, perm_r, etree, equed, R, C,
85 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
86 &gLU, memusage, stat, info);
91 struct QuerySpaceChooser<float>
93 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
95 sQuerySpace(L,U,memusage);
101#if __has_include("slu_ddefs.h")
104 struct SuperLUDenseMatChooser<double>
106 static void create(SuperMatrix *
mat,
int n,
int m,
double *dat,
int n1,
107 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
109 dCreate_Dense_Matrix(
mat, n, m, dat, n1, stype, dtype, mtype);
113 static void destroy(SuperMatrix * )
117 struct SuperLUSolveChooser<double>
119 static void solve(superlu_options_t *options, SuperMatrix *
mat,
int *perm_c,
int *perm_r,
int *etree,
120 char *equed,
double *R,
double *C, SuperMatrix *L, SuperMatrix *U,
121 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
122 double *rpg,
double *rcond,
double *ferr,
double *berr,
123 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
126 dgssvx(options,
mat, perm_c, perm_r, etree, equed, R, C,
127 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
128 &gLU, memusage, stat, info);
133 struct QuerySpaceChooser<double>
135 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
137 dQuerySpace(L,U,memusage);
142#if __has_include("slu_zdefs.h")
144 struct SuperLUDenseMatChooser<
std::complex<double> >
146 static void create(SuperMatrix *
mat,
int n,
int m, std::complex<double> *dat,
int n1,
147 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
149 zCreate_Dense_Matrix(
mat, n, m,
reinterpret_cast<doublecomplex*
>(dat), n1, stype, dtype, mtype);
153 static void destroy(SuperMatrix*)
158 struct SuperLUSolveChooser<
std::complex<double> >
160 static void solve(superlu_options_t *options, SuperMatrix *
mat,
int *perm_c,
int *perm_r,
int *etree,
161 char *equed,
double *R,
double *C, SuperMatrix *L, SuperMatrix *U,
162 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
163 double *rpg,
double *rcond,
double *ferr,
double *berr,
164 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
167 zgssvx(options,
mat, perm_c, perm_r, etree, equed, R, C,
168 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
169 &gLU, memusage, stat, info);
174 struct QuerySpaceChooser<
std::complex<double> >
176 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
178 zQuerySpace(L,U,memusage);
183#if __has_include("slu_cdefs.h")
185 struct SuperLUDenseMatChooser<
std::complex<float> >
187 static void create(SuperMatrix *
mat,
int n,
int m, std::complex<float> *dat,
int n1,
188 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
190 cCreate_Dense_Matrix(
mat, n, m,
reinterpret_cast< ::complex*
>(dat), n1, stype, dtype, mtype);
194 static void destroy(SuperMatrix* )
199 struct SuperLUSolveChooser<
std::complex<float> >
201 static void solve(superlu_options_t *options, SuperMatrix *
mat,
int *perm_c,
int *perm_r,
int *etree,
202 char *equed,
float *R,
float *C, SuperMatrix *L, SuperMatrix *U,
203 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
204 float *rpg,
float *rcond,
float *ferr,
float *berr,
205 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
208 cgssvx(options,
mat, perm_c, perm_r, etree, equed, R, C,
209 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
210 &gLU, memusage, stat, info);
215 struct QuerySpaceChooser<
std::complex<float> >
217 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
219 cQuerySpace(L,U,memusage);
227 struct SuperLUVectorChooser
230 template<
typename T,
typename A,
int n,
int m>
231 struct SuperLUVectorChooser<BCRSMatrix<FieldMatrix<T,n,m>,A > >
234 using domain_type = BlockVector<
236 typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > >;
238 using range_type = BlockVector<
240 typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > >;
243 template<
typename T,
typename A>
244 struct SuperLUVectorChooser<BCRSMatrix<T,A> >
247 using domain_type = BlockVector<T, A>;
249 using range_type = BlockVector<T, A>;
268 :
public InverseOperator<
269 typename Impl::SuperLUVectorChooser<M>::domain_type,
270 typename Impl::SuperLUVectorChooser<M>::range_type >
272 using T =
typename M::field_type;
282 using domain_type =
typename Impl::SuperLUVectorChooser<M>::domain_type;
284 using range_type =
typename Impl::SuperLUVectorChooser<M>::range_type;
307 bool reusevector=
true);
321 :
SuperLU(
mat, config.
get<bool>(
"verbose", false), config.
get<bool>(
"reuseVector", true))
350 void apply(T* x, T* b);
355 typename SuperLUMatrix::size_type
nnz()
const
357 return mat.nonzeroes();
371 const char*
name() {
return "SuperLU"; }
373 template<
class Mat,
class X,
class TM,
class TD,
class T1>
383 SuperMatrix L, U, B, X;
384 int *perm_c, *perm_r, *etree;
387 superlu_options_t options;
391 bool first, verbose, reusevector;
411 Destroy_SuperNode_Matrix(&L);
412 Destroy_CompCol_Matrix(&U);
415 if(!first && reusevector) {
416 SUPERLU_FREE(B.Store);
417 SUPERLU_FREE(X.Store);
425 : work(0), lwork(0), first(true), verbose(verbose_),
426 reusevector(reusevector_)
433 : work(0), lwork(0),verbose(false),
466 mat.setMatrix(mat_,mrs);
475 perm_c =
new int[
mat.
M()];
476 perm_r =
new int[
mat.
N()];
477 etree =
new int[
mat.
M()];
481 set_default_options(&options);
488 fakeFormat.lda=
mat.
N();
498 mem_usage_t memusage;
503 &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr,
504 &berr, &memusage, &stat, &info);
507 dinfo<<
"LU factorization: dgssvx() returns info "<< info<<std::endl;
509 auto nSuperLUCol =
static_cast<SuperMatrix&
>(
mat).ncol;
511 if ( info == 0 || info == nSuperLUCol+1 ) {
513 if ( options.PivotGrowth )
514 dinfo<<
"Recip. pivot growth = "<<rpg<<std::endl;
515 if ( options.ConditionNumber )
516 dinfo<<
"Recip. condition number = %e\n"<< rcond<<std::endl;
517 SCformat* Lstore = (SCformat *) L.Store;
518 NCformat* Ustore = (NCformat *) U.Store;
519 dinfo<<
"No of nonzeros in factor L = "<< Lstore->nnz<<std::endl;
520 dinfo<<
"No of nonzeros in factor U = "<< Ustore->nnz<<std::endl;
521 dinfo<<
"No of nonzeros in L+U = "<< Lstore->nnz + Ustore->nnz - nSuperLUCol<<std::endl;
523 dinfo<<
"L\\U MB "<<memusage.for_lu/1e6<<
" \ttotal MB needed "<<memusage.total_needed/1e6
525 std::cout<<stat.expansions<<std::endl;
527 }
else if ( info > 0 && lwork == -1 ) {
528 dinfo<<
"** Estimated memory: "<< info - nSuperLUCol<<std::endl;
530 if ( options.PrintStat ) StatPrint(&stat);
558 options.Fact = FACTORED;
565 if (
mat.
N() != b.dim())
566 DUNE_THROW(
ISTLError,
"Size of right-hand-side vector b does not match the number of matrix rows!");
567 if (
mat.
M() != x.dim())
568 DUNE_THROW(
ISTLError,
"Size of solution vector x does not match the number of matrix columns!");
570 DUNE_THROW(
ISTLError,
"Matrix of SuperLU is null!");
572 SuperMatrix* mB = &B;
573 SuperMatrix* mX = &X;
581 ((DNformat*)B.Store)->nzval=&b[0];
582 ((DNformat*)X.Store)->nzval=&x[0];
592 mem_usage_t memusage;
602 options.IterRefine=SLU_DOUBLE;
605 &L, &U, work, lwork, mB, mX, &rpg, &rcond, &ferr, &berr,
606 &memusage, &stat, &info);
626 dinfo<<
"Triangular solve: dgssvx() returns info "<< info<<std::endl;
628 auto nSuperLUCol =
static_cast<SuperMatrix&
>(
mat).ncol;
630 if ( info == 0 || info == nSuperLUCol+1 ) {
632 if ( options.IterRefine ) {
633 std::cout<<
"Iterative Refinement: steps="
634 <<stat.RefineSteps<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
636 std::cout<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
637 }
else if ( info > 0 && lwork == -1 ) {
638 std::cout<<
"** Estimated memory: "<< info - nSuperLUCol<<
" bytes"<<std::endl;
641 if ( options.PrintStat ) StatPrint(&stat);
645 SUPERLU_FREE(rB.Store);
646 SUPERLU_FREE(rX.Store);
655 DUNE_THROW(
ISTLError,
"Matrix of SuperLU is null!");
657 SuperMatrix* mB = &B;
658 SuperMatrix* mX = &X;
666 ((DNformat*) B.Store)->nzval=b;
667 ((DNformat*)X.Store)->nzval=x;
678 mem_usage_t memusage;
683 options.IterRefine=SLU_DOUBLE;
686 &L, &U, work, lwork, mB, mX, &rpg, &rcond, &ferr, &berr,
687 &memusage, &stat, &info);
690 dinfo<<
"Triangular solve: dgssvx() returns info "<< info<<std::endl;
692 auto nSuperLUCol =
static_cast<SuperMatrix&
>(
mat).ncol;
694 if ( info == 0 || info == nSuperLUCol+1 ) {
696 if ( options.IterRefine ) {
697 dinfo<<
"Iterative Refinement: steps="
698 <<stat.RefineSteps<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
700 dinfo<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
701 }
else if ( info > 0 && lwork == -1 ) {
702 dinfo<<
"** Estimated memory: "<< info - nSuperLUCol<<
" bytes"<<std::endl;
704 if ( options.PrintStat ) StatPrint(&stat);
709 SUPERLU_FREE(rB.Store);
710 SUPERLU_FREE(rX.Store);
715 template<
typename T,
typename A>
721 template<
typename T,
typename A>
730 template<
int k>
struct isValidBlock<
Dune::FieldVector<std::complex<double>,k>> : std::true_type{};
731 template<
typename TL,
typename M>
732 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
733 typename Dune::TypeListElement<2, TL>::type>>
735 std::enable_if_t<
isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
737 int verbose = config.get(
"verbose", 0);
738 return std::make_shared<Dune::SuperLU<M>>(
mat,verbose);
742 template<
typename TL,
typename M>
743 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
744 typename Dune::TypeListElement<2, TL>::type>>
746 std::enable_if_t<!
isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
749 "Unsupported Type in SuperLU (only double and std::complex<double> supported)");
759#undef FIRSTCOL_OF_SNODE
Implementation of the BCRSMatrix class.
This file implements a vector space as a tensor product of a given vector space. The number of compon...
Implementations of the inverse operator interface.
Templates characterizing the type of a solver.
#define DUNE_REGISTER_DIRECT_SOLVER(name,...)
Definition solverregistry.hh:13
void setSubMatrix(const Matrix &mat, const S &rowIndexSet)
Definition superlu.hh:457
void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition superlu.hh:563
void setVerbosity(bool v)
Definition superlu.hh:437
void free()
free allocated space.
Definition superlu.hh:403
~SuperLU()
Definition superlu.hh:396
SuperLU()
Empty default constructor.
Definition superlu.hh:432
void setMatrix(const Matrix &mat)
Initialize data from given matrix.
Definition superlu.hh:443
Matrix & mat
Definition matrixmatrix.hh:347
Definition allocator.hh:11
PropertyMapTypeSelector< Amg::VertexVisitedTag, Amg::PropertiesGraph< G, Amg::VertexProperties, EP, VM, EM > >::Type get(const Amg::VertexVisitedTag &tag, Amg::PropertiesGraph< G, Amg::VertexProperties, EP, VM, EM > &graph)
Definition dependency.hh:293
A sparse block matrix with compressed row storage.
Definition matrixutils.hh:24
derive error class from the base class in common
Definition istlexception.hh:19
Sequential overlapping Schwarz preconditioner.
Definition umfpack.hh:43
A generic dynamic dense matrix.
Definition matrixutils.hh:30
size_type M() const
Return the number of columns.
Definition matrix.hh:696
size_type N() const
Return the number of rows.
Definition matrix.hh:691
Statistics about the application of an inverse operator.
Definition solver.hh:50
int iterations
Number of iterations.
Definition solver.hh:69
bool converged
True if convergence criterion has been met.
Definition solver.hh:75
Category
Definition solvercategory.hh:23
@ sequential
Category for sequential solvers.
Definition solvercategory.hh:25
Definition solverregistry.hh:77
Definition solvertype.hh:16
@ value
Whether this is a direct solver.
Definition solvertype.hh:24
Definition solvertype.hh:30
@ value
whether the solver internally uses column compressed storage
Definition solvertype.hh:36
SuperLu Solver.
Definition supermatrix.hh:182
Dune::SuperLUMatrix< Matrix > SuperLUMatrix
The corresponding SuperLU Matrix type.
Definition superlu.hh:278
typename Impl::SuperLUVectorChooser< M >::domain_type domain_type
The type of the domain of the solver.
Definition superlu.hh:282
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition superlu.hh:287
SuperLUMatrix::size_type nnz() const
Definition superlu.hh:355
M matrix_type
Definition superlu.hh:276
SuperLU(const Matrix &mat, const ParameterTree &config)
Constructs the SuperLU solver.
Definition superlu.hh:320
const char * name()
Definition superlu.hh:371
M Matrix
The matrix type.
Definition superlu.hh:275
typename Impl::SuperLUVectorChooser< M >::range_type range_type
The type of the range of the solver.
Definition superlu.hh:284
SuperMatrixInitializer< Matrix > MatrixInitializer
Type of an associated initializer class.
Definition superlu.hh:280
void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition superlu.hh:342
Definition superlu.hh:727
std::shared_ptr< Dune::InverseOperator< typename Dune::TypeListElement< 1, TL >::type, typename Dune::TypeListElement< 2, TL >::type > > operator()(TL, const M &mat, const Dune::ParameterTree &config, std::enable_if_t< isValidBlock< typename Dune::TypeListElement< 1, TL >::type::block_type >::value, int >=0) const
Definition superlu.hh:734
Definition superlu.hh:728
Definition supermatrix.hh:132
Utility class for converting an ISTL Matrix into a SuperLU Matrix.
Definition supermatrix.hh:175
Definition supermatrix.hh:179