xref: /petsc/src/mat/utils/veccreatematdense.c (revision ed00a6c129cbc4f1b4754b50ab0c4e7536acfadd)
1 #include <petscmat.h>              /*I "petscmat.h" I*/
2 #include <petsc/private/vecimpl.h> /*I "petscvec.h" I*/
3 
4 /*@
5   VecCreateMatDense - Create a matrix that matches the type of a Vec.
6 
7   Collective
8 
9   Input Parameters:
10 + X    - the vector
11 . m    - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
12 . n    - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given)
13 . M    - number of global rows (or `PETSC_DECIDE` to have calculated if `m` is given)
14 . N    - number of global columns (or `PETSC_DECIDE` to have calculated if `n` is given)
15 - data - optional location of matrix data, which should have the same memory type as the vector. Pass `NULL` to have PETSc to control matrix.
16          memory allocation.
17 
18   Output Parameter:
19 . A - the matrix.  `A` will have the same communicator as `X` and the same `PetscMemType`.
20 
21   Level: advanced
22 
23 .seealso: [](chapter_matrices), `Mat`, `MatCreateDense()', `MatCreateDenseCUDA()`, `MatCreateDenseHIP()`, `PetscMemType`
24 @*/
25 PetscErrorCode VecCreateMatDense(Vec X, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscScalar *data, Mat *A)
26 {
27   VecType   root_type;
28   PetscBool isstd, iscuda, iship;
29   MPI_Comm  comm;
30 
31   PetscFunctionBegin;
32   PetscCall(VecGetRootType_Private(X, &root_type));
33   PetscCall(PetscObjectGetComm((PetscObject)X, &comm));
34   PetscCall(PetscStrcmp(root_type, VECSTANDARD, &isstd));
35   PetscCall(PetscStrcmp(root_type, VECCUDA, &iscuda));
36   PetscCall(PetscStrcmp(root_type, VECHIP, &iship));
37 
38   /* For performance-portable types (Kokkos, SYCL, ...) that dispatch to */
39   if (!(isstd || iscuda || iship)) {
40     const PetscScalar *array;
41     PetscMemType       memtype;
42 
43     PetscCall(VecGetArrayReadAndMemType(X, &array, &memtype));
44     PetscCall(VecRestoreArrayReadAndMemType(X, &array));
45     switch (memtype) {
46     case PETSC_MEMTYPE_HOST:
47       isstd = PETSC_TRUE;
48       break;
49     case PETSC_MEMTYPE_CUDA:
50     case PETSC_MEMTYPE_NVSHMEM:
51       iscuda = PETSC_TRUE;
52       break;
53     case PETSC_MEMTYPE_HIP:
54       iship = PETSC_TRUE;
55       break;
56     default:
57       SETERRQ(PetscObjectComm((PetscObject)X), PETSC_ERR_SUP, "Cannot figure out memory type of vector type %s", root_type);
58     }
59   }
60 
61   if (isstd) {
62     PetscCall(MatCreateDense(comm, m, n, M, N, data, A));
63   }
64 #if defined(PETSC_HAVE_CUDA)
65   else if (iscuda) {
66     PetscCall(MatCreateDenseCUDA(comm, m, n, M, N, data, A));
67   }
68 #endif
69 #if defined(PETSC_HAVE_HIP)
70   else if (iship) {
71     PetscCall(MatCreateDenseHIP(comm, m, n, M, N, data, A));
72   }
73 #endif
74   PetscFunctionReturn(PETSC_SUCCESS);
75 }
76