xref: /petsc/src/mat/utils/veccreatematdense.c (revision d16ceb7590f23315f034853cec344faff90b289e)
1*d16ceb75SStefano Zampini #include <petscmat.h>
2ad6ad4e1SHansol  Suh 
3*d16ceb75SStefano Zampini /*@C
4*d16ceb75SStefano Zampini   MatCreateDenseFromVecType - Create a matrix that matches the type of a Vec.
5ad6ad4e1SHansol  Suh 
6ad6ad4e1SHansol  Suh   Collective
7ad6ad4e1SHansol  Suh 
8ad6ad4e1SHansol  Suh   Input Parameters:
9*d16ceb75SStefano Zampini + comm  - the communicator
10*d16ceb75SStefano Zampini . vtype - the vector type
11ad6ad4e1SHansol  Suh . m     - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
12ad6ad4e1SHansol  Suh . n     - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given)
13ad6ad4e1SHansol  Suh . M     - number of global rows (or `PETSC_DECIDE` to have calculated if `m` is given)
14ad6ad4e1SHansol  Suh . N     - number of global columns (or `PETSC_DECIDE` to have calculated if `n` is given)
15*d16ceb75SStefano Zampini . lda   - optional leading dimension. Pass any non-positive number to use the default.
16*d16ceb75SStefano Zampini - data  - optional location of matrix data, which should have the same memory type as the vector. Pass `NULL` to have PETSc take care of matrix memory allocation.
17ad6ad4e1SHansol  Suh 
18ad6ad4e1SHansol  Suh   Output Parameter:
19*d16ceb75SStefano Zampini . A - the dense matrix
20ad6ad4e1SHansol  Suh 
21ad6ad4e1SHansol  Suh   Level: advanced
22ad6ad4e1SHansol  Suh 
23ad6ad4e1SHansol  Suh .seealso: [](chapter_matrices), `Mat`, `MatCreateDense()', `MatCreateDenseCUDA()`, `MatCreateDenseHIP()`, `PetscMemType`
24ad6ad4e1SHansol  Suh @*/
25*d16ceb75SStefano Zampini PetscErrorCode MatCreateDenseFromVecType(MPI_Comm comm, VecType vtype, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt lda, PetscScalar *data, Mat *A)
26ad6ad4e1SHansol  Suh {
27*d16ceb75SStefano Zampini   VecType   root_type = VECSTANDARD;
28*d16ceb75SStefano Zampini   PetscBool isstd, iscuda, iship, iskokkos;
29ad6ad4e1SHansol  Suh 
30ad6ad4e1SHansol  Suh   PetscFunctionBegin;
31*d16ceb75SStefano Zampini   PetscCall(PetscStrcmpAny(vtype, &isstd, VECSTANDARD, VECMPI, VECSEQ, ""));
32*d16ceb75SStefano Zampini   PetscCall(PetscStrcmpAny(vtype, &iscuda, VECCUDA, VECMPICUDA, VECSEQCUDA, ""));
33*d16ceb75SStefano Zampini   PetscCall(PetscStrcmpAny(vtype, &iship, VECHIP, VECMPIHIP, VECSEQHIP, ""));
34*d16ceb75SStefano Zampini   PetscCall(PetscStrcmpAny(vtype, &iskokkos, VECKOKKOS, VECMPIKOKKOS, VECSEQKOKKOS, ""));
35*d16ceb75SStefano Zampini   PetscCheck(isstd || iscuda || iship || iskokkos, comm, PETSC_ERR_SUP, "Not for type %s\n", vtype);
36*d16ceb75SStefano Zampini   if (iscuda) root_type = VECCUDA;
37*d16ceb75SStefano Zampini   else if (iship) root_type = VECHIP;
38*d16ceb75SStefano Zampini   else if (iskokkos) {
39*d16ceb75SStefano Zampini     /* We support only one type of kokkos device */
40*d16ceb75SStefano Zampini     PetscCheck(!PetscDefined(HAVE_MACRO_KOKKOS_ENABLE_SYCL), comm, PETSC_ERR_SUP, "Not for sycl backend");
41*d16ceb75SStefano Zampini     if (PetscDefined(HAVE_MACRO_KOKKOS_ENABLE_CUDA)) iscuda = PETSC_TRUE;
42*d16ceb75SStefano Zampini     else if (PetscDefined(HAVE_MACRO_KOKKOS_ENABLE_HIP)) iship = PETSC_TRUE;
43*d16ceb75SStefano Zampini     else isstd = PETSC_TRUE;
44*d16ceb75SStefano Zampini     root_type = VECKOKKOS;
45ad6ad4e1SHansol  Suh   }
46*d16ceb75SStefano Zampini   PetscCall(MatCreate(comm, A));
47*d16ceb75SStefano Zampini   PetscCall(MatSetSizes(*A, m, n, M, N));
48ad6ad4e1SHansol  Suh   if (isstd) {
49*d16ceb75SStefano Zampini     PetscCall(MatSetType(*A, MATDENSE));
50*d16ceb75SStefano Zampini     if (lda > 0) PetscCall(MatDenseSetLDA(*A, lda));
51*d16ceb75SStefano Zampini     PetscCall(MatSeqDenseSetPreallocation(*A, data));
52*d16ceb75SStefano Zampini     PetscCall(MatMPIDenseSetPreallocation(*A, data));
53*d16ceb75SStefano Zampini   } else if (iscuda) {
54*d16ceb75SStefano Zampini     PetscCheck(PetscDefined(HAVE_CUDA), comm, PETSC_ERR_SUP, "PETSc not compiled with CUDA support");
55ad6ad4e1SHansol  Suh #if defined(PETSC_HAVE_CUDA)
56*d16ceb75SStefano Zampini     PetscCall(MatSetType(*A, MATDENSECUDA));
57*d16ceb75SStefano Zampini     if (lda > 0) PetscCall(MatDenseSetLDA(*A, lda));
58*d16ceb75SStefano Zampini     PetscCall(MatDenseCUDASetPreallocation(*A, data));
59ad6ad4e1SHansol  Suh #endif
60*d16ceb75SStefano Zampini   } else if (iship) {
61*d16ceb75SStefano Zampini     PetscCheck(PetscDefined(HAVE_HIP), comm, PETSC_ERR_SUP, "PETSc not compiled with HIP support");
62ad6ad4e1SHansol  Suh #if defined(PETSC_HAVE_HIP)
63*d16ceb75SStefano Zampini     PetscCall(MatSetType(*A, MATDENSEHIP));
64*d16ceb75SStefano Zampini     if (lda > 0) PetscCall(MatDenseSetLDA(*A, lda));
65*d16ceb75SStefano Zampini     PetscCall(MatDenseHIPSetPreallocation(*A, data));
66ad6ad4e1SHansol  Suh #endif
67*d16ceb75SStefano Zampini   }
68*d16ceb75SStefano Zampini   PetscCall(MatSetVecType(*A, root_type));
69ad6ad4e1SHansol  Suh   PetscFunctionReturn(PETSC_SUCCESS);
70ad6ad4e1SHansol  Suh }
71