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