suanPan
Loading...
Searching...
No Matches
SparseMatCUDA.hpp
Go to the documentation of this file.
1/*******************************************************************************
2 * Copyright (C) 2017-2023 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 ******************************************************************************/
31// ReSharper disable CppCStyleCast
32#ifndef SPARSEMATCUDA_HPP
33#define SPARSEMATCUDA_HPP
34
35#ifdef SUANPAN_CUDA
36
37#include <cusolverSp.h>
38#include <cusparse.h>
39#include "SparseMat.hpp"
40#include "csr_form.hpp"
41
42template<sp_d T> class SparseMatCUDA final : public SparseMat<T> {
43 cusolverSpHandle_t handle = nullptr;
44 cudaStream_t stream = nullptr;
45 cusparseMatDescr_t descr = nullptr;
46
47 void* d_val_idx = nullptr;
48 void* d_col_idx = nullptr;
49 void* d_row_ptr = nullptr;
50
51 void acquire();
52 void release() const;
53
54 template<sp_d ET> void device_alloc(csr_form<ET, int>&&);
55 void device_dealloc() const;
56
57public:
58 SparseMatCUDA(uword, uword, uword = 0);
60 SparseMatCUDA(SparseMatCUDA&&) noexcept = delete;
61 SparseMatCUDA& operator=(const SparseMatCUDA&) = delete;
62 SparseMatCUDA& operator=(SparseMatCUDA&&) noexcept = delete;
63 ~SparseMatCUDA() override;
64
65 unique_ptr<MetaMat<T>> make_copy() override;
66
67 int direct_solve(Mat<T>&, const Mat<T>&) override;
68};
69
70template<sp_d T> void SparseMatCUDA<T>::acquire() {
71 cusolverSpCreate(&handle);
72 cudaStreamCreate(&stream);
73 cusolverSpSetStream(handle, stream);
74 cusparseCreateMatDescr(&descr);
75 cusparseSetMatType(descr, CUSPARSE_MATRIX_TYPE_GENERAL);
76 cusparseSetMatIndexBase(descr, CUSPARSE_INDEX_BASE_ZERO);
77 cudaMalloc(&d_row_ptr, sizeof(int) * (this->n_rows + 1));
78}
79
80template<sp_d T> void SparseMatCUDA<T>::release() const {
81 if(handle) cusolverSpDestroy(handle);
82 if(stream) cudaStreamDestroy(stream);
83 if(descr) cusparseDestroyMatDescr(descr);
84 if(d_row_ptr) cudaFree(d_row_ptr);
85}
86
87template<sp_d T> template<sp_d ET> void SparseMatCUDA<T>::device_alloc(csr_form<ET, int>&& csr_mat) {
88 const size_t n_val = sizeof(ET) * csr_mat.n_elem;
89 const size_t n_col = sizeof(int) * csr_mat.n_elem;
90
91 cudaMalloc(&d_val_idx, n_val);
92 cudaMalloc(&d_col_idx, n_col);
93
94 cudaMemcpyAsync(d_val_idx, csr_mat.val_mem(), n_val, cudaMemcpyHostToDevice, stream);
95 cudaMemcpyAsync(d_col_idx, csr_mat.col_mem(), n_col, cudaMemcpyHostToDevice, stream);
96 cudaMemcpyAsync(d_row_ptr, csr_mat.row_mem(), sizeof(int) * (csr_mat.n_rows + 1llu), cudaMemcpyHostToDevice, stream);
97}
98
99template<sp_d T> void SparseMatCUDA<T>::device_dealloc() const {
100 if(d_val_idx) cudaFree(d_val_idx);
101 if(d_col_idx) cudaFree(d_col_idx);
102}
103
104template<sp_d T> SparseMatCUDA<T>::SparseMatCUDA(const uword in_row, const uword in_col, const uword in_elem)
105 : SparseMat<T>(in_row, in_col, in_elem) { acquire(); }
106
107template<sp_d T> SparseMatCUDA<T>::SparseMatCUDA(const SparseMatCUDA& other)
108 : SparseMat<T>(other) { acquire(); }
109
110template<sp_d T> SparseMatCUDA<T>::~SparseMatCUDA() {
111 release();
112 device_dealloc();
113}
114
115template<sp_d T> unique_ptr<MetaMat<T>> SparseMatCUDA<T>::make_copy() { return std::make_unique<SparseMatCUDA<T>>(*this); }
116
117template<sp_d T> int SparseMatCUDA<T>::direct_solve(Mat<T>& X, const Mat<T>& B) {
118 if(!this->factored) {
119 // deallocate memory previously allocated for csr matrix
120 device_dealloc();
121
122 std::is_same_v<T, float> || Precision::MIXED == this->setting.precision ? device_alloc(csr_form<float, int>(this->triplet_mat)) : device_alloc(csr_form<double, int>(this->triplet_mat));
123
124 this->factored = true;
125 }
126
127 const size_t n_rhs = (std::is_same_v<T, float> || Precision::MIXED == this->setting.precision ? sizeof(float) : sizeof(double)) * B.n_elem;
128
129 void* d_b = nullptr;
130 void* d_x = nullptr;
131
132 cudaMalloc(&d_b, n_rhs);
133 cudaMalloc(&d_x, n_rhs);
134
135 int singularity;
136 auto code = 0;
137
138 if(std::is_same_v<T, float>) {
139 cudaMemcpyAsync(d_b, B.memptr(), n_rhs, cudaMemcpyHostToDevice, stream);
140
141 for(auto I = 0llu; I < B.n_elem; I += B.n_rows) code += cusolverSpScsrlsvqr(handle, int(this->n_rows), int(this->triplet_mat.n_elem), descr, (float*)d_val_idx, (int*)d_row_ptr, (int*)d_col_idx, (float*)d_b + I, float(this->setting.tolerance), 3, (float*)d_x + I, &singularity);
142
143 X.set_size(arma::size(B));
144
145 cudaMemcpyAsync(X.memptr(), d_x, n_rhs, cudaMemcpyDeviceToHost, stream);
146
147 cudaDeviceSynchronize();
148 }
149 else if(Precision::FULL == this->setting.precision) {
150 cudaMemcpyAsync(d_b, B.memptr(), n_rhs, cudaMemcpyHostToDevice, stream);
151
152 for(auto I = 0llu; I < B.n_elem; I += B.n_rows) code += cusolverSpDcsrlsvqr(handle, int(this->n_rows), int(this->triplet_mat.n_elem), descr, (double*)d_val_idx, (int*)d_row_ptr, (int*)d_col_idx, (double*)d_b + I, this->setting.tolerance, 3, (double*)d_x + I, &singularity);
153
154 X.set_size(arma::size(B));
155
156 cudaMemcpyAsync(X.memptr(), d_x, n_rhs, cudaMemcpyDeviceToHost, stream);
157
158 cudaDeviceSynchronize();
159 }
160 else {
161 X = arma::zeros(B.n_rows, B.n_cols);
162
163 mat full_residual = B;
164
165 auto multiplier = norm(full_residual);
166
167 auto counter = 0u;
168 while(counter++ < this->setting.iterative_refinement) {
169 if(multiplier < this->setting.tolerance) break;
170
171 auto residual = conv_to<fmat>::from(full_residual / multiplier);
172
173 cudaMemcpyAsync(d_b, residual.memptr(), n_rhs, cudaMemcpyHostToDevice, stream);
174
175 code = 0;
176 for(auto I = 0llu; I < B.n_elem; I += B.n_rows) code += cusolverSpScsrlsvqr(handle, int(this->n_rows), int(this->triplet_mat.n_elem), descr, (float*)d_val_idx, (int*)d_row_ptr, (int*)d_col_idx, (float*)d_b + I, float(this->setting.tolerance), 3, (float*)d_x + I, &singularity);
177 if(0 != code) break;
178
179 cudaMemcpyAsync(residual.memptr(), d_x, n_rhs, cudaMemcpyDeviceToHost, stream);
180
181 cudaDeviceSynchronize();
182
183 const mat incre = multiplier * conv_to<mat>::from(residual);
184
185 X += incre;
186
187 suanpan_debug("Mixed precision algorithm multiplier: {:.5E}.\n", multiplier = arma::norm(full_residual -= this->operator*(incre)));
188 }
189 }
190
191 if(d_b) cudaFree(d_b);
192 if(d_x) cudaFree(d_x);
193
194 return 0 == code ? SUANPAN_SUCCESS : SUANPAN_FAIL;
195}
196
197#endif
198
199#endif
200
A MetaMat class that holds matrices.
Definition: MetaMat.hpp:39
A SparseMatCUDA class that holds matrices.
A SparseMat class that holds matrices.
Definition: SparseMat.hpp:34
Definition: csr_form.hpp:25
Definition: suanPan.h:318
double norm(const vec &)
Definition: tensor.cpp:302
#define suanpan_debug(...)
Definition: suanPan.h:295
constexpr auto SUANPAN_SUCCESS
Definition: suanPan.h:162
constexpr auto SUANPAN_FAIL
Definition: suanPan.h:163