3 #ifndef DUNE_SUPERLU_HH
4 #define DUNE_SUPERLU_HH
14 #define SUPERLU_NTYPE 1
16 #ifdef SUPERLU_POST_2005_VERSION
19 #include "slu_sdefs.h"
23 #include "slu_ddefs.h"
27 #include "slu_cdefs.h"
31 #include "slu_zdefs.h"
42 #warning Support for SuperLU older than SuperLU 3.0 from August 2005 is deprecated.
61 #include <dune/common/fmatrix.hh>
62 #include <dune/common/fvector.hh>
63 #include <dune/common/stdstreams.hh>
79 template<
class Matrix>
83 template<
class M,
class T,
class TM,
class TD,
class TA>
86 template<
class T,
bool tag>
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);
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 * )
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);
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;
318 bool reusevector=
true);
339 DUNE_UNUSED_PARAMETER(reduction);
346 void apply(T* x, T* b);
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>
378 SuperMatrix L, U, B, X;
379 int *perm_c, *perm_r, *etree;
382 superlu_options_t options;
386 bool first, verbose, reusevector;
389 template<
typename T,
typename A,
int n,
int m>
397 template<
typename T,
typename A,
int n,
int m>
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>
420 : work(0), lwork(0), first(true), verbose(verbose_),
421 reusevector(reusevector_)
426 template<
typename T,
typename A,
int n,
int m>
428 : work(0), lwork(0),verbose(false),
431 template<
typename T,
typename A,
int n,
int m>
437 template<
typename T,
typename A,
int n,
int m>
450 template<
typename T,
typename A,
int n,
int m>
461 mat.setMatrix(mat_,mrs);
465 template<
typename T,
typename A,
int n,
int m>
470 perm_c =
new int[
mat.
M()];
471 perm_r =
new int[
mat.
N()];
472 etree =
new int[
mat.
M()];
476 set_default_options(&options);
483 fakeFormat.lda=
mat.
N();
493 mem_usage_t memusage;
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;
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> >
563 DUNE_THROW(
ISTLError,
"Matrix of SuperLU is null!");
565 SuperMatrix* mB = &B;
566 SuperMatrix* mX = &X;
574 ((DNformat*) B.Store)->nzval=&b[0];
575 ((DNformat*)X.Store)->nzval=&x[0];
585 mem_usage_t memusage;
595 #ifdef SUPERLU_MIN_VERSION_4_3
596 options.IterRefine=SLU_DOUBLE;
598 options.IterRefine=DOUBLE;
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>
650 DUNE_THROW(
ISTLError,
"Matrix of SuperLU is null!");
652 SuperMatrix* mB = &B;
653 SuperMatrix* mX = &X;
661 ((DNformat*) B.Store)->nzval=b;
662 ((DNformat*)X.Store)->nzval=x;
673 mem_usage_t memusage;
678 #ifdef SUPERLU_MIN_VERSION_4_3
679 options.IterRefine=SLU_DOUBLE;
681 options.IterRefine=DOUBLE;
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>
725 #endif // HAVE_SUPERLU
726 #endif // DUNE_SUPERLU_HH