suanPan
🧮 An Open Source, Parallel and Heterogeneous Finite Element Analysis Framework
Loading...
Searching...
No Matches
SparseMatSuperLU.hpp
Go to the documentation of this file.
1/*******************************************************************************
2 * Copyright (C) 2017-2026 Theodore Chang
3 *
4 * This program is free software: you can redistribute it and/or modify
5 * it under the terms of the GNU General Public License as published by
6 * the Free Software Foundation, either version 3 of the License, or
7 * (at your option) any later version.
8 *
9 * This program is distributed in the hope that it will be useful,
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 * GNU General Public License for more details.
13 *
14 * You should have received a copy of the GNU General Public License
15 * along with this program. If not, see <http://www.gnu.org/licenses/>.
16 ******************************************************************************/
29// ReSharper disable CppCStyleCast
30#ifndef SPARSEMATSUPERLU_HPP
31#define SPARSEMATSUPERLU_HPP
32
33#include "../SparseMat.hpp"
34#include "../csc_form.hpp"
35
36#include <superlu-mt/superlu-mt.h>
37
38template<sp_d T> class SparseMatSuperLU final : public SparseMat<T> {
39 SuperMatrix A{}, L{}, U{}, B{};
40
41#ifndef SUANPAN_SUPERLUMT
42 superlu_options_t options{};
43
44 SuperLUStat_t stat{};
45#else
46 const int ordering_num = 1;
47
48 Gstat_t stat{};
49#endif
50
51 std::vector<T> t_val;
52 std::vector<int> t_row, t_col, perm_r, perm_c;
53
54 bool allocated = false;
55
56 auto init_config();
57
58 template<sp_d ET> void alloc(csc_form<ET, int>&&);
59 void dealloc();
60
61 template<sp_d ET> void wrap_b(const Mat<ET>&);
62
63 int solve_trs(Mat<T>&, Mat<T>&&);
64
65protected:
66 int direct_solve(Mat<T>& out_mat, const Mat<T>& in_mat) override { return this->direct_solve(out_mat, Mat<T>(in_mat)); }
67
68 int direct_solve(Mat<T>&, Mat<T>&&) override;
69
70public:
71 SparseMatSuperLU(uword, uword, uword = 0);
76 ~SparseMatSuperLU() override;
77
78 unique_ptr<MetaMat<T>> unique_copy() override;
79};
80
81template<sp_d T> auto SparseMatSuperLU<T>::init_config() {
82#ifndef SUANPAN_SUPERLUMT
83 set_default_options(&options);
84 options.IterRefine = std::is_same_v<T, float> ? superlu::IterRefine_t::SLU_SINGLE : superlu::IterRefine_t::SLU_DOUBLE;
85 options.Equil = superlu::yes_no_t::NO;
86
87 StatInit(&stat);
88#else
89 StatAlloc(static_cast<int>(this->n_cols), SUANPAN_NUM_THREADS, sp_ienv(1), sp_ienv(2), &stat);
90 StatInit(static_cast<int>(this->n_cols), SUANPAN_NUM_THREADS, &stat);
91#endif
92
93 this->factored = false;
94}
95
96template<sp_d T> template<sp_d ET> void SparseMatSuperLU<T>::alloc(csc_form<ET, int>&& in) {
97 dealloc();
98
99 t_row = std::vector<int>(in.row_mem(), in.row_mem() + in.n_elem);
100 t_col = std::vector<int>(in.col_mem(), in.col_mem() + in.n_cols + 1);
101 t_val = std::vector<ET>(in.val_mem(), in.val_mem() + in.n_elem);
102
103 if constexpr(std::is_same_v<ET, double>) dCreate_CompCol_Matrix(&A, in.n_rows, in.n_cols, in.n_elem, t_val.data(), t_row.data(), t_col.data(), Stype_t::SLU_NC, Dtype_t::SLU_D, Mtype_t::SLU_GE);
104 else sCreate_CompCol_Matrix(&A, in.n_rows, in.n_cols, in.n_elem, t_val.data(), t_row.data(), t_col.data(), Stype_t::SLU_NC, Dtype_t::SLU_S, Mtype_t::SLU_GE);
105
106 perm_r.resize(this->n_rows + 1);
107 perm_c.resize(this->n_cols + 1);
108
109 allocated = true;
110}
111
112template<sp_d T> void SparseMatSuperLU<T>::dealloc() {
113 if(!allocated) return;
114
115 Destroy_SuperMatrix_Store(&A);
116#ifdef SUANPAN_SUPERLUMT
117 Destroy_SuperNode_SCP(&L);
118 Destroy_CompCol_NCP(&U);
119#else
120 Destroy_SuperNode_Matrix(&L);
121 Destroy_CompCol_Matrix(&U);
122#endif
123
124 allocated = false;
125}
126
127template<sp_d T> template<sp_d ET> void SparseMatSuperLU<T>::wrap_b(const Mat<ET>& in_mat) {
128 if constexpr(std::is_same_v<ET, float>) sCreate_Dense_Matrix(&B, (int)in_mat.n_rows, (int)in_mat.n_cols, (ET*)in_mat.memptr(), (int)in_mat.n_rows, Stype_t::SLU_DN, Dtype_t::SLU_S, Mtype_t::SLU_GE);
129 else dCreate_Dense_Matrix(&B, (int)in_mat.n_rows, (int)in_mat.n_cols, (ET*)in_mat.memptr(), (int)in_mat.n_rows, Stype_t::SLU_DN, Dtype_t::SLU_D, Mtype_t::SLU_GE);
130}
131
132template<sp_d T> SparseMatSuperLU<T>::SparseMatSuperLU(const uword in_row, const uword in_col, const uword in_elem)
133 : SparseMat<T>(in_row, in_col, in_elem) { init_config(); }
134
136 : SparseMat<T>(other) { init_config(); }
137
139 dealloc();
140 StatFree(&stat);
141}
142
143template<sp_d T> unique_ptr<MetaMat<T>> SparseMatSuperLU<T>::unique_copy() { return std::make_unique<SparseMatSuperLU>(*this); }
144
145template<sp_d T> int SparseMatSuperLU<T>::direct_solve(Mat<T>& out_mat, Mat<T>&& in_mat) {
146 if(this->factored) return solve_trs(out_mat, std::move(in_mat));
147
148 this->factored = true;
149
150 alloc(csc_form<T, int>(this->triplet_mat));
151
152 wrap_b(in_mat);
153
154 int flag{0};
155
156#ifdef SUANPAN_SUPERLUMT
157 get_perm_c(ordering_num, &A, perm_c.data());
158 if(std::is_same_v<T, float>) psgssv(SUANPAN_NUM_THREADS, &A, perm_c.data(), perm_r.data(), &L, &U, &B, &flag);
159 else pdgssv(SUANPAN_NUM_THREADS, &A, perm_c.data(), perm_r.data(), &L, &U, &B, &flag);
160#else
161 superlu::gssv<T>(&options, &A, perm_c.data(), perm_r.data(), &L, &U, &B, &stat, &flag);
162#endif
163
164 Destroy_SuperMatrix_Store(&B);
165
166 out_mat = std::move(in_mat);
167
168 return flag;
169}
170
171template<sp_d T> int SparseMatSuperLU<T>::solve_trs(Mat<T>& out_mat, Mat<T>&& in_mat) {
172 wrap_b(in_mat);
173
174 int flag{0};
175
176#ifdef SUANPAN_SUPERLUMT
177 if(std::is_same_v<T, float>) sgstrs(NOTRANS, &L, &U, perm_c.data(), perm_r.data(), &B, &stat, &flag);
178 else dgstrs(NOTRANS, &L, &U, perm_c.data(), perm_r.data(), &B, &stat, &flag);
179#else
180 superlu::gstrs<T>(options.Trans, &L, &U, perm_c.data(), perm_r.data(), &B, &stat, &flag);
181#endif
182
183 Destroy_SuperMatrix_Store(&B);
184
185 out_mat = std::move(in_mat);
186
187 return flag;
188}
189#endif
190
A SparseMat class that holds matrices.
Definition SparseMat.hpp:34
A SparseMatSuperLU class that holds matrices.
Definition SparseMatSuperLU.hpp:38
SparseMatSuperLU & operator=(const SparseMatSuperLU &)=delete
SparseMatSuperLU(SparseMatSuperLU &&)=delete
SparseMatSuperLU & operator=(SparseMatSuperLU &&)=delete
int direct_solve(Mat< T > &out_mat, const Mat< T > &in_mat) override
Definition SparseMatSuperLU.hpp:66
Definition csc_form.hpp:25
int SUANPAN_NUM_THREADS
Definition command.cpp:70
~SparseMatSuperLU() override
Definition SparseMatSuperLU.hpp:138
unique_ptr< MetaMat< T > > unique_copy() override
Definition SparseMatSuperLU.hpp:143
SparseMatSuperLU(uword, uword, uword=0)
Definition SparseMatSuperLU.hpp:132