Loading...
Searching...
No Matches
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
46namespace 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 void set_lower_triangle_matrix(const CSRMatrix<refine_t,integer_t>& A);
115 void set_lower_triangle_matrix(const CSRMatrix<factor_t,integer_t>& A);
116
117 void update_matrix_values(const CSRMatrix<refine_t,integer_t>& A);
118 void update_matrix_values(const CSRMatrix<factor_t,integer_t>& A);
119
120 ReturnCode factor();
121 ReturnCode reorder(int nx=1, int ny=1, int nz=1);
122
123 ReturnCode solve(const DenseMatrix<refine_t>& b,
125 bool use_initial_guess=false);
126 ReturnCode solve(const refine_t* b, refine_t* x,
127 bool use_initial_guess=false);
128 ReturnCode solve(int nrhs, const refine_t* b, int ldb,
129 refine_t* x, int ldx,
130 bool use_initial_guess=false);
131
132 ReturnCode solve(const DenseMatrix<factor_t>& b,
134 bool use_initial_guess=false);
135 ReturnCode solve(const factor_t* b, factor_t* x,
136 bool use_initial_guess=false);
137 ReturnCode solve(int nrhs, const factor_t* b, int ldb,
138 factor_t* x, int ldx,
139 bool use_initial_guess=false);
140
141 SPOptions<refine_t>& options() { return opts_; }
142 const SPOptions<refine_t>& options() const { return opts_; }
143
144 SparseSolver<factor_t,integer_t>& solver() { return solver_; }
145 const SparseSolver<factor_t,integer_t>& solver() const { return solver_; }
146
151 int Krylov_iterations() const { return Krylov_its_; }
152
153 private:
157 int Krylov_its_ = 0;
158 };
159
160 template<typename factor_t,typename refine_t,typename integer_t>
161 using StrumpackSparseSolverMixedPrecision =
162 SparseSolverMixedPrecision<factor_t,refine_t,integer_t>;
163
164} //end namespace strumpack
165
166#endif // STRUMPACK_SPARSE_SOLVER_MIXED_PRECISION_HPP
Holds options for the sparse solver.
Class for storing a compressed sparse row matrix (single node).
Definition CSRMatrix.hpp:55
This class represents a matrix, stored in column major format, to allow direct use of BLAS/LAPACK rou...
Definition DenseMatrix.hpp:139
Options for the sparse solver.
Definition StrumpackOptions.hpp:217
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:151
SparseSolver is the main sequential or multithreaded sparse solver class.
Definition StrumpackSparseSolver.hpp:75
Definition StrumpackOptions.hpp:44
ReturnCode
Enumeration for the possible return codes.
Definition StrumpackParameters.hpp:50