3#ifndef DUNE_ISTL_SUPERLU_HH
4#define DUNE_ISTL_SUPERLU_HH
14#define SUPERLU_NTYPE 1
16#ifdef SUPERLU_POST_2005_VERSION
42#warning Support for SuperLU older than SuperLU 3.0 from August 2005 is deprecated.
55#include "supermatrix.hh"
60#include "istlexception.hh"
79 template<
class Matrix>
83 template<
class M,
class T,
class TM,
class TD,
class TA>
84 class SeqOverlappingSchwarz;
86 template<
class T,
bool tag>
87 struct SeqOverlappingSchwarzAssemblerHelper;
90 struct SuperLUSolveChooser
94 struct SuperLUDenseMatChooser
98 struct SuperLUQueryChooser
102 struct QuerySpaceChooser
107 struct SuperLUDenseMatChooser<float>
109 static void create(SuperMatrix *mat,
int n,
int m,
float *dat,
int n1,
110 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
112 sCreate_Dense_Matrix(mat, n, m, dat, n1, stype, dtype, mtype);
116 static void destroy(SuperMatrix*)
121 struct SuperLUSolveChooser<float>
123 static void solve(superlu_options_t *options, SuperMatrix *mat,
int *perm_c,
int *perm_r,
int *etree,
124 char *equed,
float *R,
float *C, SuperMatrix *L, SuperMatrix *U,
125 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
126 float *rpg,
float *rcond,
float *ferr,
float *berr,
127 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
129#if SUPERLU_MIN_VERSION_5
131 sgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
132 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
133 &gLU, memusage, stat, info);
135 sgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
136 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
137 memusage, stat, info);
143 struct QuerySpaceChooser<float>
145 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
147 sQuerySpace(L,U,memusage);
156 struct SuperLUDenseMatChooser<double>
158 static void create(SuperMatrix *mat,
int n,
int m,
double *dat,
int n1,
159 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
161 dCreate_Dense_Matrix(mat, n, m, dat, n1, stype, dtype, mtype);
165 static void destroy(SuperMatrix * )
169 struct SuperLUSolveChooser<double>
171 static void solve(superlu_options_t *options, SuperMatrix *mat,
int *perm_c,
int *perm_r,
int *etree,
172 char *equed,
double *R,
double *C, SuperMatrix *L, SuperMatrix *U,
173 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
174 double *rpg,
double *rcond,
double *ferr,
double *berr,
175 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
177#if SUPERLU_MIN_VERSION_5
179 dgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
180 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
181 &gLU, memusage, stat, info);
183 dgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
184 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
185 memusage, stat, info);
191 struct QuerySpaceChooser<double>
193 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
195 dQuerySpace(L,U,memusage);
202 struct SuperLUDenseMatChooser<
std::complex<double> >
204 static void create(SuperMatrix *mat,
int n,
int m, std::complex<double> *dat,
int n1,
205 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
207 zCreate_Dense_Matrix(mat, n, m,
reinterpret_cast<doublecomplex*
>(dat), n1, stype, dtype, mtype);
211 static void destroy(SuperMatrix*)
216 struct SuperLUSolveChooser<
std::complex<double> >
218 static void solve(superlu_options_t *options, SuperMatrix *mat,
int *perm_c,
int *perm_r,
int *etree,
219 char *equed,
double *R,
double *C, SuperMatrix *L, SuperMatrix *U,
220 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
221 double *rpg,
double *rcond,
double *ferr,
double *berr,
222 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
224#if SUPERLU_MIN_VERSION_5
226 zgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
227 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
228 &gLU, memusage, stat, info);
230 zgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
231 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
232 memusage, stat, info);
238 struct QuerySpaceChooser<
std::complex<double> >
240 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
242 zQuerySpace(L,U,memusage);
249 struct SuperLUDenseMatChooser<
std::complex<float> >
251 static void create(SuperMatrix *mat,
int n,
int m, std::complex<float> *dat,
int n1,
252 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
254 cCreate_Dense_Matrix(mat, n, m,
reinterpret_cast< ::complex*
>(dat), n1, stype, dtype, mtype);
258 static void destroy(SuperMatrix* )
263 struct SuperLUSolveChooser<
std::complex<float> >
265 static void solve(superlu_options_t *options, SuperMatrix *mat,
int *perm_c,
int *perm_r,
int *etree,
266 char *equed,
float *R,
float *C, SuperMatrix *L, SuperMatrix *U,
267 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
268 float *rpg,
float *rcond,
float *ferr,
float *berr,
269 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
271#if SUPERLU_MIN_VERSION_5
273 cgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
274 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
275 &gLU, memusage, stat, info);
277 cgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
278 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
279 memusage, stat, info);
285 struct QuerySpaceChooser<
std::complex<float> >
287 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
289 cQuerySpace(L,U,memusage);
307 template<
typename T,
typename A,
int n,
int m>
310 BlockVector<FieldVector<T,m>,
311 typename A::template rebind<FieldVector<T,m> >::other>,
312 BlockVector<FieldVector<T,n>,
313 typename A::template rebind<FieldVector<T,n> >::other> >
326 typename A::template rebind<FieldVector<T,m> >::other>
domain_type;
330 typename A::template rebind<FieldVector<T,n> >::other>
range_type;
345 explicit SuperLU(
const Matrix& mat,
bool verbose=
false,
346 bool reusevector=
true);
374 void apply(T* x, T* b);
377 void setMatrix(
const Matrix& mat);
379 typename SuperLUMatrix::size_type nnz()
const
385 void setSubMatrix(
const Matrix& mat,
const S& rowIndexSet);
387 void setVerbosity(
bool v);
395 const char* name() {
return "SuperLU"; }
397 friend class std::mem_fun_ref_t<void,SuperLU>;
398 template<
class M,
class X,
class TM,
class TD,
class T1>
399 friend class SeqOverlappingSchwarz;
400 friend struct SeqOverlappingSchwarzAssemblerHelper<SuperLU<Matrix>,true>;
402 SuperLUMatrix& getInternalMatrix() {
return mat; }
408 SuperMatrix L, U, B, X;
409 int *perm_c, *perm_r, *etree;
410 typename GetSuperLUType<T>::float_type *R, *C;
412 superlu_options_t options;
416 bool first, verbose, reusevector;
419 template<
typename T,
typename A,
int n,
int m>
420 SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
423 if(mat.
N()+mat.
M()>0)
427 template<
typename T,
typename A,
int n,
int m>
428 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::free()
436 Destroy_SuperNode_Matrix(&L);
437 Destroy_CompCol_Matrix(&U);
440 if(!first && reusevector) {
441 SUPERLU_FREE(B.Store);
442 SUPERLU_FREE(X.Store);
447 template<
typename T,
typename A,
int n,
int m>
448 SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
449 ::SuperLU(
const Matrix& mat_,
bool verbose_,
bool reusevector_)
450 : work(0), lwork(0), first(true), verbose(verbose_),
451 reusevector(reusevector_)
456 template<
typename T,
typename A,
int n,
int m>
457 SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::SuperLU()
458 : work(0), lwork(0),verbose(false),
461 template<
typename T,
typename A,
int n,
int m>
462 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setVerbosity(
bool v)
467 template<
typename T,
typename A,
int n,
int m>
468 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setMatrix(
const Matrix& mat_)
470 if(mat.
N()+mat.
M()>0) {
480 template<
typename T,
typename A,
int n,
int m>
482 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setSubMatrix(
const Matrix& mat_,
485 if(mat.
N()+mat.
M()>0) {
491 mat.setMatrix(mat_,mrs);
495 template<
typename T,
typename A,
int n,
int m>
496 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::decompose()
500 perm_c =
new int[mat.
M()];
501 perm_r =
new int[mat.
N()];
502 etree =
new int[mat.
M()];
503 R =
new typename GetSuperLUType<T>::float_type[mat.
N()];
504 C =
new typename GetSuperLUType<T>::float_type[mat.
M()];
506 set_default_options(&options);
510 B.Dtype=GetSuperLUType<T>::type;
513 fakeFormat.lda=mat.
N();
516 X.Dtype=GetSuperLUType<T>::type;
521 typename GetSuperLUType<T>::float_type rpg, rcond, ferr=1e10, berr=1e10;
523 mem_usage_t memusage;
527 SuperLUSolveChooser<T>::solve(&options, &
static_cast<SuperMatrix&
>(mat), perm_c, perm_r, etree, &equed, R, C,
528 &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr,
529 &berr, &memusage, &stat, &info);
532 dinfo<<
"LU factorization: dgssvx() returns info "<< info<<std::endl;
534 if ( info == 0 || info == n+1 ) {
536 if ( options.PivotGrowth )
537 dinfo<<
"Recip. pivot growth = "<<rpg<<std::endl;
538 if ( options.ConditionNumber )
539 dinfo<<
"Recip. condition number = %e\n"<< rcond<<std::endl;
540 SCformat* Lstore = (SCformat *) L.Store;
541 NCformat* Ustore = (NCformat *) U.Store;
542 dinfo<<
"No of nonzeros in factor L = "<< Lstore->nnz<<std::endl;
543 dinfo<<
"No of nonzeros in factor U = "<< Ustore->nnz<<std::endl;
544 dinfo<<
"No of nonzeros in L+U = "<< Lstore->nnz + Ustore->nnz - n<<std::endl;
545 QuerySpaceChooser<T>::querySpace(&L, &U, &memusage);
546 dinfo<<
"L\\U MB "<<memusage.for_lu/1e6<<
" \ttotal MB needed "<<memusage.total_needed/1e6
549#ifdef HAVE_MEM_USAGE_T_EXPANSIONS
550 std::cout<<memusage.expansions<<std::endl;
552 std::cout<<stat.expansions<<std::endl;
554 }
else if ( info > 0 && lwork == -1 ) {
555 dinfo<<
"** Estimated memory: "<< info - n<<std::endl;
557 if ( options.PrintStat ) StatPrint(&stat);
585 options.Fact = FACTORED;
588 template<
typename T,
typename A,
int n,
int m>
589 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
592 if (mat.
N() != b.
dim())
593 DUNE_THROW(
ISTLError,
"Size of right-hand-side vector b does not match the number of matrix rows!");
594 if (mat.
M() != x.
dim())
595 DUNE_THROW(
ISTLError,
"Size of solution vector x does not match the number of matrix columns!");
596 if (mat.
M()+mat.
N()==0)
599 SuperMatrix* mB = &B;
600 SuperMatrix* mX = &X;
604 SuperLUDenseMatChooser<T>::create(&B, (
int)mat.
N(), 1,
reinterpret_cast<T*
>(&b[0]), (
int)mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
605 SuperLUDenseMatChooser<T>::create(&X, (
int)mat.
N(), 1,
reinterpret_cast<T*
>(&x[0]), (
int)mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
608 ((DNformat*)B.Store)->nzval=&b[0];
609 ((DNformat*)X.Store)->nzval=&x[0];
612 SuperLUDenseMatChooser<T>::create(&rB, (
int)mat.
N(), 1,
reinterpret_cast<T*
>(&b[0]), (
int)mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
613 SuperLUDenseMatChooser<T>::create(&rX, (
int)mat.
N(), 1,
reinterpret_cast<T*
>(&x[0]), (
int)mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
617 typename GetSuperLUType<T>::float_type rpg, rcond, ferr=1e10, berr;
619 mem_usage_t memusage;
629#ifdef SUPERLU_MIN_VERSION_4_3
630 options.IterRefine=SLU_DOUBLE;
632 options.IterRefine=DOUBLE;
635 SuperLUSolveChooser<T>::solve(&options, &
static_cast<SuperMatrix&
>(mat), perm_c, perm_r, etree, &equed, R, C,
636 &L, &U, work, lwork, mB, mX, &rpg, &rcond, &ferr, &berr,
637 &memusage, &stat, &info);
657 dinfo<<
"Triangular solve: dgssvx() returns info "<< info<<std::endl;
659 if ( info == 0 || info == n+1 ) {
661 if ( options.IterRefine ) {
662 std::cout<<
"Iterative Refinement: steps="
663 <<stat.RefineSteps<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
665 std::cout<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
666 }
else if ( info > 0 && lwork == -1 ) {
667 std::cout<<
"** Estimated memory: "<< info - n<<
" bytes"<<std::endl;
670 if ( options.PrintStat ) StatPrint(&stat);
674 SUPERLU_FREE(rB.Store);
675 SUPERLU_FREE(rX.Store);
679 template<
typename T,
typename A,
int n,
int m>
680 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
683 if(mat.
N()+mat.
M()==0)
686 SuperMatrix* mB = &B;
687 SuperMatrix* mX = &X;
691 SuperLUDenseMatChooser<T>::create(&B, mat.
N(), 1, b, mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
692 SuperLUDenseMatChooser<T>::create(&X, mat.
N(), 1, x, mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
695 ((DNformat*) B.Store)->nzval=b;
696 ((DNformat*)X.Store)->nzval=x;
699 SuperLUDenseMatChooser<T>::create(&rB, mat.
N(), 1, b, mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
700 SuperLUDenseMatChooser<T>::create(&rX, mat.
N(), 1, x, mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
705 typename GetSuperLUType<T>::float_type rpg, rcond, ferr=1e10, berr;
707 mem_usage_t memusage;
712#ifdef SUPERLU_MIN_VERSION_4_3
713 options.IterRefine=SLU_DOUBLE;
715 options.IterRefine=DOUBLE;
718 SuperLUSolveChooser<T>::solve(&options, &
static_cast<SuperMatrix&
>(mat), perm_c, perm_r, etree, &equed, R, C,
719 &L, &U, work, lwork, mB, mX, &rpg, &rcond, &ferr, &berr,
720 &memusage, &stat, &info);
723 dinfo<<
"Triangular solve: dgssvx() returns info "<< info<<std::endl;
725 if ( info == 0 || info == n+1 ) {
727 if ( options.IterRefine ) {
728 dinfo<<
"Iterative Refinement: steps="
729 <<stat.RefineSteps<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
731 dinfo<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
732 }
else if ( info > 0 && lwork == -1 ) {
733 dinfo<<
"** Estimated memory: "<< info - n<<
" bytes"<<std::endl;
735 if ( options.PrintStat ) StatPrint(&stat);
740 SUPERLU_FREE(rB.Store);
741 SUPERLU_FREE(rX.Store);
746 template<
typename T,
typename A,
int n,
int m>
752 template<
typename T,
typename A,
int n,
int m>
753 struct StoresColumnCompressed<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
755 enum { value =
true };
760#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...
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:413
A vector of blocks with memory management.
Definition: bvector.hh:253
A dense n x m matrix.
Definition: fmatrix.hh:67
vector space out of a tensor product of fields.
Definition: fvector.hh:94
derive error class from the base class in common
Definition: istlexception.hh:16
Abstract base class for all solvers.
Definition: solver.hh:79
A generic dynamic dense matrix.
Definition: matrix.hh:25
size_type M() const
Return the number of columns.
Definition: matrix.hh:165
size_type N() const
Return the number of rows.
Definition: matrix.hh:160
Dune::BlockVector< FieldVector< T, n >, typename A::template rebind< FieldVector< T, n > >::other > range_type
The type of the range of the solver.
Definition: superlu.hh:330
Dune::SuperLUMatrix< Matrix > SuperLUMatrix
The corresponding SuperLU Matrix type.
Definition: superlu.hh:320
Dune::BlockVector< FieldVector< T, m >, typename A::template rebind< FieldVector< T, m > >::other > domain_type
The type of the domain of the solver.
Definition: superlu.hh:326
void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: superlu.hh:365
SuperMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A > > MatrixInitializer
Type of an associated initializer class.
Definition: superlu.hh:322
size_type dim() const
dimension of the vector space
Definition: bvector.hh:223
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:243
DInfoType dinfo(std::cout)
Stream for informative output.
Definition: stdstreams.hh:138
Dune namespace.
Definition: alignment.hh:10
Implementations of the inverse operator interface.
Templates characterizing the type of a solver.
Standard Dune debug streams.
Statistics about the application of an inverse operator.
Definition: solver.hh:32
int iterations
Number of iterations.
Definition: solver.hh:50
bool converged
True if convergence criterion has been met.
Definition: solver.hh:56
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentional unused function parameters with.
Definition: unused.hh:18