53 #ifndef AMESOS2_SUPERLUMT_FUNCTIONMAP_HPP
54 #define AMESOS2_SUPERLUMT_FUNCTIONMAP_HPP
56 #ifdef HAVE_TEUCHOS_COMPLEX
61 #include "Amesos2_MatrixAdapter.hpp"
71 #include "slu_mt_util.h"
72 #include "pxgstrf_synch.h"
75 #include "pssp_defs.h"
79 #include "pdsp_defs.h"
82 #ifdef HAVE_TEUCHOS_COMPLEX
84 #include "pcsp_defs.h"
88 #include "pzsp_defs.h"
90 #endif // HAVE_TEUCHOS_COMPLEX
98 template <
class Matrix,
class Vector>
class Superlumt;
135 struct FunctionMap<Superlumt,float>
137 typedef TypeMap<Superlumt,float> type_map;
142 static void gssvx(SLUMT::superlumt_options_t* options, SLUMT::SuperMatrix* A,
143 int* perm_c,
int* perm_r,
int* etree, SLUMT::equed_t* equed,
float* R,
float* C,
144 SLUMT::SuperMatrix* L, SLUMT::SuperMatrix* U,
void* work,
int lwork,
145 SLUMT::SuperMatrix* B, SLUMT::SuperMatrix* X,
float* recip_pivot_growth,
146 float* rcond,
float* ferr,
float* berr, SLUMT::superlu_memusage_t* mem_usage,
147 SLUMT::Gstat_t* stat,
int* info)
149 options->etree = etree;
150 options->perm_c = perm_c;
151 options->perm_r = perm_r;
153 options->work = work;
154 options->lwork = lwork;
156 SLUMT::S::psgssvx(options->nprocs, options, A, perm_c, perm_r,
157 equed, R, C, L, U, B, X, recip_pivot_growth, rcond, ferr,
158 berr, mem_usage, info);
164 static void gstrs(SLUMT::trans_t trans, SLUMT::SuperMatrix* L,
165 SLUMT::SuperMatrix* U,
int* perm_r,
int* perm_c,
166 SLUMT::SuperMatrix* B, SLUMT::Gstat_t* Gstat,
int* info)
168 SLUMT::S::sgstrs(trans, L, U, perm_r, perm_c, B, Gstat, info);
188 static void gstrf(SLUMT::superlumt_options_t* options, SLUMT::SuperMatrix* A,
189 int* perm_r, SLUMT::SuperMatrix* L, SLUMT::SuperMatrix* U,
190 SLUMT::Gstat_t* stat,
int* info)
192 SLUMT::S::psgstrf(options, A, perm_r, L, U, stat, info);
208 static void create_CompCol_Matrix(SLUMT::SuperMatrix* A,
int m,
int n,
int nnz,
209 type_map::type* nzval,
int* rowind,
int* colptr,
210 SLUMT::Stype_t stype, SLUMT::Dtype_t dtype, SLUMT::Mtype_t mtype)
212 SLUMT::S::sCreate_CompCol_Matrix(A, m, n, nnz, nzval, rowind, colptr,
213 stype, dtype, mtype);
227 static void create_Dense_Matrix(SLUMT::SuperMatrix* X,
int m,
int n,
228 type_map::type* x,
int ldx, SLUMT::Stype_t stype,
229 SLUMT::Dtype_t dtype, SLUMT::Mtype_t mtype)
231 SLUMT::S::sCreate_Dense_Matrix(X, m, n, x, ldx, stype, dtype, mtype);
240 static void gsequ(SLUMT::SuperMatrix* A,
241 type_map::magnitude_type* r,
242 type_map::magnitude_type* c,
243 type_map::magnitude_type* rowcnd,
244 type_map::magnitude_type* colcnd,
245 type_map::magnitude_type* amax,
248 SLUMT::S::sgsequ(A, r, c, rowcnd, colcnd, amax, info);
259 static void laqgs(SLUMT::SuperMatrix* A,
260 type_map::magnitude_type* r,
261 type_map::magnitude_type* c,
262 type_map::magnitude_type rowcnd,
263 type_map::magnitude_type colcnd,
264 type_map::magnitude_type amax,
265 SLUMT::equed_t* equed)
267 SLUMT::S::slaqgs(A, r, c, rowcnd, colcnd, amax, equed);
273 struct FunctionMap<Superlumt,double>
275 typedef TypeMap<Superlumt,double> type_map;
277 static void gssvx(SLUMT::superlumt_options_t* options, SLUMT::SuperMatrix* A,
278 int* perm_c,
int* perm_r,
int* etree, SLUMT::equed_t* equed,
double* R,
double* C,
279 SLUMT::SuperMatrix* L, SLUMT::SuperMatrix* U,
void* work,
int lwork,
280 SLUMT::SuperMatrix* B, SLUMT::SuperMatrix* X,
double* recip_pivot_growth,
281 double* rcond,
double* ferr,
double* berr, SLUMT::superlu_memusage_t* mem_usage,
282 SLUMT::Gstat_t* stat,
int* info)
284 options->etree = etree;
285 options->perm_c = perm_c;
286 options->perm_r = perm_r;
288 options->work = work;
289 options->lwork = lwork;
291 SLUMT::D::pdgssvx(options->nprocs, options, A, perm_c, perm_r,
292 equed, R, C, L, U, B, X, recip_pivot_growth, rcond, ferr,
293 berr, mem_usage, info);
296 static void gstrs(SLUMT::trans_t trans, SLUMT::SuperMatrix* L,
297 SLUMT::SuperMatrix* U,
int* perm_r,
int* perm_c,
298 SLUMT::SuperMatrix* B, SLUMT::Gstat_t* Gstat,
int* info)
300 SLUMT::D::dgstrs(trans, L, U, perm_r, perm_c, B, Gstat, info);
303 static void gstrf(SLUMT::superlumt_options_t* options, SLUMT::SuperMatrix* A,
304 int* perm_r, SLUMT::SuperMatrix* L, SLUMT::SuperMatrix* U,
305 SLUMT::Gstat_t* stat,
int* info)
307 SLUMT::D::pdgstrf(options, A, perm_r, L, U, stat, info);
310 static void create_CompCol_Matrix(SLUMT::SuperMatrix* A,
int m,
int n,
int nnz,
311 type_map::type* nzval,
int* rowind,
int* colptr,
312 SLUMT::Stype_t stype, SLUMT::Dtype_t dtype, SLUMT::Mtype_t mtype)
314 SLUMT::D::dCreate_CompCol_Matrix(A, m, n, nnz, nzval, rowind, colptr,
315 stype, dtype, mtype);
318 static void create_Dense_Matrix(SLUMT::SuperMatrix* X,
int m,
int n,
319 type_map::type* x,
int ldx, SLUMT::Stype_t stype,
320 SLUMT::Dtype_t dtype, SLUMT::Mtype_t mtype)
322 SLUMT::D::dCreate_Dense_Matrix(X, m, n, x, ldx, stype, dtype, mtype);
325 static void gsequ(SLUMT::SuperMatrix* A,
326 type_map::magnitude_type* r,
327 type_map::magnitude_type* c,
328 type_map::magnitude_type* rowcnd,
329 type_map::magnitude_type* colcnd,
330 type_map::magnitude_type* amax,
333 SLUMT::D::dgsequ(A, r, c, rowcnd, colcnd, amax, info);
336 static void laqgs(SLUMT::SuperMatrix* A,
337 type_map::magnitude_type* r,
338 type_map::magnitude_type* c,
339 type_map::magnitude_type rowcnd,
340 type_map::magnitude_type colcnd,
341 type_map::magnitude_type amax,
342 SLUMT::equed_t* equed)
344 SLUMT::D::dlaqgs(A, r, c, rowcnd, colcnd, amax, equed);
349 #ifdef HAVE_TEUCHOS_COMPLEX
354 struct FunctionMap<Superlumt,SLUMT::C::complex>
356 typedef TypeMap<Superlumt,SLUMT::C::complex> type_map;
358 static void gssvx(SLUMT::superlumt_options_t* options, SLUMT::SuperMatrix* A,
359 int* perm_c,
int* perm_r,
int* etree, SLUMT::equed_t* equed,
float* R,
float* C,
360 SLUMT::SuperMatrix* L, SLUMT::SuperMatrix* U,
void* work,
int lwork,
361 SLUMT::SuperMatrix* B, SLUMT::SuperMatrix* X,
float* recip_pivot_growth,
362 float* rcond,
float* ferr,
float* berr, SLUMT::superlu_memusage_t* mem_usage,
363 SLUMT::Gstat_t* stat,
int* info)
365 options->etree = etree;
366 options->perm_c = perm_c;
367 options->perm_r = perm_r;
369 options->work = work;
370 options->lwork = lwork;
372 SLUMT::C::pcgssvx(options->nprocs, options, A, perm_c, perm_r,
373 equed, R, C, L, U, B, X, recip_pivot_growth, rcond, ferr,
374 berr, mem_usage, info);
377 static void gstrs(SLUMT::trans_t trans, SLUMT::SuperMatrix* L,
378 SLUMT::SuperMatrix* U,
int* perm_r,
int* perm_c,
379 SLUMT::SuperMatrix* B, SLUMT::Gstat_t* Gstat,
int* info)
381 SLUMT::C::cgstrs(trans, L, U, perm_r, perm_c, B, Gstat, info);
384 static void gstrf(SLUMT::superlumt_options_t* options, SLUMT::SuperMatrix* A,
385 int* perm_r, SLUMT::SuperMatrix* L, SLUMT::SuperMatrix* U,
386 SLUMT::Gstat_t* stat,
int* info)
388 SLUMT::C::pcgstrf(options, A, perm_r, L, U, stat, info);
391 static void create_CompCol_Matrix(SLUMT::SuperMatrix* A,
int m,
int n,
int nnz,
392 type_map::type* nzval,
int* rowind,
int* colptr,
393 SLUMT::Stype_t stype, SLUMT::Dtype_t dtype, SLUMT::Mtype_t mtype)
395 SLUMT::C::cCreate_CompCol_Matrix(A, m, n, nnz, nzval, rowind, colptr,
396 stype, dtype, mtype);
399 static void create_Dense_Matrix(SLUMT::SuperMatrix* X,
int m,
int n,
400 type_map::type* x,
int ldx, SLUMT::Stype_t stype,
401 SLUMT::Dtype_t dtype, SLUMT::Mtype_t mtype)
403 SLUMT::C::cCreate_Dense_Matrix(X, m, n, x, ldx, stype, dtype, mtype);
406 static void gsequ(SLUMT::SuperMatrix* A,
float* r,
float* c,
407 float* rowcnd,
float* colcnd,
float* amax,
int* info)
409 SLUMT::C::cgsequ(A, r, c, rowcnd, colcnd, amax, info);
412 static void laqgs(SLUMT::SuperMatrix* A,
float* r,
float* c,
float rowcnd,
413 float colcnd,
float amax, SLUMT::equed_t* equed)
415 SLUMT::C::claqgs(A, r, c, rowcnd, colcnd, amax, equed);
421 struct FunctionMap<Superlumt,SLUMT::Z::doublecomplex>
423 typedef TypeMap<Superlumt,SLUMT::Z::doublecomplex> type_map;
425 static void gssvx(SLUMT::superlumt_options_t* options, SLUMT::SuperMatrix* A,
426 int* perm_c,
int* perm_r,
int* etree, SLUMT::equed_t* equed,
double* R,
double* C,
427 SLUMT::SuperMatrix* L, SLUMT::SuperMatrix* U,
void* work,
int lwork,
428 SLUMT::SuperMatrix* B, SLUMT::SuperMatrix* X,
double* recip_pivot_growth,
429 double* rcond,
double* ferr,
double* berr, SLUMT::superlu_memusage_t* mem_usage,
430 SLUMT::Gstat_t* stat,
int* info)
432 options->etree = etree;
433 options->perm_c = perm_c;
434 options->perm_r = perm_r;
436 options->work = work;
437 options->lwork = lwork;
439 SLUMT::Z::pzgssvx(options->nprocs, options, A, perm_c, perm_r,
440 equed, R, C, L, U, B, X, recip_pivot_growth, rcond, ferr,
441 berr, mem_usage, info);
444 static void gstrs(SLUMT::trans_t trans, SLUMT::SuperMatrix* L,
445 SLUMT::SuperMatrix* U,
int* perm_r,
int* perm_c,
446 SLUMT::SuperMatrix* B, SLUMT::Gstat_t* Gstat,
int* info)
448 SLUMT::Z::zgstrs(trans, L, U, perm_r, perm_c, B, Gstat, info);
451 static void gstrf(SLUMT::superlumt_options_t* options, SLUMT::SuperMatrix* A,
452 int* perm_r, SLUMT::SuperMatrix* L, SLUMT::SuperMatrix* U,
453 SLUMT::Gstat_t* stat,
int* info)
455 SLUMT::Z::pzgstrf(options, A, perm_r, L, U, stat, info);
458 static void create_CompCol_Matrix(SLUMT::SuperMatrix* A,
int m,
int n,
int nnz,
459 type_map::type* nzval,
int* rowind,
int* colptr,
460 SLUMT::Stype_t stype, SLUMT::Dtype_t dtype, SLUMT::Mtype_t mtype)
462 SLUMT::Z::zCreate_CompCol_Matrix(A, m, n, nnz, nzval, rowind, colptr,
463 stype, dtype, mtype);
466 static void create_Dense_Matrix(SLUMT::SuperMatrix* X,
int m,
int n,
467 type_map::type* x,
int ldx, SLUMT::Stype_t stype,
468 SLUMT::Dtype_t dtype, SLUMT::Mtype_t mtype)
470 SLUMT::Z::zCreate_Dense_Matrix(X, m, n, x, ldx, stype, dtype, mtype);
473 static void gsequ(SLUMT::SuperMatrix* A,
double* r,
double* c,
474 double* rowcnd,
double* colcnd,
double* amax,
int* info)
476 SLUMT::Z::zgsequ(A, r, c, rowcnd, colcnd, amax, info);
479 static void laqgs(SLUMT::SuperMatrix* A,
double* r,
double* c,
double rowcnd,
480 double colcnd,
double amax, SLUMT::equed_t* equed)
482 SLUMT::Z::zlaqgs(A, r, c, rowcnd, colcnd, amax, equed);
485 #endif // HAVE_TEUCHOS_COMPLEX
492 #endif // AMESOS2_SUPERLUMT_FUNCTIONMAP_HPP
Declaration of Function mapping class for Amesos2.
Provides definition of SuperLU_MT types as well as conversions and type traits.