StrumpackSparseSolverMixedPrecision.hpp
1 /*
2  * STRUMPACK -- STRUctured Matrices PACKage, Copyright (c) 2014, The
3  * Regents of the University of California, through Lawrence Berkeley
4  * National Laboratory (subject to receipt of any required approvals
5  * from the U.S. Dept. of Energy). All rights reserved.
6  *
7  * If you have questions about your rights to use or distribute this
8  * software, please contact Berkeley Lab's Technology Transfer
9  * Department at TTD@lbl.gov.
10  *
11  * NOTICE. This software is owned by the U.S. Department of Energy. As
12  * such, the U.S. Government has been granted for itself and others
13  * acting on its behalf a paid-up, nonexclusive, irrevocable,
14  * worldwide license in the Software to reproduce, prepare derivative
15  * works, and perform publicly and display publicly. Beginning five
16  * (5) years after the date permission to assert copyright is obtained
17  * from the U.S. Department of Energy, and subject to any subsequent
18  * five (5) year renewals, the U.S. Government is granted for itself
19  * and others acting on its behalf a paid-up, nonexclusive,
20  * irrevocable, worldwide license in the Software to reproduce,
21  * prepare derivative works, distribute copies to the public, perform
22  * publicly and display publicly, and to permit others to do so.
23  *
24  * Developers: Pieter Ghysels, Francois-Henry Rouet, Xiaoye S. Li.
25  * (Lawrence Berkeley National Lab, Computational Research
26  * Division).
27  */
33 #ifndef SPARSE_SOLVER_MIXED_PRECISION_HPP
34 #define SPARSE_SOLVER_MIXED_PRECISION_HPP
35 
36 #include <new>
37 #include <memory>
38 #include <vector>
39 #include <string>
40 
41 #include "iterative/IterativeSolvers.hpp"
42 #include "StrumpackOptions.hpp"
43 #include "StrumpackSparseSolver.hpp"
44 
45 
46 namespace strumpack {
47 
84  template<typename factor_t,typename refine_t,typename integer_t>
86 
87  public:
91  SparseSolverMixedPrecision(bool verbose=true, bool root=true);
92 
103  SparseSolverMixedPrecision(int argc, char* argv[],
104  bool verbose=true, bool root=true);
105 
110 
111  void set_matrix(const CSRMatrix<refine_t,integer_t>& A);
112  void set_matrix(const CSRMatrix<factor_t,integer_t>& A);
113 
114  ReturnCode factor();
115  ReturnCode reorder(int nx=1, int ny=1, int nz=1);
116 
117  ReturnCode solve(const DenseMatrix<refine_t>& b,
119  bool use_initial_guess=false);
120  ReturnCode solve(const refine_t* b, refine_t* x,
121  bool use_initial_guess=false);
122 
123  ReturnCode solve(const DenseMatrix<factor_t>& b,
125  bool use_initial_guess=false);
126  ReturnCode solve(const factor_t* b, factor_t* x,
127  bool use_initial_guess=false);
128 
129  SPOptions<refine_t>& options() { return opts_; }
130  const SPOptions<refine_t>& options() const { return opts_; }
131 
132  SparseSolver<factor_t,integer_t>& solver() { return solver_; }
133  const SparseSolver<factor_t,integer_t>& solver() const { return solver_; }
134 
139  int Krylov_iterations() const { return Krylov_its_; }
140 
141  private:
144  SPOptions<refine_t> opts_;
145  int Krylov_its_ = 0;
146  };
147 
148  template<typename factor_t,typename refine_t,typename integer_t>
149  using StrumpackSparseSolverMixedPrecision =
150  SparseSolverMixedPrecision<factor_t,refine_t,integer_t>;
151 
152 } //end namespace strumpack
153 
154 #endif // STRUMPACK_SPARSE_SOLVER_MIXED_PRECISION_HPP
Holds options for the sparse solver.
This class represents a matrix, stored in column major format, to allow direct use of BLAS/LAPACK rou...
Definition: DenseMatrix.hpp:138
SparseSolverMixedPrecision Allows to use lower precision (float) for the factorization/preconditioner...
Definition: StrumpackSparseSolverMixedPrecision.hpp:85
SparseSolverMixedPrecision(int argc, char *argv[], bool verbose=true, bool root=true)
SparseSolverMixedPrecision(bool verbose=true, bool root=true)
int Krylov_iterations() const
Definition: StrumpackSparseSolverMixedPrecision.hpp:139
Definition: StrumpackOptions.hpp:42
ReturnCode
Enumeration for the possible return codes.
Definition: StrumpackParameters.hpp:60