petsc-3.15.0 2021-03-30
MatDenseCUDAPlaceArray
Allows one to replace the GPU array in a dense matrix with an array provided by the user. This is useful to avoid copying an array into a matrix
Synopsis
#include "petscmat.h"
PetscErrorCode MatDenseCUDAPlaceArray(Mat mat,const PetscScalar *array)
Not Collective
Input Parameters
| mat | - the matrix
|
| array | - the array in column major order
|
Notes
You can return to the original array with a call to MatDenseCUDAResetArray(). The user is responsible for freeing this array; it will not be
freed when the matrix is destroyed. The array must have been allocated with cudaMalloc().
See Also
MatDenseCUDAGetArray(), MatDenseCUDAResetArray()
Level
developer
Location
src/mat/impls/dense/mpi/mpidense.c
Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages