xref: /petsc/src/mat/impls/dense/mpi/cupm/cuda/matmpidensecuda.cu (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
14742e46bSJacob Faibussowitsch #include "../matmpidensecupm.hpp"
24742e46bSJacob Faibussowitsch 
34742e46bSJacob Faibussowitsch using namespace Petsc::mat::cupm;
44742e46bSJacob Faibussowitsch using Petsc::device::cupm::DeviceType;
54742e46bSJacob Faibussowitsch 
64742e46bSJacob Faibussowitsch static constexpr impl::MatDense_MPI_CUPM<DeviceType::CUDA> mat_cupm{};
74742e46bSJacob Faibussowitsch 
84742e46bSJacob Faibussowitsch /*MC
94742e46bSJacob Faibussowitsch   MATDENSECUDA - "densecuda" - A matrix type to be used for dense matrices on GPUs.
104742e46bSJacob Faibussowitsch 
114742e46bSJacob Faibussowitsch   This matrix type is identical to `MATSEQDENSECUDA` when constructed with a single process
124742e46bSJacob Faibussowitsch   communicator, and `MATMPIDENSECUDA` otherwise.
134742e46bSJacob Faibussowitsch 
144742e46bSJacob Faibussowitsch   Options Database Key:
154742e46bSJacob Faibussowitsch . -mat_type densecuda - sets the matrix type to `MATDENSECUDA` during a call to
164742e46bSJacob Faibussowitsch                         `MatSetFromOptions()`
174742e46bSJacob Faibussowitsch 
184742e46bSJacob Faibussowitsch   Level: beginner
194742e46bSJacob Faibussowitsch 
20*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATSEQDENSECUDA`, `MATMPIDENSECUDA`, `MATSEQDENSEHIP`,
214742e46bSJacob Faibussowitsch `MATMPIDENSEHIP`, `MATDENSE`
224742e46bSJacob Faibussowitsch M*/
234742e46bSJacob Faibussowitsch 
244742e46bSJacob Faibussowitsch /*MC
254742e46bSJacob Faibussowitsch   MATMPIDENSECUDA - "mpidensecuda" - A matrix type to be used for distributed dense matrices on
264742e46bSJacob Faibussowitsch   GPUs.
274742e46bSJacob Faibussowitsch 
284742e46bSJacob Faibussowitsch   Options Database Key:
294742e46bSJacob Faibussowitsch . -mat_type mpidensecuda - sets the matrix type to `MATMPIDENSECUDA` during a call to
304742e46bSJacob Faibussowitsch                            `MatSetFromOptions()`
314742e46bSJacob Faibussowitsch 
324742e46bSJacob Faibussowitsch   Level: beginner
334742e46bSJacob Faibussowitsch 
34*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATDENSECUDA`, `MATMPIDENSE`, `MATSEQDENSE`,
354742e46bSJacob Faibussowitsch `MATSEQDENSECUDA`, `MATSEQDENSEHIP`
364742e46bSJacob Faibussowitsch M*/
374742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatCreate_MPIDenseCUDA(Mat A)
384742e46bSJacob Faibussowitsch {
394742e46bSJacob Faibussowitsch   PetscFunctionBegin;
404742e46bSJacob Faibussowitsch   PetscCall(mat_cupm.Create(A));
414742e46bSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
424742e46bSJacob Faibussowitsch }
434742e46bSJacob Faibussowitsch 
444742e46bSJacob Faibussowitsch PetscErrorCode MatConvert_MPIDense_MPIDenseCUDA(Mat A, MatType type, MatReuse reuse, Mat *ret)
454742e46bSJacob Faibussowitsch {
464742e46bSJacob Faibussowitsch   PetscFunctionBegin;
474742e46bSJacob Faibussowitsch   PetscCall(mat_cupm.Convert_MPIDense_MPIDenseCUPM(A, type, reuse, ret));
484742e46bSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
494742e46bSJacob Faibussowitsch }
504742e46bSJacob Faibussowitsch 
514742e46bSJacob Faibussowitsch /*@C
524742e46bSJacob Faibussowitsch   MatCreateDenseCUDA - Creates a matrix in `MATDENSECUDA` format using CUDA.
534742e46bSJacob Faibussowitsch 
544742e46bSJacob Faibussowitsch   Collective
554742e46bSJacob Faibussowitsch 
564742e46bSJacob Faibussowitsch   Input Parameters:
574742e46bSJacob Faibussowitsch + comm - MPI communicator
584742e46bSJacob Faibussowitsch . m    - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
594742e46bSJacob Faibussowitsch . n    - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given)
604742e46bSJacob Faibussowitsch . M    - number of global rows (or `PETSC_DECIDE` to have calculated if `m` is given)
614742e46bSJacob Faibussowitsch . N    - number of global columns (or `PETSC_DECIDE` to have calculated if `n` is given)
622fe279fdSBarry Smith - data - optional location of GPU matrix data. Pass `NULL` to have PETSc to control matrix memory allocation.
634742e46bSJacob Faibussowitsch 
644742e46bSJacob Faibussowitsch   Output Parameter:
654742e46bSJacob Faibussowitsch . A - the matrix
664742e46bSJacob Faibussowitsch 
674742e46bSJacob Faibussowitsch   Level: intermediate
684742e46bSJacob Faibussowitsch 
694742e46bSJacob Faibussowitsch .seealso: `MATDENSECUDA`, `MatCreate()`, `MatCreateDense()`
704742e46bSJacob Faibussowitsch @*/
714742e46bSJacob Faibussowitsch PetscErrorCode MatCreateDenseCUDA(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscScalar *data, Mat *A)
724742e46bSJacob Faibussowitsch {
734742e46bSJacob Faibussowitsch   PetscFunctionBegin;
744742e46bSJacob Faibussowitsch   PetscCall(MatCreateDenseCUPM<DeviceType::CUDA>(comm, m, n, M, N, data, A));
754742e46bSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
764742e46bSJacob Faibussowitsch }
774742e46bSJacob Faibussowitsch 
784742e46bSJacob Faibussowitsch /*@C
794742e46bSJacob Faibussowitsch   MatDenseCUDAPlaceArray - Allows one to replace the GPU array in a `MATDENSECUDA` matrix with an
804742e46bSJacob Faibussowitsch   array provided by the user. This is useful to avoid copying an array into a matrix.
814742e46bSJacob Faibussowitsch 
824742e46bSJacob Faibussowitsch   Not Collective
834742e46bSJacob Faibussowitsch 
844742e46bSJacob Faibussowitsch   Input Parameters:
854742e46bSJacob Faibussowitsch + mat   - the matrix
864742e46bSJacob Faibussowitsch - array - the array in column major order
874742e46bSJacob Faibussowitsch 
884742e46bSJacob Faibussowitsch   Level: developer
894742e46bSJacob Faibussowitsch 
904742e46bSJacob Faibussowitsch   Note:
914742e46bSJacob Faibussowitsch   You can return to the original array with a call to `MatDenseCUDAResetArray()`. The user is
924742e46bSJacob Faibussowitsch   responsible for freeing this array; it will not be freed when the matrix is destroyed. The
934742e46bSJacob Faibussowitsch   array must have been allocated with `cudaMalloc()`.
944742e46bSJacob Faibussowitsch 
954742e46bSJacob Faibussowitsch .seealso: `MATDENSECUDA`, `MatDenseCUDAGetArray()`, `MatDenseCUDAResetArray()`,
964742e46bSJacob Faibussowitsch           `MatDenseCUDAReplaceArray()`
974742e46bSJacob Faibussowitsch @*/
984742e46bSJacob Faibussowitsch PetscErrorCode MatDenseCUDAPlaceArray(Mat mat, const PetscScalar *array)
994742e46bSJacob Faibussowitsch {
1004742e46bSJacob Faibussowitsch   PetscFunctionBegin;
1014742e46bSJacob Faibussowitsch   PetscCall(MatDenseCUPMPlaceArray<DeviceType::CUDA>(mat, array));
1024742e46bSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1034742e46bSJacob Faibussowitsch }
1044742e46bSJacob Faibussowitsch 
1054742e46bSJacob Faibussowitsch /*@C
1064742e46bSJacob Faibussowitsch   MatDenseCUDAResetArray - Resets the matrix array to that it previously had before the call to
1074742e46bSJacob Faibussowitsch   `MatDenseCUDAPlaceArray()`
1084742e46bSJacob Faibussowitsch 
1094742e46bSJacob Faibussowitsch   Not Collective
1104742e46bSJacob Faibussowitsch 
1112fe279fdSBarry Smith   Input Parameter:
1124742e46bSJacob Faibussowitsch . mat - the matrix
1134742e46bSJacob Faibussowitsch 
1144742e46bSJacob Faibussowitsch   Level: developer
1154742e46bSJacob Faibussowitsch 
1164742e46bSJacob Faibussowitsch   Note:
1174742e46bSJacob Faibussowitsch   You can only call this after a call to `MatDenseCUDAPlaceArray()`
1184742e46bSJacob Faibussowitsch 
1194742e46bSJacob Faibussowitsch .seealso: `MATDENSECUDA`, `MatDenseCUDAGetArray()`, `MatDenseCUDAPlaceArray()`
1204742e46bSJacob Faibussowitsch @*/
1214742e46bSJacob Faibussowitsch PetscErrorCode MatDenseCUDAResetArray(Mat mat)
1224742e46bSJacob Faibussowitsch {
1234742e46bSJacob Faibussowitsch   PetscFunctionBegin;
1244742e46bSJacob Faibussowitsch   PetscCall(MatDenseCUPMResetArray<DeviceType::CUDA>(mat));
1254742e46bSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1264742e46bSJacob Faibussowitsch }
1274742e46bSJacob Faibussowitsch 
1284742e46bSJacob Faibussowitsch /*@C
1294742e46bSJacob Faibussowitsch   MatDenseCUDAReplaceArray - Allows one to replace the GPU array in a `MATDENSECUDA` matrix
1304742e46bSJacob Faibussowitsch   with an array provided by the user. This is useful to avoid copying an array into a matrix.
1314742e46bSJacob Faibussowitsch 
1324742e46bSJacob Faibussowitsch   Not Collective
1334742e46bSJacob Faibussowitsch 
1344742e46bSJacob Faibussowitsch   Input Parameters:
1354742e46bSJacob Faibussowitsch + mat   - the matrix
1364742e46bSJacob Faibussowitsch - array - the array in column major order
1374742e46bSJacob Faibussowitsch 
1384742e46bSJacob Faibussowitsch   Level: developer
1394742e46bSJacob Faibussowitsch 
1404742e46bSJacob Faibussowitsch   Note:
1414742e46bSJacob Faibussowitsch   This permanently replaces the GPU array and frees the memory associated with the old GPU
1424742e46bSJacob Faibussowitsch   array. The memory passed in CANNOT be freed by the user. It will be freed when the matrix is
1434742e46bSJacob Faibussowitsch   destroyed. The array should respect the matrix leading dimension.
1444742e46bSJacob Faibussowitsch 
1454742e46bSJacob Faibussowitsch .seealso: `MatDenseCUDAGetArray()`, `MatDenseCUDAPlaceArray()`, `MatDenseCUDAResetArray()`
1464742e46bSJacob Faibussowitsch @*/
1474742e46bSJacob Faibussowitsch PetscErrorCode MatDenseCUDAReplaceArray(Mat mat, const PetscScalar *array)
1484742e46bSJacob Faibussowitsch {
1494742e46bSJacob Faibussowitsch   PetscFunctionBegin;
1504742e46bSJacob Faibussowitsch   PetscCall(MatDenseCUPMReplaceArray<DeviceType::CUDA>(mat, array));
1514742e46bSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1524742e46bSJacob Faibussowitsch }
1534742e46bSJacob Faibussowitsch 
1544742e46bSJacob Faibussowitsch /*@C
1554742e46bSJacob Faibussowitsch   MatDenseCUDAGetArrayWrite - Provides write access to the CUDA buffer inside a `MATDENSECUDA`
1564742e46bSJacob Faibussowitsch   matrix.
1574742e46bSJacob Faibussowitsch 
1584742e46bSJacob Faibussowitsch   Not Collective
1594742e46bSJacob Faibussowitsch 
1602fe279fdSBarry Smith   Input Parameter:
1614742e46bSJacob Faibussowitsch . A - the matrix
1624742e46bSJacob Faibussowitsch 
1632fe279fdSBarry Smith   Output Parameter:
1642fe279fdSBarry Smith . a - the GPU array in column major order
1654742e46bSJacob Faibussowitsch 
1664742e46bSJacob Faibussowitsch   Level: developer
1674742e46bSJacob Faibussowitsch 
1684742e46bSJacob Faibussowitsch   Notes:
1694742e46bSJacob Faibussowitsch   The data on the GPU may not be updated due to operations done on the CPU. If you need updated
1704742e46bSJacob Faibussowitsch   data, use `MatDenseCUDAGetArray()`.
1714742e46bSJacob Faibussowitsch 
1724742e46bSJacob Faibussowitsch   The array must be restored with `MatDenseCUDARestoreArrayWrite()` when no longer needed.
1734742e46bSJacob Faibussowitsch 
1744742e46bSJacob Faibussowitsch .seealso: `MATDENSECUDA`, `MatDenseCUDAGetArray()`, `MatDenseCUDARestoreArray()`,
1754742e46bSJacob Faibussowitsch           `MatDenseCUDARestoreArrayWrite()`, `MatDenseCUDAGetArrayRead()`,
1764742e46bSJacob Faibussowitsch           `MatDenseCUDARestoreArrayRead()`
1774742e46bSJacob Faibussowitsch @*/
1784742e46bSJacob Faibussowitsch PetscErrorCode MatDenseCUDAGetArrayWrite(Mat A, PetscScalar **a)
1794742e46bSJacob Faibussowitsch {
1804742e46bSJacob Faibussowitsch   PetscFunctionBegin;
1814742e46bSJacob Faibussowitsch   PetscCall(MatDenseCUPMGetArrayWrite<DeviceType::CUDA>(A, a));
1824742e46bSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1834742e46bSJacob Faibussowitsch }
1844742e46bSJacob Faibussowitsch 
1854742e46bSJacob Faibussowitsch /*@C
1864742e46bSJacob Faibussowitsch   MatDenseCUDARestoreArrayWrite - Restore write access to the CUDA buffer inside a
1874742e46bSJacob Faibussowitsch   `MATDENSECUDA` matrix previously obtained with `MatDenseCUDAGetArrayWrite()`.
1884742e46bSJacob Faibussowitsch 
1894742e46bSJacob Faibussowitsch   Not Collective
1904742e46bSJacob Faibussowitsch 
1914742e46bSJacob Faibussowitsch   Input Parameters:
1924742e46bSJacob Faibussowitsch + A     - the matrix
1932fe279fdSBarry Smith - a - the GPU array in column major order
1944742e46bSJacob Faibussowitsch 
1954742e46bSJacob Faibussowitsch   Level: developer
1964742e46bSJacob Faibussowitsch 
1974742e46bSJacob Faibussowitsch .seealso: `MATDENSECUDA`, `MatDenseCUDAGetArray()`, `MatDenseCUDARestoreArray()`,
1984742e46bSJacob Faibussowitsch `MatDenseCUDAGetArrayWrite()`, `MatDenseCUDARestoreArrayRead()`, `MatDenseCUDAGetArrayRead()`
1994742e46bSJacob Faibussowitsch @*/
2004742e46bSJacob Faibussowitsch PetscErrorCode MatDenseCUDARestoreArrayWrite(Mat A, PetscScalar **a)
2014742e46bSJacob Faibussowitsch {
2024742e46bSJacob Faibussowitsch   PetscFunctionBegin;
2034742e46bSJacob Faibussowitsch   PetscCall(MatDenseCUPMRestoreArrayWrite<DeviceType::CUDA>(A, a));
2044742e46bSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2054742e46bSJacob Faibussowitsch }
2064742e46bSJacob Faibussowitsch 
2074742e46bSJacob Faibussowitsch /*@C
2084742e46bSJacob Faibussowitsch   MatDenseCUDAGetArrayRead - Provides read-only access to the CUDA buffer inside a
2094742e46bSJacob Faibussowitsch   `MATDENSECUDA` matrix. The array must be restored with `MatDenseCUDARestoreArrayRead()` when
2104742e46bSJacob Faibussowitsch   no longer needed.
2114742e46bSJacob Faibussowitsch 
2124742e46bSJacob Faibussowitsch   Not Collective
2134742e46bSJacob Faibussowitsch 
2142fe279fdSBarry Smith   Input Parameter:
2154742e46bSJacob Faibussowitsch . A - the matrix
2164742e46bSJacob Faibussowitsch 
2172fe279fdSBarry Smith   Output Parameter:
2182fe279fdSBarry Smith . a - the GPU array in column major order
2194742e46bSJacob Faibussowitsch 
2204742e46bSJacob Faibussowitsch   Level: developer
2214742e46bSJacob Faibussowitsch 
2224742e46bSJacob Faibussowitsch   Note:
2234742e46bSJacob Faibussowitsch   Data may be copied to the GPU due to operations done on the CPU. If you need write only
2244742e46bSJacob Faibussowitsch   access, use `MatDenseCUDAGetArrayWrite()`.
2254742e46bSJacob Faibussowitsch 
2264742e46bSJacob Faibussowitsch .seealso: `MATDENSECUDA`, `MatDenseCUDAGetArray()`, `MatDenseCUDARestoreArray()`,
2274742e46bSJacob Faibussowitsch           `MatDenseCUDARestoreArrayWrite()`, `MatDenseCUDAGetArrayWrite()`,
2284742e46bSJacob Faibussowitsch           `MatDenseCUDARestoreArrayRead()`
2294742e46bSJacob Faibussowitsch @*/
2304742e46bSJacob Faibussowitsch PetscErrorCode MatDenseCUDAGetArrayRead(Mat A, const PetscScalar **a)
2314742e46bSJacob Faibussowitsch {
2324742e46bSJacob Faibussowitsch   PetscFunctionBegin;
2334742e46bSJacob Faibussowitsch   PetscCall(MatDenseCUPMGetArrayRead<DeviceType::CUDA>(A, a));
2344742e46bSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2354742e46bSJacob Faibussowitsch }
2364742e46bSJacob Faibussowitsch 
2374742e46bSJacob Faibussowitsch /*@C
2384742e46bSJacob Faibussowitsch   MatDenseCUDARestoreArrayRead - Restore read-only access to the CUDA buffer inside a
2394742e46bSJacob Faibussowitsch   `MATDENSECUDA` matrix previously obtained with a call to `MatDenseCUDAGetArrayRead()`.
2404742e46bSJacob Faibussowitsch 
2414742e46bSJacob Faibussowitsch   Not Collective
2424742e46bSJacob Faibussowitsch 
2434742e46bSJacob Faibussowitsch   Input Parameters:
2444742e46bSJacob Faibussowitsch + A     - the matrix
2452fe279fdSBarry Smith - a - the GPU array in column major order
2464742e46bSJacob Faibussowitsch 
2474742e46bSJacob Faibussowitsch   Level: developer
2484742e46bSJacob Faibussowitsch 
2494742e46bSJacob Faibussowitsch   Note:
2504742e46bSJacob Faibussowitsch   Data can be copied to the GPU due to operations done on the CPU. If you need write only
2514742e46bSJacob Faibussowitsch   access, use `MatDenseCUDAGetArrayWrite()`.
2524742e46bSJacob Faibussowitsch 
2534742e46bSJacob Faibussowitsch .seealso: `MATDENSECUDA`, `MatDenseCUDAGetArray()`, `MatDenseCUDARestoreArray()`,
2544742e46bSJacob Faibussowitsch           `MatDenseCUDARestoreArrayWrite()`, `MatDenseCUDAGetArrayWrite()`, `MatDenseCUDAGetArrayRead()`
2554742e46bSJacob Faibussowitsch @*/
2564742e46bSJacob Faibussowitsch PetscErrorCode MatDenseCUDARestoreArrayRead(Mat A, const PetscScalar **a)
2574742e46bSJacob Faibussowitsch {
2584742e46bSJacob Faibussowitsch   PetscFunctionBegin;
2594742e46bSJacob Faibussowitsch   PetscCall(MatDenseCUPMRestoreArrayRead<DeviceType::CUDA>(A, a));
2604742e46bSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2614742e46bSJacob Faibussowitsch }
2624742e46bSJacob Faibussowitsch 
2634742e46bSJacob Faibussowitsch /*@C
2644742e46bSJacob Faibussowitsch   MatDenseCUDAGetArray - Provides access to the CUDA buffer inside a `MATDENSECUDA` matrix. The
2654742e46bSJacob Faibussowitsch   array must be restored with `MatDenseCUDARestoreArray()` when no longer needed.
2664742e46bSJacob Faibussowitsch 
2674742e46bSJacob Faibussowitsch   Not Collective
2684742e46bSJacob Faibussowitsch 
2692fe279fdSBarry Smith   Input Parameter:
2704742e46bSJacob Faibussowitsch . A - the matrix
2714742e46bSJacob Faibussowitsch 
2722fe279fdSBarry Smith   Output Parameter:
2732fe279fdSBarry Smith . a - the GPU array in column major order
2744742e46bSJacob Faibussowitsch 
2754742e46bSJacob Faibussowitsch   Level: developer
2764742e46bSJacob Faibussowitsch 
2774742e46bSJacob Faibussowitsch   Note:
2784742e46bSJacob Faibussowitsch   Data can be copied to the GPU due to operations done on the CPU. If you need write only
2794742e46bSJacob Faibussowitsch   access, use `MatDenseCUDAGetArrayWrite()`. For read-only access, use
2804742e46bSJacob Faibussowitsch   `MatDenseCUDAGetArrayRead()`.
2814742e46bSJacob Faibussowitsch 
2824742e46bSJacob Faibussowitsch .seealso: `MATDENSECUDA`, `MatDenseCUDAGetArrayRead()`, `MatDenseCUDARestoreArray()`,
2834742e46bSJacob Faibussowitsch           `MatDenseCUDARestoreArrayWrite()`, `MatDenseCUDAGetArrayWrite()`,
2844742e46bSJacob Faibussowitsch           `MatDenseCUDARestoreArrayRead()`
2854742e46bSJacob Faibussowitsch @*/
2864742e46bSJacob Faibussowitsch PetscErrorCode MatDenseCUDAGetArray(Mat A, PetscScalar **a)
2874742e46bSJacob Faibussowitsch {
2884742e46bSJacob Faibussowitsch   PetscFunctionBegin;
2894742e46bSJacob Faibussowitsch   PetscCall(MatDenseCUPMGetArray<DeviceType::CUDA>(A, a));
2904742e46bSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2914742e46bSJacob Faibussowitsch }
2924742e46bSJacob Faibussowitsch 
2934742e46bSJacob Faibussowitsch /*@C
2944742e46bSJacob Faibussowitsch   MatDenseCUDARestoreArray - Restore access to the CUDA buffer inside a `MATDENSECUDA` matrix
2954742e46bSJacob Faibussowitsch   previously obtained with `MatDenseCUDAGetArray()`.
2964742e46bSJacob Faibussowitsch 
2974742e46bSJacob Faibussowitsch   Not Collective
2984742e46bSJacob Faibussowitsch 
2994742e46bSJacob Faibussowitsch   Level: developer
3004742e46bSJacob Faibussowitsch 
3014742e46bSJacob Faibussowitsch   Input Parameters:
3024742e46bSJacob Faibussowitsch + A - the matrix
3032fe279fdSBarry Smith - a - the GPU array in column major order
3044742e46bSJacob Faibussowitsch 
3054742e46bSJacob Faibussowitsch .seealso: `MATDENSECUDA`, `MatDenseCUDAGetArray()`, `MatDenseCUDARestoreArrayWrite()`,
3064742e46bSJacob Faibussowitsch           `MatDenseCUDAGetArrayWrite()`, `MatDenseCUDARestoreArrayRead()`, `MatDenseCUDAGetArrayRead()`
3074742e46bSJacob Faibussowitsch @*/
3084742e46bSJacob Faibussowitsch PetscErrorCode MatDenseCUDARestoreArray(Mat A, PetscScalar **a)
3094742e46bSJacob Faibussowitsch {
3104742e46bSJacob Faibussowitsch   PetscFunctionBegin;
3114742e46bSJacob Faibussowitsch   PetscCall(MatDenseCUPMRestoreArray<DeviceType::CUDA>(A, a));
3124742e46bSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3134742e46bSJacob Faibussowitsch }
314