3 #ifndef DUNE_ISTL_UMFPACK_HH 4 #define DUNE_ISTL_UMFPACK_HH 6 #if HAVE_SUITESPARSE_UMFPACK || defined DOXYGEN 13 #include<dune/common/exceptions.hh> 14 #include<dune/common/fmatrix.hh> 15 #include<dune/common/fvector.hh> 16 #include<dune/common/unused.hh> 37 template<
class M,
class T,
class TM,
class TD,
class TA>
38 class SeqOverlappingSchwarz;
40 template<
class T,
bool tag>
41 struct SeqOverlappingSchwarzAssemblerHelper;
48 template<
class Matrix>
61 template<
typename... A>
64 umfpack_di_defaults(args...);
66 template<
typename... A>
69 umfpack_di_free_numeric(args...);
71 template<
typename... A>
74 umfpack_di_free_symbolic(args...);
76 template<
typename... A>
79 return umfpack_di_load_numeric(args...);
81 template<
typename... A>
84 umfpack_di_numeric(args...);
86 template<
typename... A>
89 umfpack_di_report_info(args...);
91 template<
typename... A>
94 umfpack_di_report_status(args...);
96 template<
typename... A>
99 return umfpack_di_save_numeric(args...);
101 template<
typename... A>
104 umfpack_di_solve(args...);
106 template<
typename... A>
109 umfpack_di_symbolic(args...);
116 template<
typename... A>
119 umfpack_zi_defaults(args...);
121 template<
typename... A>
124 umfpack_zi_free_numeric(args...);
126 template<
typename... A>
129 umfpack_zi_free_symbolic(args...);
131 template<
typename... A>
134 return umfpack_zi_load_numeric(args...);
136 template<
typename... A>
137 static void numeric(
const int* cs,
const int* ri,
const double* val, A... args)
139 umfpack_zi_numeric(cs,ri,val,NULL,args...);
141 template<
typename... A>
144 umfpack_zi_report_info(args...);
146 template<
typename... A>
149 umfpack_zi_report_status(args...);
151 template<
typename... A>
154 return umfpack_zi_save_numeric(args...);
156 template<
typename... A>
157 static void solve(
int m,
const int* cs,
const int* ri, std::complex<double>* val,
double* x,
const double* b,A... args)
159 const double* cval =
reinterpret_cast<const double*
>(val);
160 umfpack_zi_solve(m,cs,ri,cval,NULL,x,NULL,b,NULL,args...);
162 template<
typename... A>
163 static void symbolic(
int m,
int n,
const int* cs,
const int* ri,
const double* val, A... args)
165 umfpack_zi_symbolic(m,n,cs,ri,val,NULL,args...);
182 template<
typename T,
typename A,
int n,
int m>
185 BlockVector<FieldVector<T,m>,
186 typename A::template rebind<FieldVector<T,m> >::other>,
187 BlockVector<FieldVector<T,n>,
188 typename A::template rebind<FieldVector<T,n> >::other> >
201 typename A::template rebind<FieldVector<T,m> >::other>
domain_type;
205 typename A::template rebind<FieldVector<T,n> >::other>
range_type;
215 UMFPack(
const Matrix& matrix,
int verbose=0) : matrixIsLoaded_(false)
218 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
219 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
220 Caller::defaults(UMF_Control);
221 setVerbosity(verbose);
233 UMFPack(
const Matrix& matrix,
int verbose,
bool) : matrixIsLoaded_(false)
236 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
237 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
238 Caller::defaults(UMF_Control);
239 setVerbosity(verbose);
245 UMFPack() : matrixIsLoaded_(false), verbosity_(0)
248 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
249 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
250 Caller::defaults(UMF_Control);
263 UMFPack(
const Matrix& mat_,
const char* file,
int verbose=0)
266 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
267 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
268 Caller::defaults(UMF_Control);
269 setVerbosity(verbose);
270 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
271 if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
273 matrixIsLoaded_ =
false;
275 saveDecomposition(file);
279 matrixIsLoaded_ =
true;
280 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
293 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
294 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
295 Caller::defaults(UMF_Control);
296 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
297 if (errcode == UMFPACK_ERROR_out_of_memory)
298 DUNE_THROW(Dune::Exception,
"ran out of memory while loading UMFPack decomposition");
299 if (errcode == UMFPACK_ERROR_file_IO)
300 DUNE_THROW(Dune::Exception,
"IO error while loading UMFPack decomposition");
301 matrixIsLoaded_ =
true;
302 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
303 setVerbosity(verbose);
308 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
317 if (umfpackMatrix_.N() != b.
dim())
318 DUNE_THROW(
Dune::ISTLError,
"Size of right-hand-side vector b does not match the number of matrix rows!");
319 if (umfpackMatrix_.M() != x.
dim())
320 DUNE_THROW(
Dune::ISTLError,
"Size of solution vector x does not match the number of matrix columns!");
322 double UMF_Apply_Info[UMFPACK_INFO];
323 Caller::solve(UMFPACK_A,
324 umfpackMatrix_.getColStart(),
325 umfpackMatrix_.getRowIndex(),
326 umfpackMatrix_.getValues(),
327 reinterpret_cast<double*
>(&x[0]),
328 reinterpret_cast<double*>(&b[0]),
336 res.
elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];
338 printOnApply(UMF_Apply_Info);
346 DUNE_UNUSED_PARAMETER(reduction);
357 double UMF_Apply_Info[UMFPACK_INFO];
358 Caller::solve(UMFPACK_A,
359 umfpackMatrix_.getColStart(),
360 umfpackMatrix_.getRowIndex(),
361 umfpackMatrix_.getValues(),
367 printOnApply(UMF_Apply_Info);
383 if (option >= UMFPACK_CONTROL)
384 DUNE_THROW(RangeError,
"Requested non-existing UMFPack option");
386 UMF_Control[option] = value;
394 int errcode = Caller::save_numeric(UMF_Numeric, const_cast<char*>(file));
395 if (errcode != UMFPACK_OK)
396 DUNE_THROW(Dune::Exception,
"IO ERROR while trying to save UMFPack decomposition");
402 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
404 umfpackMatrix_ = matrix;
411 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
413 umfpackMatrix_.setMatrix(_mat,rowIndexSet);
429 UMF_Control[UMFPACK_PRL] = 1;
431 UMF_Control[UMFPACK_PRL] = 2;
433 UMF_Control[UMFPACK_PRL] = 4;
451 return umfpackMatrix_;
460 if (!matrixIsLoaded_)
462 Caller::free_symbolic(&UMF_Symbolic);
463 umfpackMatrix_.free();
465 Caller::free_numeric(&UMF_Numeric);
466 matrixIsLoaded_ =
false;
469 const char*
name() {
return "UMFPACK"; }
474 template<
class M,
class X,
class TM,
class TD,
class T1>
481 double UMF_Decomposition_Info[UMFPACK_INFO];
482 Caller::symbolic(static_cast<int>(umfpackMatrix_.N()),
483 static_cast<int>(umfpackMatrix_.N()),
484 umfpackMatrix_.getColStart(),
485 umfpackMatrix_.getRowIndex(),
486 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
489 UMF_Decomposition_Info);
490 Caller::numeric(umfpackMatrix_.getColStart(),
491 umfpackMatrix_.getRowIndex(),
492 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
496 UMF_Decomposition_Info);
497 Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
500 std::cout <<
"[UMFPack Decomposition]" << std::endl;
501 std::cout <<
"Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] <<
" (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] <<
")" << std::endl;
502 std::cout <<
"Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
503 std::cout <<
"Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] <<
" bytes" << std::endl;
504 std::cout <<
"Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
505 std::cout <<
"Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] <<
" U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
509 Caller::report_info(UMF_Control,UMF_Decomposition_Info);
513 void printOnApply(
double* UMF_Info)
515 Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
518 std::cout <<
"[UMFPack Solve]" << std::endl;
519 std::cout <<
"Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] <<
" (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] <<
")" << std::endl;
520 std::cout <<
"Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
521 std::cout <<
"Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
522 std::cout <<
"Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] <<
" resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
526 UMFPackMatrix umfpackMatrix_;
527 bool matrixIsLoaded_;
531 double UMF_Control[UMFPACK_CONTROL];
534 template<
typename T,
typename A,
int n,
int m>
540 template<
typename T,
typename A,
int n,
int m>
543 enum { value =
true };
547 #endif // HAVE_SUITESPARSE_UMFPACK 549 #endif //DUNE_ISTL_UMFPACK_HH Definition: solvertype.hh:27
UMFPack()
default constructor
Definition: umfpack.hh:245
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: umfpack.hh:315
Use the UMFPack package to directly solve linear systems – empty default class.
Definition: umfpack.hh:49
UMFPack(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: umfpack.hh:215
void apply(T *x, T *b)
additional apply method with c-arrays in analogy to superlu
Definition: umfpack.hh:355
UMFPackMatrix & getInternalMatrix()
Return the column compress matrix from UMFPack.
Definition: umfpack.hh:449
static int load_numeric(A... args)
Definition: umfpack.hh:77
Statistics about the application of an inverse operator.
Definition: solver.hh:31
Dune::BlockVector< FieldVector< T, m >, typename A::template rebind< FieldVector< T, m > >::other > domain_type
The type of the domain of the solver.
Definition: umfpack.hh:201
static void free_numeric(A... args)
Definition: umfpack.hh:67
static void free_symbolic(A... args)
Definition: umfpack.hh:72
double elapsed
Elapsed time in seconds.
Definition: solver.hh:62
Implementations of the inverse operator interface.
Inititializer for the ColCompMatrix as needed by OverlappingSchwarz.
Definition: colcompmatrix.hh:153
static int save_numeric(A... args)
Definition: umfpack.hh:152
static void report_info(A... args)
Definition: umfpack.hh:142
static void report_status(A... args)
Definition: umfpack.hh:92
size_type dim() const
dimension of the vector space
Definition: bvector.hh:279
Templates characterizing the type of a solver.
static void numeric(const int *cs, const int *ri, const double *val, A... args)
Definition: umfpack.hh:137
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:344
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:392
Definition: colcompmatrix.hh:160
Definition: solvertype.hh:13
const char * name()
Definition: umfpack.hh:469
Dune::BlockVector< FieldVector< T, n >, typename A::template rebind< FieldVector< T, n > >::other > range_type
The type of the range of the solver.
Definition: umfpack.hh:205
static void free_numeric(A... args)
Definition: umfpack.hh:122
virtual ~UMFPack()
Definition: umfpack.hh:306
Definition: matrixutils.hh:25
UMFPack(const Matrix &mat_, const char *file, int verbose=0)
Try loading a decomposition from file and do a decomposition if unsuccessful.
Definition: umfpack.hh:263
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:411
static void numeric(A... args)
Definition: umfpack.hh:82
int iterations
Number of iterations.
Definition: solver.hh:50
Dune::ColCompMatrix< Matrix > UMFPackMatrix
The corresponding SuperLU Matrix type.
Definition: umfpack.hh:195
Definition: basearray.hh:19
static void report_info(A... args)
Definition: umfpack.hh:87
static void defaults(A... args)
Definition: umfpack.hh:62
bool converged
True if convergence criterion has been met.
Definition: solver.hh:56
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: umfpack.hh:400
static void solve(A... args)
Definition: umfpack.hh:102
Implementation of the BCRSMatrix class.
A vector of blocks with memory management.
Definition: bvector.hh:308
void * getFactorization()
Return the matrix factorization.
Definition: umfpack.hh:440
static int save_numeric(A... args)
Definition: umfpack.hh:97
Definition: umfpack.hh:55
void setSubMatrix(const Matrix &_mat, const S &rowIndexSet)
Definition: umfpack.hh:409
static void solve(int m, const int *cs, const int *ri, std::complex< double > *val, double *x, const double *b, A... args)
Definition: umfpack.hh:157
void free()
free allocated space.
Definition: umfpack.hh:458
static void report_status(A... args)
Definition: umfpack.hh:147
static void symbolic(int m, int n, const int *cs, const int *ri, const double *val, A... args)
Definition: umfpack.hh:163
Abstract base class for all solvers.
Definition: solver.hh:79
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:290
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:424
UMFPack(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:233
Definition: colcompmatrix.hh:250
derive error class from the base class in common
Definition: istlexception.hh:16
Sequential overlapping Schwarz preconditioner.
Definition: colcompmatrix.hh:157
static int load_numeric(A... args)
Definition: umfpack.hh:132
static void symbolic(A... args)
Definition: umfpack.hh:107
static void defaults(A... args)
Definition: umfpack.hh:117
static void free_symbolic(A... args)
Definition: umfpack.hh:127
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:381