37 #if defined(GMM_USES_SUPERLU) && !defined(GETFEM_VERSION)
39 #ifndef GMM_SUPERLU_INTERFACE_H
40 #define GMM_SUPERLU_INTERFACE_H
54 #include "superlu/slu_Cnames.h"
55 #include "superlu/supermatrix.h"
56 #include "superlu/slu_util.h"
59 #include "superlu/slu_sdefs.h"
62 #include "superlu/slu_ddefs.h"
65 #include "superlu/slu_cdefs.h"
68 #include "superlu/slu_zdefs.h"
77 inline void Create_CompCol_Matrix(SuperMatrix *A,
int m,
int n,
int nnz,
78 float *a,
int *ir,
int *jc) {
79 SuperLU_S::sCreate_CompCol_Matrix(A, m, n, nnz, a, ir, jc,
80 SLU_NC, SLU_S, SLU_GE);
83 inline void Create_CompCol_Matrix(SuperMatrix *A,
int m,
int n,
int nnz,
84 double *a,
int *ir,
int *jc) {
85 SuperLU_D::dCreate_CompCol_Matrix(A, m, n, nnz, a, ir, jc,
86 SLU_NC, SLU_D, SLU_GE);
89 inline void Create_CompCol_Matrix(SuperMatrix *A,
int m,
int n,
int nnz,
90 std::complex<float> *a,
int *ir,
int *jc) {
91 SuperLU_C::cCreate_CompCol_Matrix(A, m, n, nnz, (SuperLU_C::complex *)(a),
92 ir, jc, SLU_NC, SLU_C, SLU_GE);
95 inline void Create_CompCol_Matrix(SuperMatrix *A,
int m,
int n,
int nnz,
96 std::complex<double> *a,
int *ir,
int *jc) {
97 SuperLU_Z::zCreate_CompCol_Matrix(A, m, n, nnz,
98 (SuperLU_Z::doublecomplex *)(a), ir, jc,
99 SLU_NC, SLU_Z, SLU_GE);
104 inline void Create_Dense_Matrix(SuperMatrix *A,
int m,
int n,
float *a,
int k)
105 { SuperLU_S::sCreate_Dense_Matrix(A, m, n, a, k, SLU_DN, SLU_S, SLU_GE); }
106 inline void Create_Dense_Matrix(SuperMatrix *A,
int m,
int n,
double *a,
int k)
107 { SuperLU_D::dCreate_Dense_Matrix(A, m, n, a, k, SLU_DN, SLU_D, SLU_GE); }
108 inline void Create_Dense_Matrix(SuperMatrix *A,
int m,
int n,
109 std::complex<float> *a,
int k) {
110 SuperLU_C::cCreate_Dense_Matrix(A, m, n, (SuperLU_C::complex *)(a),
111 k, SLU_DN, SLU_C, SLU_GE);
113 inline void Create_Dense_Matrix(SuperMatrix *A,
int m,
int n,
114 std::complex<double> *a,
int k) {
115 SuperLU_Z::zCreate_Dense_Matrix(A, m, n, (SuperLU_Z::doublecomplex *)(a),
116 k, SLU_DN, SLU_Z, SLU_GE);
121 #define DECL_GSSV(NAMESPACE,FNAME,FLOATTYPE,KEYTYPE) \
122 inline void SuperLU_gssv(superlu_options_t *options, SuperMatrix *A, int *p, \
123 int *q, SuperMatrix *L, SuperMatrix *U, SuperMatrix *B, \
124 SuperLUStat_t *stats, int *info, KEYTYPE) { \
125 NAMESPACE::FNAME(options, A, p, q, L, U, B, stats, info); \
128 DECL_GSSV(SuperLU_S,sgssv,
float,
float)
129 DECL_GSSV(SuperLU_C,cgssv,
float,std::complex<float>)
130 DECL_GSSV(SuperLU_D,dgssv,
double,
double)
131 DECL_GSSV(SuperLU_Z,zgssv,
double,std::complex<double>)
135 #define DECL_GSSVX(NAMESPACE,FNAME,FLOATTYPE,KEYTYPE) \
136 inline float SuperLU_gssvx(superlu_options_t *options, SuperMatrix *A, \
137 int *perm_c, int *perm_r, int *etree, char *equed, \
138 FLOATTYPE *R, FLOATTYPE *C, SuperMatrix *L, \
139 SuperMatrix *U, void *work, int lwork, \
140 SuperMatrix *B, SuperMatrix *X, \
141 FLOATTYPE *recip_pivot_growth, \
142 FLOATTYPE *rcond, FLOATTYPE *ferr, FLOATTYPE *berr, \
143 SuperLUStat_t *stats, int *info, KEYTYPE) { \
144 NAMESPACE::mem_usage_t mem_usage; \
145 NAMESPACE::FNAME(options, A, perm_c, perm_r, etree, equed, R, C, L, \
146 U, work, lwork, B, X, recip_pivot_growth, rcond, \
147 ferr, berr, &mem_usage, stats, info); \
148 return mem_usage.for_lu; \
151 DECL_GSSVX(SuperLU_S,sgssvx,
float,
float)
152 DECL_GSSVX(SuperLU_C,cgssvx,
float,std::complex<float>)
153 DECL_GSSVX(SuperLU_D,dgssvx,
double,
double)
154 DECL_GSSVX(SuperLU_Z,zgssvx,
double,std::complex<double>)
160 template <
typename MAT,
typename VECTX,
typename VECTB>
161 int SuperLU_solve(
const MAT &A,
const VECTX &X_,
const VECTB &B,
162 double& rcond_,
int permc_spec = 3) {
163 VECTX &X =
const_cast<VECTX &
>(X_);
171 typedef typename linalg_traits<MAT>::value_type T;
172 typedef typename number_traits<T>::magnitude_type R;
174 int m = mat_nrows(A), n = mat_ncols(A), nrhs = 1, info = 0;
176 csc_matrix<T> csc_A(m, n); gmm::copy(A, csc_A);
177 std::vector<T> rhs(m), sol(m);
181 if ((2 * nz / n) >= m)
182 GMM_WARNING2(
"CAUTION : it seems that SuperLU has a problem"
183 " for nearly dense sparse matrices");
185 superlu_options_t options;
186 set_default_options(&options);
187 options.ColPerm = NATURAL;
188 options.PrintStat = NO;
189 options.ConditionNumber = YES;
190 switch (permc_spec) {
191 case 1 : options.ColPerm = MMD_ATA;
break;
192 case 2 : options.ColPerm = MMD_AT_PLUS_A;
break;
193 case 3 : options.ColPerm = COLAMD;
break;
198 SuperMatrix SA, SL, SU, SB, SX;
199 Create_CompCol_Matrix(&SA, m, n, nz, (
double *)(&(csc_A.pr[0])),
200 (
int *)(&(csc_A.ir[0])), (
int *)(&(csc_A.jc[0])));
201 Create_Dense_Matrix(&SB, m, nrhs, &rhs[0], m);
202 Create_Dense_Matrix(&SX, m, nrhs, &sol[0], m);
203 memset(&SL,0,
sizeof SL);
204 memset(&SU,0,
sizeof SU);
206 std::vector<int> etree(n);
208 std::vector<R> Rscale(m),Cscale(n);
209 std::vector<R> ferr(nrhs), berr(nrhs);
210 R recip_pivot_gross, rcond;
211 std::vector<int> perm_r(m), perm_c(n);
213 SuperLU_gssvx(&options, &SA, &perm_c[0], &perm_r[0],
229 Destroy_SuperMatrix_Store(&SB);
230 Destroy_SuperMatrix_Store(&SX);
231 Destroy_SuperMatrix_Store(&SA);
232 Destroy_SuperNode_Matrix(&SL);
233 Destroy_CompCol_Matrix(&SU);
235 GMM_ASSERT1(info >= 0,
"SuperLU solve failed: info =" << info);
236 if (info > 0) GMM_WARNING1(
"SuperLU solve failed: info =" << info);
241 template <
class T>
class SuperLU_factor {
242 typedef typename number_traits<T>::magnitude_type R;
245 mutable SuperMatrix SA, SL, SB, SU, SX;
246 mutable SuperLUStat_t stat;
247 mutable superlu_options_t options;
249 mutable std::vector<int> etree, perm_r, perm_c;
250 mutable std::vector<R> Rscale, Cscale;
251 mutable std::vector<R> ferr, berr;
252 mutable std::vector<T> rhs;
253 mutable std::vector<T> sol;
254 mutable bool is_init;
258 enum { LU_NOTRANSP, LU_TRANSP, LU_CONJUGATED };
259 void free_supermatrix(
void);
260 template <
class MAT>
void build_with(
const MAT &A,
int permc_spec = 3);
261 template <
typename VECTX,
typename VECTB>
265 void solve(
const VECTX &X_,
const VECTB &B,
int transp=LU_NOTRANSP)
const;
266 SuperLU_factor(
void) { is_init =
false; }
267 SuperLU_factor(
const SuperLU_factor& other) {
268 GMM_ASSERT2(!(other.is_init),
269 "copy of initialized SuperLU_factor is forbidden");
272 SuperLU_factor& operator=(
const SuperLU_factor& other) {
273 GMM_ASSERT2(!(other.is_init) && !is_init,
274 "assignment of initialized SuperLU_factor is forbidden");
277 ~SuperLU_factor() { free_supermatrix(); }
278 float memsize() {
return memory_used; }
282 template <
class T>
void SuperLU_factor<T>::free_supermatrix(
void) {
284 if (SB.Store) Destroy_SuperMatrix_Store(&SB);
285 if (SX.Store) Destroy_SuperMatrix_Store(&SX);
286 if (SA.Store) Destroy_SuperMatrix_Store(&SA);
287 if (SL.Store) Destroy_SuperNode_Matrix(&SL);
288 if (SU.Store) Destroy_CompCol_Matrix(&SU);
293 template <
class T>
template <
class MAT>
303 int n = mat_nrows(A), m = mat_ncols(A), info = 0;
306 rhs.resize(m); sol.resize(m);
310 set_default_options(&options);
311 options.ColPerm = NATURAL;
312 options.PrintStat = NO;
313 options.ConditionNumber = NO;
314 switch (permc_spec) {
315 case 1 : options.ColPerm = MMD_ATA;
break;
316 case 2 : options.ColPerm = MMD_AT_PLUS_A;
break;
317 case 3 : options.ColPerm = COLAMD;
break;
321 Create_CompCol_Matrix(&SA, m, n, nz, (
double *)(&(csc_A.pr[0])),
322 (
int *)(&(csc_A.ir[0])), (
int *)(&(csc_A.jc[0])));
324 Create_Dense_Matrix(&SB, m, 0, &rhs[0], m);
325 Create_Dense_Matrix(&SX, m, 0, &sol[0], m);
326 memset(&SL,0,
sizeof SL);
327 memset(&SU,0,
sizeof SU);
329 Rscale.resize(m); Cscale.resize(n); etree.resize(n);
330 ferr.resize(1); berr.resize(1);
331 R recip_pivot_gross, rcond;
332 perm_r.resize(m); perm_c.resize(n);
333 memory_used = SuperLU_gssvx(&options, &SA, &perm_c[0], &perm_r[0],
349 Destroy_SuperMatrix_Store(&SB);
350 Destroy_SuperMatrix_Store(&SX);
351 Create_Dense_Matrix(&SB, m, 1, &rhs[0], m);
352 Create_Dense_Matrix(&SX, m, 1, &sol[0], m);
355 GMM_ASSERT1(info == 0,
"SuperLU solve failed: info=" << info);
359 template <
class T>
template <
typename VECTX,
typename VECTB>
362 VECTX &X =
const_cast<VECTX &
>(X_);
364 options.Fact = FACTORED;
365 options.IterRefine = NOREFINE;
367 case LU_NOTRANSP: options.Trans = NOTRANS;
break;
368 case LU_TRANSP: options.Trans = TRANS;
break;
369 case LU_CONJUGATED: options.Trans = CONJ;
break;
370 default: GMM_ASSERT1(
false,
"invalid value for transposition option");
374 R recip_pivot_gross, rcond;
375 SuperLU_gssvx(&options, &SA, &perm_c[0], &perm_r[0],
391 GMM_ASSERT1(info == 0,
"SuperLU solve failed: info=" << info);
395 template <
typename T,
typename V1,
typename V2>
inline
396 void mult(
const SuperLU_factor<T>& P,
const V1 &v1,
const V2 &v2) {
400 template <
typename T,
typename V1,
typename V2>
inline
401 void transposed_mult(
const SuperLU_factor<T>& P,
const V1 &v1,
const V2 &v2) {
402 P.solve(v2, v1, SuperLU_factor<T>::LU_TRANSP);
408 #endif // GMM_SUPERLU_INTERFACE_H
410 #endif // GMM_USES_SUPERLU