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 *m)
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 sgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
130 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
131 memusage, stat, info);
136 struct QuerySpaceChooser<float>
138 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
140 sQuerySpace(L,U,memusage);
149 struct SuperLUDenseMatChooser<double>
151 static void create(SuperMatrix *mat,
int n,
int m,
double *dat,
int n1,
152 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
154 dCreate_Dense_Matrix(mat, n, m, dat, n1, stype, dtype, mtype);
158 static void destroy(SuperMatrix * )
162 struct SuperLUSolveChooser<double>
164 static void solve(superlu_options_t *options, SuperMatrix *mat,
int *perm_c,
int *perm_r,
int *etree,
165 char *equed,
double *R,
double *C, SuperMatrix *L, SuperMatrix *U,
166 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
167 double *rpg,
double *rcond,
double *ferr,
double *berr,
168 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
170 dgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
171 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
172 memusage, stat, info);
177 struct QuerySpaceChooser<double>
179 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
181 dQuerySpace(L,U,memusage);
188 struct SuperLUDenseMatChooser<
std::complex<double> >
190 static void create(SuperMatrix *mat,
int n,
int m, std::complex<double> *dat,
int n1,
191 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
193 zCreate_Dense_Matrix(mat, n, m,
reinterpret_cast<doublecomplex*
>(dat), n1, stype, dtype, mtype);
197 static void destroy(SuperMatrix *mat)
202 struct SuperLUSolveChooser<
std::complex<double> >
204 static void solve(superlu_options_t *options, SuperMatrix *mat,
int *perm_c,
int *perm_r,
int *etree,
205 char *equed,
double *R,
double *C, SuperMatrix *L, SuperMatrix *U,
206 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
207 double *rpg,
double *rcond,
double *ferr,
double *berr,
208 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
210 zgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
211 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
212 memusage, stat, info);
217 struct QuerySpaceChooser<
std::complex<double> >
219 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
221 zQuerySpace(L,U,memusage);
228 struct SuperLUDenseMatChooser<
std::complex<float> >
230 static void create(SuperMatrix *mat,
int n,
int m, std::complex<float> *dat,
int n1,
231 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
233 cCreate_Dense_Matrix(mat, n, m,
reinterpret_cast< ::complex*
>(dat), n1, stype, dtype, mtype);
237 static void destroy(SuperMatrix *mat)
242 struct SuperLUSolveChooser<
std::complex<float> >
244 static void solve(superlu_options_t *options, SuperMatrix *mat,
int *perm_c,
int *perm_r,
int *etree,
245 char *equed,
float *R,
float *C, SuperMatrix *L, SuperMatrix *U,
246 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
247 float *rpg,
float *rcond,
float *ferr,
float *berr,
248 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
250 cgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
251 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
252 memusage, stat, info);
257 struct QuerySpaceChooser<
std::complex<float> >
259 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
261 cQuerySpace(L,U,memusage);
279 template<
typename T,
typename A,
int n,
int m>
282 BlockVector<FieldVector<T,m>,
283 typename A::template rebind<FieldVector<T,m> >::other>,
284 BlockVector<FieldVector<T,n>,
285 typename A::template rebind<FieldVector<T,n> >::other> >
298 typename A::template rebind<FieldVector<T,m> >::other>
domain_type;
302 typename A::template rebind<FieldVector<T,n> >::other>
range_type;
317 explicit SuperLU(
const Matrix& mat,
bool verbose=
false,
318 bool reusevector=
true);
346 void apply(T* x, T* b);
349 void setMatrix(
const Matrix& mat);
351 typename SuperLUMatrix::size_type nnz()
const
357 void setSubMatrix(
const Matrix& mat,
const S& rowIndexSet);
359 void setVerbosity(
bool v);
367 const char* name() {
return "SuperLU"; }
369 friend class std::mem_fun_ref_t<void,SuperLU>;
370 template<
class M,
class X,
class TM,
class TD,
class T1>
371 friend class SeqOverlappingSchwarz;
372 friend struct SeqOverlappingSchwarzAssemblerHelper<SuperLU<Matrix>,true>;
378 SuperMatrix L, U, B, X;
379 int *perm_c, *perm_r, *etree;
380 typename GetSuperLUType<T>::float_type *R, *C;
382 superlu_options_t options;
386 bool first, verbose, reusevector;
389 template<
typename T,
typename A,
int n,
int m>
390 SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
393 if(mat.
N()+mat.
M()>0)
397 template<
typename T,
typename A,
int n,
int m>
398 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::free()
406 Destroy_SuperNode_Matrix(&L);
407 Destroy_CompCol_Matrix(&U);
410 if(!first && reusevector) {
411 SUPERLU_FREE(B.Store);
412 SUPERLU_FREE(X.Store);
417 template<
typename T,
typename A,
int n,
int m>
418 SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
419 ::SuperLU(
const Matrix& mat_,
bool verbose_,
bool reusevector_)
420 : work(0), lwork(0), first(true), verbose(verbose_),
421 reusevector(reusevector_)
426 template<
typename T,
typename A,
int n,
int m>
427 SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::SuperLU()
428 : work(0), lwork(0),verbose(false),
431 template<
typename T,
typename A,
int n,
int m>
432 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setVerbosity(
bool v)
437 template<
typename T,
typename A,
int n,
int m>
438 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setMatrix(
const Matrix& mat_)
440 if(mat.
N()+mat.
M()>0) {
450 template<
typename T,
typename A,
int n,
int m>
452 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::setSubMatrix(
const Matrix& mat_,
455 if(mat.
N()+mat.
M()>0) {
461 mat.setMatrix(mat_,mrs);
465 template<
typename T,
typename A,
int n,
int m>
466 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >::decompose()
470 perm_c =
new int[mat.
M()];
471 perm_r =
new int[mat.
N()];
472 etree =
new int[mat.
M()];
473 R =
new typename GetSuperLUType<T>::float_type[mat.
N()];
474 C =
new typename GetSuperLUType<T>::float_type[mat.
M()];
476 set_default_options(&options);
480 B.Dtype=GetSuperLUType<T>::type;
483 fakeFormat.lda=mat.
N();
486 X.Dtype=GetSuperLUType<T>::type;
491 typename GetSuperLUType<T>::float_type rpg, rcond, ferr=1e10, berr=1e10;
493 mem_usage_t memusage;
497 SuperLUSolveChooser<T>::solve(&options, &
static_cast<SuperMatrix&
>(mat), perm_c, perm_r, etree, &equed, R, C,
498 &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr,
499 &berr, &memusage, &stat, &info);
502 dinfo<<
"LU factorization: dgssvx() returns info "<< info<<std::endl;
504 if ( info == 0 || info == n+1 ) {
506 if ( options.PivotGrowth )
507 dinfo<<
"Recip. pivot growth = "<<rpg<<std::endl;
508 if ( options.ConditionNumber )
509 dinfo<<
"Recip. condition number = %e\n"<< rcond<<std::endl;
510 SCformat* Lstore = (SCformat *) L.Store;
511 NCformat* Ustore = (NCformat *) U.Store;
512 dinfo<<
"No of nonzeros in factor L = "<< Lstore->nnz<<std::endl;
513 dinfo<<
"No of nonzeros in factor U = "<< Ustore->nnz<<std::endl;
514 dinfo<<
"No of nonzeros in L+U = "<< Lstore->nnz + Ustore->nnz - n<<std::endl;
515 QuerySpaceChooser<T>::querySpace(&L, &U, &memusage);
516 dinfo<<
"L\\U MB "<<memusage.for_lu/1e6<<
" \ttotal MB needed "<<memusage.total_needed/1e6
519#ifdef HAVE_MEM_USAGE_T_EXPANSIONS
520 std::cout<<memusage.expansions<<std::endl;
522 std::cout<<stat.expansions<<std::endl;
524 }
else if ( info > 0 && lwork == -1 ) {
525 dinfo<<
"** Estimated memory: "<< info - n<<std::endl;
527 if ( options.PrintStat ) StatPrint(&stat);
555 options.Fact = FACTORED;
558 template<
typename T,
typename A,
int n,
int m>
559 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
562 if(mat.
M()+mat.
N()==0)
565 SuperMatrix* mB = &B;
566 SuperMatrix* mX = &X;
570 SuperLUDenseMatChooser<T>::create(&B, (
int)mat.
N(), 1,
reinterpret_cast<T*
>(&b[0]), (
int)mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
571 SuperLUDenseMatChooser<T>::create(&X, (
int)mat.
N(), 1,
reinterpret_cast<T*
>(&x[0]), (
int)mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
574 ((DNformat*) B.Store)->nzval=&b[0];
575 ((DNformat*)X.Store)->nzval=&x[0];
578 SuperLUDenseMatChooser<T>::create(&rB, (
int)mat.
N(), 1,
reinterpret_cast<T*
>(&b[0]), (
int)mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
579 SuperLUDenseMatChooser<T>::create(&rX, (
int)mat.
N(), 1,
reinterpret_cast<T*
>(&x[0]), (
int)mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
583 typename GetSuperLUType<T>::float_type rpg, rcond, ferr=1e10, berr;
585 mem_usage_t memusage;
595#ifdef SUPERLU_MIN_VERSION_4_3
596 options.IterRefine=SLU_DOUBLE;
598 options.IterRefine=DOUBLE;
601 SuperLUSolveChooser<T>::solve(&options, &
static_cast<SuperMatrix&
>(mat), perm_c, perm_r, etree, &equed, R, C,
602 &L, &U, work, lwork, mB, mX, &rpg, &rcond, &ferr, &berr,
603 &memusage, &stat, &info);
623 dinfo<<
"Triangular solve: dgssvx() returns info "<< info<<std::endl;
625 if ( info == 0 || info == n+1 ) {
627 if ( options.IterRefine ) {
628 std::cout<<
"Iterative Refinement: steps="
629 <<stat.RefineSteps<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
631 std::cout<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
632 }
else if ( info > 0 && lwork == -1 ) {
633 std::cout<<
"** Estimated memory: "<< info - n<<
" bytes"<<std::endl;
636 if ( options.PrintStat ) StatPrint(&stat);
640 SUPERLU_FREE(rB.Store);
641 SUPERLU_FREE(rX.Store);
645 template<
typename T,
typename A,
int n,
int m>
646 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> >
649 if(mat.
N()+mat.
M()==0)
652 SuperMatrix* mB = &B;
653 SuperMatrix* mX = &X;
657 SuperLUDenseMatChooser<T>::create(&B, mat.
N(), 1, b, mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
658 SuperLUDenseMatChooser<T>::create(&X, mat.
N(), 1, x, mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
661 ((DNformat*) B.Store)->nzval=b;
662 ((DNformat*)X.Store)->nzval=x;
665 SuperLUDenseMatChooser<T>::create(&rB, mat.
N(), 1, b, mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
666 SuperLUDenseMatChooser<T>::create(&rX, mat.
N(), 1, x, mat.
N(), SLU_DN, GetSuperLUType<T>::type, SLU_GE);
671 typename GetSuperLUType<T>::float_type rpg, rcond, ferr=1e10, berr;
673 mem_usage_t memusage;
678#ifdef SUPERLU_MIN_VERSION_4_3
679 options.IterRefine=SLU_DOUBLE;
681 options.IterRefine=DOUBLE;
684 SuperLUSolveChooser<T>::solve(&options, &
static_cast<SuperMatrix&
>(mat), perm_c, perm_r, etree, &equed, R, C,
685 &L, &U, work, lwork, mB, mX, &rpg, &rcond, &ferr, &berr,
686 &memusage, &stat, &info);
689 dinfo<<
"Triangular solve: dgssvx() returns info "<< info<<std::endl;
691 if ( info == 0 || info == n+1 ) {
693 if ( options.IterRefine ) {
694 dinfo<<
"Iterative Refinement: steps="
695 <<stat.RefineSteps<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
697 dinfo<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
698 }
else if ( info > 0 && lwork == -1 ) {
699 dinfo<<
"** Estimated memory: "<< info - n<<
" bytes"<<std::endl;
701 if ( options.PrintStat ) StatPrint(&stat);
706 SUPERLU_FREE(rB.Store);
707 SUPERLU_FREE(rX.Store);
712 template<
typename T,
typename A,
int n,
int m>
718 template<
typename T,
typename A,
int n,
int m>
719 struct StoresColumnCompressed<SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
721 enum { value =
true };
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:414
A vector of blocks with memory management.
Definition: bvector.hh:254
A dense n x m matrix.
Definition: fmatrix.hh:67
vector space out of a tensor product of fields.
Definition: fvector.hh:92
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:302
Dune::SuperLUMatrix< Matrix > SuperLUMatrix
The corresponding SuperLU Matrix type.
Definition: superlu.hh:292
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:298
void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: superlu.hh:337
SuperMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A > > MatrixInitializer
Type of an associated initializer class.
Definition: superlu.hh:294
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:244
DInfoType dinfo(std::cout)
Stream for informative output.
Definition: stdstreams.hh:139
Dune namespace.
Definition: alignment.hh:14
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