xref: /petsc/src/mat/impls/dense/mpi/cupm/hip/matmpidensehip.hip.cxx (revision daba9d70159ea2f6905738fcbec7404635487b2b)
1 #include "../matmpidensecupm.hpp"
2 
3 using namespace Petsc::mat::cupm;
4 using Petsc::device::cupm::DeviceType;
5 
6 static constexpr impl::MatDense_MPI_CUPM<DeviceType::HIP> mat_cupm{};
7 
8 /*MC
9   MATDENSEHIP - "densehip" - A matrix type to be used for dense matrices on GPUs.
10 
11   This matrix type is identical to `MATSEQDENSEHIP` when constructed with a single process
12   communicator, and `MATMPIDENSEHIP` otherwise.
13 
14   Options Database Key:
15 . -mat_type densehip - sets the matrix type to `MATDENSEHIP` during a call to
16                         `MatSetFromOptions()`
17 
18   Level: beginner
19 
20 .seealso: [](ch_matrices), `Mat`, `MATSEQDENSEHIP`, `MATMPIDENSEHIP`, `MATSEQDENSECUDA`,
21 `MATMPIDENSECUDA`, `MATDENSE`
22 M*/
23 
24 /*MC
25   MATMPIDENSEHIP - "mpidensehip" - A matrix type to be used for distributed dense matrices on
26   GPUs.
27 
28   Options Database Key:
29 . -mat_type mpidensehip - sets the matrix type to `MATMPIDENSEHIP` during a call to
30                            `MatSetFromOptions()`
31 
32   Level: beginner
33 
34 .seealso: [](ch_matrices), `Mat`, `MATDENSEHIP`, `MATMPIDENSE`, `MATSEQDENSE`,
35 `MATSEQDENSEHIP`, `MATSEQDENSECUDA`
36 M*/
MatCreate_MPIDenseHIP(Mat A)37 PETSC_INTERN PetscErrorCode MatCreate_MPIDenseHIP(Mat A)
38 {
39   PetscFunctionBegin;
40   PetscCall(mat_cupm.Create(A));
41   PetscFunctionReturn(PETSC_SUCCESS);
42 }
43 
MatConvert_MPIDense_MPIDenseHIP(Mat A,MatType type,MatReuse reuse,Mat * ret)44 PetscErrorCode MatConvert_MPIDense_MPIDenseHIP(Mat A, MatType type, MatReuse reuse, Mat *ret)
45 {
46   PetscFunctionBegin;
47   PetscCall(mat_cupm.Convert_MPIDense_MPIDenseCUPM(A, type, reuse, ret));
48   PetscFunctionReturn(PETSC_SUCCESS);
49 }
50 
51 /*@C
52   MatCreateDenseHIP - Creates a matrix in `MATDENSEHIP` format using HIP.
53 
54   Collective
55 
56   Input Parameters:
57 + comm - MPI communicator
58 . m    - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
59 . n    - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given)
60 . M    - number of global rows (or `PETSC_DECIDE` to have calculated if `m` is given)
61 . N    - number of global columns (or `PETSC_DECIDE` to have calculated if `n` is given)
62 - data - optional location of GPU matrix data. Pass `NULL` to have PETSc to control matrix
63          memory allocation.
64 
65   Output Parameter:
66 . A - the matrix
67 
68   Level: intermediate
69 
70 .seealso: `MATDENSEHIP`, `MatCreate()`, `MatCreateDense()`
71 @*/
MatCreateDenseHIP(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar * data,Mat * A)72 PetscErrorCode MatCreateDenseHIP(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscScalar *data, Mat *A)
73 {
74   PetscFunctionBegin;
75   PetscCall(MatCreateDenseCUPM<DeviceType::HIP>(comm, m, n, M, N, data, A));
76   PetscFunctionReturn(PETSC_SUCCESS);
77 }
78 
79 /*@C
80   MatDenseHIPPlaceArray - Allows one to replace the GPU array in a `MATDENSEHIP` matrix with an
81   array provided by the user. This is useful to avoid copying an array into a matrix.
82 
83   Not Collective
84 
85   Input Parameters:
86 + mat   - the matrix
87 - array - the array in column major order
88 
89   Level: developer
90 
91   Note:
92   Adding `const` to `array` was an oversight, see notes in `VecPlaceArray()`.
93 
94   You can return to the original array with a call to `MatDenseHIPResetArray()`. The user is
95   responsible for freeing this array; it will not be freed when the matrix is destroyed. The
96   array must have been allocated with `hipMalloc()`.
97 
98 .seealso: `MATDENSEHIP`, `MatDenseHIPGetArray()`, `MatDenseHIPResetArray()`,
99 `MatDenseHIPReplaceArray()`
100 @*/
MatDenseHIPPlaceArray(Mat mat,const PetscScalar * array)101 PetscErrorCode MatDenseHIPPlaceArray(Mat mat, const PetscScalar *array)
102 {
103   PetscFunctionBegin;
104   PetscCall(MatDenseCUPMPlaceArray<DeviceType::HIP>(mat, array));
105   PetscFunctionReturn(PETSC_SUCCESS);
106 }
107 
108 /*@C
109   MatDenseHIPResetArray - Resets the matrix array to that it previously had before the call to
110   `MatDenseHIPPlaceArray()`
111 
112   Not Collective
113 
114   Input Parameter:
115 . mat - the matrix
116 
117   Level: developer
118 
119   Note:
120   You can only call this after a call to `MatDenseHIPPlaceArray()`
121 
122 .seealso: `MATDENSEHIP`, `MatDenseHIPGetArray()`, `MatDenseHIPPlaceArray()`
123 @*/
MatDenseHIPResetArray(Mat mat)124 PetscErrorCode MatDenseHIPResetArray(Mat mat)
125 {
126   PetscFunctionBegin;
127   PetscCall(MatDenseCUPMResetArray<DeviceType::HIP>(mat));
128   PetscFunctionReturn(PETSC_SUCCESS);
129 }
130 
131 /*@C
132   MatDenseHIPReplaceArray - Allows one to replace the GPU array in a `MATDENSEHIP` matrix
133   with an array provided by the user. This is useful to avoid copying an array into a matrix.
134 
135   Not Collective
136 
137   Input Parameters:
138 + mat   - the matrix
139 - array - the array in column major order
140 
141   Level: developer
142 
143   Note:
144   Adding `const` to `array` was an oversight, see notes in `VecPlaceArray()`.
145 
146   This permanently replaces the GPU array and frees the memory associated with the old GPU
147   array. The memory passed in CANNOT be freed by the user. It will be freed when the matrix is
148   destroyed. The array should respect the matrix leading dimension.
149 
150 .seealso: `MatDenseHIPGetArray()`, `MatDenseHIPPlaceArray()`, `MatDenseHIPResetArray()`
151 @*/
MatDenseHIPReplaceArray(Mat mat,const PetscScalar * array)152 PetscErrorCode MatDenseHIPReplaceArray(Mat mat, const PetscScalar *array)
153 {
154   PetscFunctionBegin;
155   PetscCall(MatDenseCUPMReplaceArray<DeviceType::HIP>(mat, array));
156   PetscFunctionReturn(PETSC_SUCCESS);
157 }
158 
159 /*@C
160   MatDenseHIPGetArrayWrite - Provides write access to the HIP buffer inside a `MATDENSEHIP`
161   matrix.
162 
163   Not Collective
164 
165   Input Parameter:
166 . A - the matrix
167 
168   Output Parameter:
169 . a - the GPU array in column major order
170 
171   Level: developer
172 
173   Notes:
174   The data on the GPU may not be updated due to operations done on the CPU. If you need updated
175   data, use `MatDenseHIPGetArray()`.
176 
177   The array must be restored with `MatDenseHIPRestoreArrayWrite()` when no longer needed.
178 
179 .seealso: `MATDENSEHIP`, `MatDenseHIPGetArray()`, `MatDenseHIPRestoreArray()`,
180 `MatDenseHIPRestoreArrayWrite()`, `MatDenseHIPGetArrayRead()`,
181 `MatDenseHIPRestoreArrayRead()`
182 @*/
MatDenseHIPGetArrayWrite(Mat A,PetscScalar ** a)183 PetscErrorCode MatDenseHIPGetArrayWrite(Mat A, PetscScalar **a)
184 {
185   PetscFunctionBegin;
186   PetscCall(MatDenseCUPMGetArrayWrite<DeviceType::HIP>(A, a));
187   PetscFunctionReturn(PETSC_SUCCESS);
188 }
189 
190 /*@C
191   MatDenseHIPRestoreArrayWrite - Restore write access to the HIP buffer inside a
192   `MATDENSEHIP` matrix previously obtained with `MatDenseHIPGetArrayWrite()`.
193 
194   Not Collective
195 
196   Input Parameters:
197 + A - the matrix
198 - a - the GPU array in column major order
199 
200   Level: developer
201 
202 .seealso: `MATDENSEHIP`, `MatDenseHIPGetArray()`, `MatDenseHIPRestoreArray()`,
203 `MatDenseHIPGetArrayWrite()`, `MatDenseHIPRestoreArrayRead()`, `MatDenseHIPGetArrayRead()`
204 @*/
MatDenseHIPRestoreArrayWrite(Mat A,PetscScalar ** a)205 PetscErrorCode MatDenseHIPRestoreArrayWrite(Mat A, PetscScalar **a)
206 {
207   PetscFunctionBegin;
208   PetscCall(MatDenseCUPMRestoreArrayWrite<DeviceType::HIP>(A, a));
209   PetscFunctionReturn(PETSC_SUCCESS);
210 }
211 
212 /*@C
213   MatDenseHIPGetArrayRead - Provides read-only access to the HIP buffer inside a
214   `MATDENSEHIP` matrix. The array must be restored with `MatDenseHIPRestoreArrayRead()` when
215   no longer needed.
216 
217   Not Collective
218 
219   Input Parameter:
220 . A - the matrix
221 
222   Output Parameter:
223 . a - the GPU array in column major order
224 
225   Level: developer
226 
227   Note:
228   Data may be copied to the GPU due to operations done on the CPU. If you need write only
229   access, use `MatDenseHIPGetArrayWrite()`.
230 
231 .seealso: `MATDENSEHIP`, `MatDenseHIPGetArray()`, `MatDenseHIPRestoreArray()`,
232 `MatDenseHIPRestoreArrayWrite()`, `MatDenseHIPGetArrayWrite()`,
233 `MatDenseHIPRestoreArrayRead()`
234 @*/
MatDenseHIPGetArrayRead(Mat A,const PetscScalar ** a)235 PetscErrorCode MatDenseHIPGetArrayRead(Mat A, const PetscScalar **a)
236 {
237   PetscFunctionBegin;
238   PetscCall(MatDenseCUPMGetArrayRead<DeviceType::HIP>(A, a));
239   PetscFunctionReturn(PETSC_SUCCESS);
240 }
241 
242 /*@C
243   MatDenseHIPRestoreArrayRead - Restore read-only access to the HIP buffer inside a
244   `MATDENSEHIP` matrix previously obtained with a call to `MatDenseHIPGetArrayRead()`.
245 
246   Not Collective
247 
248   Input Parameters:
249 + A - the matrix
250 - a - the GPU array in column major order
251 
252   Level: developer
253 
254   Note:
255   Data can be copied to the GPU due to operations done on the CPU. If you need write only
256   access, use `MatDenseHIPGetArrayWrite()`.
257 
258 .seealso: `MATDENSEHIP`, `MatDenseHIPGetArray()`, `MatDenseHIPRestoreArray()`,
259 `MatDenseHIPRestoreArrayWrite()`, `MatDenseHIPGetArrayWrite()`, `MatDenseHIPGetArrayRead()`
260 @*/
MatDenseHIPRestoreArrayRead(Mat A,const PetscScalar ** a)261 PetscErrorCode MatDenseHIPRestoreArrayRead(Mat A, const PetscScalar **a)
262 {
263   PetscFunctionBegin;
264   PetscCall(MatDenseCUPMRestoreArrayRead<DeviceType::HIP>(A, a));
265   PetscFunctionReturn(PETSC_SUCCESS);
266 }
267 
268 /*@C
269   MatDenseHIPGetArray - Provides access to the HIP buffer inside a `MATDENSEHIP` matrix. The
270   array must be restored with `MatDenseHIPRestoreArray()` when no longer needed.
271 
272   Not Collective
273 
274   Input Parameter:
275 . A - the matrix
276 
277   Output Parameter:
278 . a - the GPU array in column major order
279 
280   Level: developer
281 
282   Note:
283   Data can be copied to the GPU due to operations done on the CPU. If you need write only
284   access, use `MatDenseHIPGetArrayWrite()`. For read-only access, use
285   `MatDenseHIPGetArrayRead()`.
286 
287 .seealso: `MATDENSEHIP`, `MatDenseHIPGetArrayRead()`, `MatDenseHIPRestoreArray()`,
288 `MatDenseHIPRestoreArrayWrite()`, `MatDenseHIPGetArrayWrite()`,
289 `MatDenseHIPRestoreArrayRead()`
290 @*/
MatDenseHIPGetArray(Mat A,PetscScalar ** a)291 PetscErrorCode MatDenseHIPGetArray(Mat A, PetscScalar **a)
292 {
293   PetscFunctionBegin;
294   PetscCall(MatDenseCUPMGetArray<DeviceType::HIP>(A, a));
295   PetscFunctionReturn(PETSC_SUCCESS);
296 }
297 
298 /*@C
299   MatDenseHIPRestoreArray - Restore access to the HIP buffer inside a `MATDENSEHIP` matrix
300   previously obtained with `MatDenseHIPGetArray()`.
301 
302   Not Collective
303 
304   Level: developer
305 
306   Input Parameters:
307 + A - the matrix
308 - a - the GPU array in column major order
309 
310 .seealso: `MATDENSEHIP`, `MatDenseHIPGetArray()`, `MatDenseHIPRestoreArrayWrite()`,
311 `MatDenseHIPGetArrayWrite()`, `MatDenseHIPRestoreArrayRead()`, `MatDenseHIPGetArrayRead()`
312 @*/
MatDenseHIPRestoreArray(Mat A,PetscScalar ** a)313 PetscErrorCode MatDenseHIPRestoreArray(Mat A, PetscScalar **a)
314 {
315   PetscFunctionBegin;
316   PetscCall(MatDenseCUPMRestoreArray<DeviceType::HIP>(A, a));
317   PetscFunctionReturn(PETSC_SUCCESS);
318 }
319 
320 /*@C
321   MatDenseHIPSetPreallocation - Set the device array used for storing the matrix elements of a
322   `MATDENSEHIP` matrix
323 
324   Collective
325 
326   Input Parameters:
327 + A            - the matrix
328 - device_array - the array (or `NULL`)
329 
330   Level: intermediate
331 
332 .seealso: [](ch_matrices), `Mat`, `MATDENSEHIP`, `MatCreate()`, `MatCreateDenseHIP()`,
333 `MatSetValues()`, `MatDenseSetLDA()`
334 @*/
MatDenseHIPSetPreallocation(Mat A,PetscScalar * device_array)335 PetscErrorCode MatDenseHIPSetPreallocation(Mat A, PetscScalar *device_array)
336 {
337   PetscFunctionBegin;
338   PetscCall(MatDenseCUPMSetPreallocation<DeviceType::HIP>(A, device_array));
339   PetscFunctionReturn(PETSC_SUCCESS);
340 }
341