1 #ifndef DUNE_SUPERLU_HH
2 #define DUNE_SUPERLU_HH
12 #define SUPERLU_NTYPE 1
14 #ifdef SUPERLU_POST_2005_VERSION
17 #include "slu_sdefs.h"
21 #include "slu_ddefs.h"
25 #include "slu_cdefs.h"
29 #include "slu_zdefs.h"
40 #warning Support for SuperLU older than SuperLU 3.0 from August 2005 is deprecated.
59 #include<dune/common/fmatrix.hh>
60 #include<dune/common/fvector.hh>
61 #include<dune/common/stdstreams.hh>
77 template<
class Matrix>
82 template<
class M,
class T,
class TM,
class TD,
class TA>
112 static void create(SuperMatrix *
mat,
int n,
int m,
float *dat,
int n1,
113 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
115 sCreate_Dense_Matrix(mat, n, m, dat, n1, stype, dtype, mtype);
119 static void destroy(SuperMatrix *m)
125 struct SuperLUSolveChooser<float>
127 static void solve(superlu_options_t *options, SuperMatrix *
mat,
int *perm_c,
int *perm_r,
int *etree,
128 char *equed,
float *R,
float *C, SuperMatrix *L, SuperMatrix *U,
129 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
130 float *rpg,
float *rcond,
float *ferr,
float *berr,
131 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
133 sgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
134 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
135 memusage, stat, info);
140 struct QuerySpaceChooser<float>
142 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
144 sQuerySpace(L,U,memusage);
155 static void create(SuperMatrix *
mat,
int n,
int m,
double *dat,
int n1,
156 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
158 dCreate_Dense_Matrix(mat, n, m, dat, n1, stype, dtype, mtype);
162 static void destroy(SuperMatrix *
mat)
169 static void solve(superlu_options_t *options, SuperMatrix *
mat,
int *perm_c,
int *perm_r,
int *etree,
170 char *equed,
double *R,
double *C, SuperMatrix *L, SuperMatrix *U,
171 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
172 double *rpg,
double *rcond,
double *ferr,
double *berr,
173 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
175 dgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
176 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
177 memusage, stat, info);
184 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
186 dQuerySpace(L,U,memusage);
193 struct SuperLUDenseMatChooser<std::complex<double> >
195 static void create(SuperMatrix *
mat,
int n,
int m, std::complex<double> *dat,
int n1,
196 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
198 zCreate_Dense_Matrix(mat, n, m, reinterpret_cast<doublecomplex*>(dat), n1, stype, dtype, mtype);
202 static void destroy(SuperMatrix *mat)
208 struct SuperLUSolveChooser<std::complex<double> >
210 static void solve(superlu_options_t *options, SuperMatrix *
mat,
int *perm_c,
int *perm_r,
int *etree,
211 char *equed,
double *R,
double *C, SuperMatrix *L, SuperMatrix *U,
212 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
213 double *rpg,
double *rcond,
double *ferr,
double *berr,
214 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
216 zgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
217 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
218 memusage, stat, info);
223 struct QuerySpaceChooser<std::complex<double> >
225 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
227 zQuerySpace(L,U,memusage);
234 struct SuperLUDenseMatChooser<std::complex<float> >
236 static void create(SuperMatrix *
mat,
int n,
int m, std::complex<float> *dat,
int n1,
237 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
239 cCreate_Dense_Matrix(mat, n, m, reinterpret_cast< ::complex*>(dat), n1, stype, dtype, mtype);
243 static void destroy(SuperMatrix *mat)
249 struct SuperLUSolveChooser<std::complex<float> >
251 static void solve(superlu_options_t *options, SuperMatrix *
mat,
int *perm_c,
int *perm_r,
int *etree,
252 char *equed,
float *R,
float *C, SuperMatrix *L, SuperMatrix *U,
253 void *work,
int lwork, SuperMatrix *B, SuperMatrix *X,
254 float *rpg,
float *rcond,
float *ferr,
float *berr,
255 mem_usage_t *memusage, SuperLUStat_t *stat,
int *info)
257 cgssvx(options, mat, perm_c, perm_r, etree, equed, R, C,
258 L, U, work, lwork, B, X, rpg, rcond, ferr, berr,
259 memusage, stat, info);
264 struct QuerySpaceChooser<std::complex<float> >
266 static void querySpace(SuperMatrix* L, SuperMatrix* U, mem_usage_t* memusage)
268 cQuerySpace(L,U,memusage);
286 template<
typename T,
typename A,
int n,
int m>
289 BlockVector<FieldVector<T,m>,
290 typename A::template rebind<FieldVector<T,m> >::other>,
291 BlockVector<FieldVector<T,n>,
292 typename A::template rebind<FieldVector<T,n> >::other> >
302 typename A::template rebind<FieldVector<T,m> >::other>
domain_type;
306 typename A::template rebind<FieldVector<T,n> >::other>
range_type;
342 void apply(T* x, T* b);
347 typename SuperLUMatrix::size_type nnz()
const
353 void setSubMatrix(
const Matrix&
mat,
const S& rowIndexSet);
355 void setVerbosity(
bool v);
363 friend class std::mem_fun_ref_t<void,
SuperLU>;
364 template<
class M,
class X,
class TM,
class TD,
class T1>
372 SuperMatrix L, U, B, X;
373 int *perm_c, *perm_r, *etree;
376 superlu_options_t options;
383 template<
typename T,
typename A,
int n,
int m>
391 template<
typename T,
typename A,
int n,
int m>
400 Destroy_SuperNode_Matrix(&L);
401 Destroy_CompCol_Matrix(&U);
405 SUPERLU_FREE(B.Store);
406 SUPERLU_FREE(X.Store);
411 template<
typename T,
typename A,
int n,
int m>
414 : work(0), lwork(0), first(true), verbose(verbose_)
419 template<
typename T,
typename A,
int n,
int m>
421 : work(0), lwork(0),verbose(false)
423 template<
typename T,
typename A,
int n,
int m>
429 template<
typename T,
typename A,
int n,
int m>
442 template<
typename T,
typename A,
int n,
int m>
453 mat.setMatrix(mat_,mrs);
457 template<
typename T,
typename A,
int n,
int m>
462 perm_c =
new int[
mat.
M()];
463 perm_r =
new int[
mat.
N()];
464 etree =
new int[
mat.
M()];
468 set_default_options(&options);
475 fakeFormat.lda=
mat.
N();
485 mem_usage_t memusage;
490 &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr,
491 &berr, &memusage, &stat, &info);
494 dinfo<<
"LU factorization: dgssvx() returns info "<< info<<std::endl;
496 if ( info == 0 || info == n+1 ) {
498 if ( options.PivotGrowth )
499 dinfo<<
"Recip. pivot growth = "<<rpg<<std::endl;
500 if ( options.ConditionNumber )
501 dinfo<<
"Recip. condition number = %e\n"<< rcond<<std::endl;
502 SCformat* Lstore = (SCformat *) L.Store;
503 NCformat* Ustore = (NCformat *) U.Store;
504 dinfo<<
"No of nonzeros in factor L = "<< Lstore->nnz<<std::endl;
505 dinfo<<
"No of nonzeros in factor U = "<< Ustore->nnz<<std::endl;
506 dinfo<<
"No of nonzeros in L+U = "<< Lstore->nnz + Ustore->nnz - n<<std::endl;
508 dinfo<<
"L\\U MB "<<memusage.for_lu/1e6<<
" \ttotal MB needed "<<memusage.total_needed/1e6
511 #ifdef HAVE_MEM_USAGE_T_EXPANSIONS
512 std::cout<<memusage.expansions<<std::endl;
514 std::cout<<stat.expansions<<std::endl;
516 }
else if ( info > 0 && lwork == -1 ) {
517 dinfo<<
"** Estimated memory: "<< info - n<<std::endl;
519 if ( options.PrintStat ) StatPrint(&stat);
547 options.Fact = FACTORED;
550 template<
typename T,
typename A,
int n,
int m>
551 void SuperLU<BCRSMatrix<FieldMatrix<T,n,m>,
A> >
555 DUNE_THROW(
ISTLError,
"Matrix of SuperLU is null!");
562 ((DNformat*) B.Store)->nzval=&b[0];
563 ((DNformat*)X.Store)->nzval=&x[0];
567 mem_usage_t memusage;
577 #ifdef SUPERLU_MIN_VERSION_4_3
578 options.IterRefine=SLU_DOUBLE;
580 options.IterRefine=DOUBLE;
584 &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr, &berr,
585 &memusage, &stat, &info);
605 dinfo<<
"Triangular solve: dgssvx() returns info "<< info<<std::endl;
607 if ( info == 0 || info == n+1 ) {
609 if ( options.IterRefine ) {
610 std::cout<<
"Iterative Refinement: steps="
611 <<stat.RefineSteps<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
613 std::cout<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
614 }
else if ( info > 0 && lwork == -1 ) {
615 std::cout<<
"** Estimated memory: "<< info - n<<
" bytes"<<std::endl;
618 if ( options.PrintStat ) StatPrint(&stat);
623 template<
typename T,
typename A,
int n,
int m>
628 DUNE_THROW(
ISTLError,
"Matrix of SuperLU is null!");
635 ((DNformat*) B.Store)->nzval=b;
636 ((DNformat*)X.Store)->nzval=x;
641 mem_usage_t memusage;
646 #ifdef SUPERLU_MIN_VERSION_4_3
647 options.IterRefine=SLU_DOUBLE;
649 options.IterRefine=DOUBLE;
653 &L, &U, work, lwork, &B, &X, &rpg, &rcond, &ferr, &berr,
654 &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 dinfo<<
"Iterative Refinement: steps="
663 <<stat.RefineSteps<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
665 dinfo<<
" FERR="<<ferr<<
" BERR="<<berr<<std::endl;
666 }
else if ( info > 0 && lwork == -1 ) {
667 dinfo<<
"** Estimated memory: "<< info - n<<
" bytes"<<std::endl;
669 if ( options.PrintStat ) StatPrint(&stat);
676 template<
typename T,
typename A,
int n,
int m>
683 #endif // HAVE_SUPERLU
684 #endif // DUNE_SUPERLU_HH