1*5d83a8b1SBarry Smith #include <petscmat.h> /*I "petscmat.h" I*/ 2ad6ad4e1SHansol Suh 3*5d83a8b1SBarry Smith /*@ 4d16ceb75SStefano Zampini MatCreateDenseFromVecType - Create a matrix that matches the type of a Vec. 5ad6ad4e1SHansol Suh 6ad6ad4e1SHansol Suh Collective 7ad6ad4e1SHansol Suh 8ad6ad4e1SHansol Suh Input Parameters: 9d16ceb75SStefano Zampini + comm - the communicator 10d16ceb75SStefano 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) 15d16ceb75SStefano Zampini . lda - optional leading dimension. Pass any non-positive number to use the default. 16d16ceb75SStefano 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: 19d16ceb75SStefano Zampini . A - the dense matrix 20ad6ad4e1SHansol Suh 21ad6ad4e1SHansol Suh Level: advanced 22ad6ad4e1SHansol Suh 232913281cSPierre Jolivet .seealso: [](ch_matrices), `Mat`, `MatCreateDense()`, `MatCreateDenseCUDA()`, `MatCreateDenseHIP()`, `PetscMemType` 24ad6ad4e1SHansol Suh @*/ 25d16ceb75SStefano Zampini PetscErrorCode MatCreateDenseFromVecType(MPI_Comm comm, VecType vtype, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt lda, PetscScalar *data, Mat *A) 26ad6ad4e1SHansol Suh { 27d16ceb75SStefano Zampini VecType root_type = VECSTANDARD; 28d16ceb75SStefano Zampini PetscBool isstd, iscuda, iship, iskokkos; 29ad6ad4e1SHansol Suh 30ad6ad4e1SHansol Suh PetscFunctionBegin; 31d16ceb75SStefano Zampini PetscCall(PetscStrcmpAny(vtype, &isstd, VECSTANDARD, VECMPI, VECSEQ, "")); 32d16ceb75SStefano Zampini PetscCall(PetscStrcmpAny(vtype, &iscuda, VECCUDA, VECMPICUDA, VECSEQCUDA, "")); 33d16ceb75SStefano Zampini PetscCall(PetscStrcmpAny(vtype, &iship, VECHIP, VECMPIHIP, VECSEQHIP, "")); 34d16ceb75SStefano Zampini PetscCall(PetscStrcmpAny(vtype, &iskokkos, VECKOKKOS, VECMPIKOKKOS, VECSEQKOKKOS, "")); 3500045ab3SPierre Jolivet PetscCheck(isstd || iscuda || iship || iskokkos, comm, PETSC_ERR_SUP, "Not for type %s", vtype); 36d16ceb75SStefano Zampini if (iscuda) root_type = VECCUDA; 37d16ceb75SStefano Zampini else if (iship) root_type = VECHIP; 38d16ceb75SStefano Zampini else if (iskokkos) { 39d16ceb75SStefano Zampini /* We support only one type of kokkos device */ 40d16ceb75SStefano Zampini PetscCheck(!PetscDefined(HAVE_MACRO_KOKKOS_ENABLE_SYCL), comm, PETSC_ERR_SUP, "Not for sycl backend"); 41d16ceb75SStefano Zampini if (PetscDefined(HAVE_MACRO_KOKKOS_ENABLE_CUDA)) iscuda = PETSC_TRUE; 42d16ceb75SStefano Zampini else if (PetscDefined(HAVE_MACRO_KOKKOS_ENABLE_HIP)) iship = PETSC_TRUE; 43d16ceb75SStefano Zampini else isstd = PETSC_TRUE; 44d16ceb75SStefano Zampini root_type = VECKOKKOS; 45ad6ad4e1SHansol Suh } 46d16ceb75SStefano Zampini PetscCall(MatCreate(comm, A)); 47d16ceb75SStefano Zampini PetscCall(MatSetSizes(*A, m, n, M, N)); 48ad6ad4e1SHansol Suh if (isstd) { 49d16ceb75SStefano Zampini PetscCall(MatSetType(*A, MATDENSE)); 50d16ceb75SStefano Zampini if (lda > 0) PetscCall(MatDenseSetLDA(*A, lda)); 51d16ceb75SStefano Zampini PetscCall(MatSeqDenseSetPreallocation(*A, data)); 52d16ceb75SStefano Zampini PetscCall(MatMPIDenseSetPreallocation(*A, data)); 53d16ceb75SStefano Zampini } else if (iscuda) { 54d16ceb75SStefano Zampini PetscCheck(PetscDefined(HAVE_CUDA), comm, PETSC_ERR_SUP, "PETSc not compiled with CUDA support"); 55ad6ad4e1SHansol Suh #if defined(PETSC_HAVE_CUDA) 56d16ceb75SStefano Zampini PetscCall(MatSetType(*A, MATDENSECUDA)); 57d16ceb75SStefano Zampini if (lda > 0) PetscCall(MatDenseSetLDA(*A, lda)); 58d16ceb75SStefano Zampini PetscCall(MatDenseCUDASetPreallocation(*A, data)); 59ad6ad4e1SHansol Suh #endif 60d16ceb75SStefano Zampini } else if (iship) { 61d16ceb75SStefano Zampini PetscCheck(PetscDefined(HAVE_HIP), comm, PETSC_ERR_SUP, "PETSc not compiled with HIP support"); 62ad6ad4e1SHansol Suh #if defined(PETSC_HAVE_HIP) 63d16ceb75SStefano Zampini PetscCall(MatSetType(*A, MATDENSEHIP)); 64d16ceb75SStefano Zampini if (lda > 0) PetscCall(MatDenseSetLDA(*A, lda)); 65d16ceb75SStefano Zampini PetscCall(MatDenseHIPSetPreallocation(*A, data)); 66ad6ad4e1SHansol Suh #endif 67d16ceb75SStefano Zampini } 68d16ceb75SStefano Zampini PetscCall(MatSetVecType(*A, root_type)); 69ad6ad4e1SHansol Suh PetscFunctionReturn(PETSC_SUCCESS); 70ad6ad4e1SHansol Suh } 71