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