xref: /honee/include/petsc-ceed-utils.h (revision ae2b091fac884a554e48acc4b4c187524c2a2818)
1*ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2*ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
340d80af1SJames Wright #pragma once
440d80af1SJames Wright 
540d80af1SJames Wright #include <ceed.h>
640d80af1SJames Wright #include <petscdm.h>
740d80af1SJames Wright 
840d80af1SJames Wright /**
940d80af1SJames Wright   @brief Copy the reference to a `Vec`.
1040d80af1SJames Wright          Note: If `vec_copy` is non-null, it is assumed to be a valid pointer to a `Vec` and `VecDestroy()` will be called.
1140d80af1SJames Wright 
1240d80af1SJames Wright   Collective across MPI processes.
1340d80af1SJames Wright 
1440d80af1SJames Wright   @param[in]   vec       `Vec` to reference
1540d80af1SJames Wright   @param[out]  vec_copy  Copy of reference
1640d80af1SJames Wright 
1740d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
1840d80af1SJames Wright **/
1940d80af1SJames Wright static inline PetscErrorCode VecReferenceCopy(Vec vec, Vec *vec_copy) {
2040d80af1SJames Wright   PetscFunctionBeginUser;
2140d80af1SJames Wright   PetscCall(PetscObjectReference((PetscObject)vec));
2240d80af1SJames Wright   PetscCall(VecDestroy(vec_copy));
2340d80af1SJames Wright   *vec_copy = vec;
2440d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2540d80af1SJames Wright }
2640d80af1SJames Wright 
2740d80af1SJames Wright /**
2840d80af1SJames Wright   @brief Copy the reference to a `DM`.
2940d80af1SJames Wright          Note: If `dm_copy` is non-null, it is assumed to be a valid pointer to a `DM` and `DMDestroy()` will be called.
3040d80af1SJames Wright 
3140d80af1SJames Wright   Collective across MPI processes.
3240d80af1SJames Wright 
3340d80af1SJames Wright   @param[in]   dm       `DM` to reference
3440d80af1SJames Wright   @param[out]  dm_copy  Copy of reference
3540d80af1SJames Wright 
3640d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
3740d80af1SJames Wright **/
3840d80af1SJames Wright static inline PetscErrorCode DMReferenceCopy(DM dm, DM *dm_copy) {
3940d80af1SJames Wright   PetscFunctionBeginUser;
4040d80af1SJames Wright   PetscCall(PetscObjectReference((PetscObject)dm));
4140d80af1SJames Wright   PetscCall(DMDestroy(dm_copy));
4240d80af1SJames Wright   *dm_copy = dm;
4340d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4440d80af1SJames Wright }
4540d80af1SJames Wright 
4640d80af1SJames Wright /**
4740d80af1SJames Wright   @brief Translate PetscMemType to CeedMemType
4840d80af1SJames Wright 
4940d80af1SJames Wright   @param[in]  mem_type  PetscMemType
5040d80af1SJames Wright 
5140d80af1SJames Wright   @return Equivalent CeedMemType
5240d80af1SJames Wright **/
5340d80af1SJames Wright static inline CeedMemType MemTypePetscToCeed(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; }
5440d80af1SJames Wright 
5540d80af1SJames Wright /**
5640d80af1SJames Wright   @brief Translate array of `PetscInt` to `CeedInt`.
5740d80af1SJames Wright     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
5840d80af1SJames Wright     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
5940d80af1SJames Wright 
6040d80af1SJames Wright   Not collective across MPI processes.
6140d80af1SJames Wright 
6240d80af1SJames Wright   @param[in]      num_entries  Number of array entries
6340d80af1SJames Wright   @param[in,out]  array_petsc  Array of `PetscInt`
6440d80af1SJames Wright   @param[out]     array_ceed   Array of `CeedInt`
6540d80af1SJames Wright 
6640d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
6740d80af1SJames Wright **/
6840d80af1SJames Wright static inline PetscErrorCode IntArrayCeedToPetsc(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) {
6940d80af1SJames Wright   const CeedInt  int_c = 0;
7040d80af1SJames Wright   const PetscInt int_p = 0;
7140d80af1SJames Wright 
7240d80af1SJames Wright   PetscFunctionBeginUser;
7340d80af1SJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
7440d80af1SJames Wright     *array_petsc = (PetscInt *)*array_ceed;
7540d80af1SJames Wright   } else {
7640d80af1SJames Wright     *array_petsc = malloc(num_entries * sizeof(PetscInt));
7740d80af1SJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i];
7840d80af1SJames Wright     free(*array_ceed);
7940d80af1SJames Wright   }
8040d80af1SJames Wright   *array_ceed = NULL;
8140d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
8240d80af1SJames Wright }
8340d80af1SJames Wright 
8440d80af1SJames Wright /**
8540d80af1SJames Wright   @brief Translate array of `PetscInt` to `CeedInt`.
8640d80af1SJames Wright     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
8740d80af1SJames Wright     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
8840d80af1SJames Wright 
8940d80af1SJames Wright   Not collective across MPI processes.
9040d80af1SJames Wright 
9140d80af1SJames Wright   @param[in]      num_entries  Number of array entries
9240d80af1SJames Wright   @param[in,out]  array_petsc  Array of `PetscInt`
9340d80af1SJames Wright   @param[out]     array_ceed   Array of `CeedInt`
9440d80af1SJames Wright 
9540d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
9640d80af1SJames Wright **/
9740d80af1SJames Wright static inline PetscErrorCode IntArrayPetscToCeed(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) {
9840d80af1SJames Wright   const CeedInt  int_c = 0;
9940d80af1SJames Wright   const PetscInt int_p = 0;
10040d80af1SJames Wright 
10140d80af1SJames Wright   PetscFunctionBeginUser;
10240d80af1SJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
10340d80af1SJames Wright     *array_ceed = (CeedInt *)*array_petsc;
10440d80af1SJames Wright   } else {
10540d80af1SJames Wright     PetscCall(PetscMalloc1(num_entries, array_ceed));
10640d80af1SJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i];
10740d80af1SJames Wright     PetscCall(PetscFree(*array_petsc));
10840d80af1SJames Wright   }
10940d80af1SJames Wright   *array_petsc = NULL;
11040d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
11140d80af1SJames Wright }
11240d80af1SJames Wright 
11340d80af1SJames Wright /**
11440d80af1SJames Wright   @brief Transfer array from PETSc `Vec` to `CeedVector`.
11540d80af1SJames Wright 
11640d80af1SJames Wright   Collective across MPI processes.
11740d80af1SJames Wright 
11840d80af1SJames Wright   @param[in]   X_petsc   PETSc `Vec`
11940d80af1SJames Wright   @param[out]  mem_type  PETSc `MemType`
12040d80af1SJames Wright   @param[out]  x_ceed    `CeedVector`
12140d80af1SJames Wright 
12240d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
12340d80af1SJames Wright **/
12440d80af1SJames Wright static inline PetscErrorCode VecPetscToCeed(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed) {
12540d80af1SJames Wright   PetscScalar *x;
12640d80af1SJames Wright 
12740d80af1SJames Wright   PetscFunctionBeginUser;
12840d80af1SJames Wright   PetscCall(VecGetArrayAndMemType(X_petsc, &x, mem_type));
12940d80af1SJames Wright   PetscCallCeed(CeedVectorReturnCeed(x_ceed), CeedVectorSetArray(x_ceed, MemTypePetscToCeed(*mem_type), CEED_USE_POINTER, x));
13040d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
13140d80af1SJames Wright }
13240d80af1SJames Wright 
13340d80af1SJames Wright /**
13440d80af1SJames Wright   @brief Transfer array from `CeedVector` to PETSc `Vec`.
13540d80af1SJames Wright 
13640d80af1SJames Wright   Collective across MPI processes.
13740d80af1SJames Wright 
13840d80af1SJames Wright   @param[in]   x_ceed    `CeedVector`
13940d80af1SJames Wright   @param[in]   mem_type  PETSc `MemType`
14040d80af1SJames Wright   @param[out]  X_petsc   PETSc `Vec`
14140d80af1SJames Wright 
14240d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
14340d80af1SJames Wright **/
14440d80af1SJames Wright static inline PetscErrorCode VecCeedToPetsc(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc) {
14540d80af1SJames Wright   PetscScalar *x;
14640d80af1SJames Wright 
14740d80af1SJames Wright   PetscFunctionBeginUser;
14840d80af1SJames Wright   PetscCallCeed(CeedVectorReturnCeed(x_ceed), CeedVectorTakeArray(x_ceed, MemTypePetscToCeed(mem_type), &x));
14940d80af1SJames Wright   PetscCall(VecRestoreArrayAndMemType(X_petsc, &x));
15040d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
15140d80af1SJames Wright }
15240d80af1SJames Wright 
15340d80af1SJames Wright /**
15440d80af1SJames Wright   @brief Transfer read only array from PETSc `Vec` to `CeedVector`.
15540d80af1SJames Wright 
15640d80af1SJames Wright   Collective across MPI processes.
15740d80af1SJames Wright 
15840d80af1SJames Wright   @param[in]   X_petsc   PETSc `Vec`
15940d80af1SJames Wright   @param[out]  mem_type  PETSc `MemType`
16040d80af1SJames Wright   @param[out]  x_ceed    `CeedVector`
16140d80af1SJames Wright 
16240d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
16340d80af1SJames Wright **/
16440d80af1SJames Wright static inline PetscErrorCode VecReadPetscToCeed(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed) {
16540d80af1SJames Wright   PetscScalar *x;
16640d80af1SJames Wright 
16740d80af1SJames Wright   PetscFunctionBeginUser;
16840d80af1SJames Wright   PetscCall(VecGetArrayReadAndMemType(X_petsc, (const PetscScalar **)&x, mem_type));
16940d80af1SJames Wright   PetscCallCeed(CeedVectorReturnCeed(x_ceed), CeedVectorSetArray(x_ceed, MemTypePetscToCeed(*mem_type), CEED_USE_POINTER, x));
17040d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
17140d80af1SJames Wright }
17240d80af1SJames Wright 
17340d80af1SJames Wright /**
17440d80af1SJames Wright   @brief Transfer read only array from `CeedVector` to PETSc `Vec`.
17540d80af1SJames Wright 
17640d80af1SJames Wright   Collective across MPI processes.
17740d80af1SJames Wright 
17840d80af1SJames Wright   @param[in]   x_ceed    `CeedVector`
17940d80af1SJames Wright   @param[in]   mem_type  PETSc `MemType`
18040d80af1SJames Wright   @param[out]  X_petsc   PETSc `Vec`
18140d80af1SJames Wright 
18240d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
18340d80af1SJames Wright **/
18440d80af1SJames Wright static inline PetscErrorCode VecReadCeedToPetsc(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc) {
18540d80af1SJames Wright   PetscScalar *x;
18640d80af1SJames Wright 
18740d80af1SJames Wright   PetscFunctionBeginUser;
18840d80af1SJames Wright   PetscCallCeed(CeedVectorReturnCeed(x_ceed), CeedVectorTakeArray(x_ceed, MemTypePetscToCeed(mem_type), &x));
18940d80af1SJames Wright   PetscCall(VecRestoreArrayReadAndMemType(X_petsc, (const PetscScalar **)&x));
19040d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
19140d80af1SJames Wright }
19240d80af1SJames Wright 
19340d80af1SJames Wright /**
19440d80af1SJames Wright   @brief Copy PETSc `Vec` data into `CeedVector`
19540d80af1SJames Wright 
19640d80af1SJames Wright   @param[in]   X_petsc PETSc `Vec`
19740d80af1SJames Wright   @param[out]  x_ceed  `CeedVector`
19840d80af1SJames Wright 
19940d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
20040d80af1SJames Wright **/
20140d80af1SJames Wright static inline PetscErrorCode VecCopyPetscToCeed(Vec X_petsc, CeedVector x_ceed) {
20240d80af1SJames Wright   PetscScalar *x;
20340d80af1SJames Wright   PetscMemType mem_type;
20440d80af1SJames Wright   PetscInt     X_size;
20540d80af1SJames Wright   CeedSize     x_size;
20640d80af1SJames Wright   Ceed         ceed;
20740d80af1SJames Wright 
20840d80af1SJames Wright   PetscFunctionBeginUser;
20940d80af1SJames Wright   PetscCall(CeedVectorGetCeed(x_ceed, &ceed));
21040d80af1SJames Wright   PetscCall(VecGetLocalSize(X_petsc, &X_size));
21140d80af1SJames Wright   PetscCallCeed(ceed, CeedVectorGetLength(x_ceed, &x_size));
21240d80af1SJames Wright   PetscCheck(X_size == x_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "X_petsc (%" PetscInt_FMT ") and x_ceed (%" CeedSize_FMT ") must be same size",
21340d80af1SJames Wright              X_size, x_size);
21440d80af1SJames Wright 
21540d80af1SJames Wright   PetscCall(VecGetArrayReadAndMemType(X_petsc, (const PetscScalar **)&x, &mem_type));
21640d80af1SJames Wright   PetscCallCeed(ceed, CeedVectorSetArray(x_ceed, MemTypePetscToCeed(mem_type), CEED_COPY_VALUES, x));
21740d80af1SJames Wright   PetscCall(VecRestoreArrayReadAndMemType(X_petsc, (const PetscScalar **)&x));
21840d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
21940d80af1SJames Wright }
22040d80af1SJames Wright 
22140d80af1SJames Wright /**
22240d80af1SJames Wright   @brief Return the quadrature size from the eval mode, dimension, and number of components
22340d80af1SJames Wright 
22440d80af1SJames Wright   @param[in]  eval_mode       The basis evaluation mode
22540d80af1SJames Wright   @param[in]  dim             The basis dimension
22640d80af1SJames Wright   @param[in]  num_components  The basis number of components
22740d80af1SJames Wright 
22840d80af1SJames Wright   @return The maximum of the two integers
22940d80af1SJames Wright 
23040d80af1SJames Wright **/
23140d80af1SJames Wright static inline CeedInt GetCeedQuadratureSize(CeedEvalMode eval_mode, CeedInt dim, CeedInt num_components) {
23240d80af1SJames Wright   switch (eval_mode) {
23340d80af1SJames Wright     case CEED_EVAL_INTERP:
23440d80af1SJames Wright       return num_components;
23540d80af1SJames Wright     case CEED_EVAL_GRAD:
23640d80af1SJames Wright       return dim * num_components;
23740d80af1SJames Wright     default:
23840d80af1SJames Wright       return -1;
23940d80af1SJames Wright   }
24040d80af1SJames Wright }
24140d80af1SJames Wright 
24240d80af1SJames Wright /**
24340d80af1SJames Wright   @brief Convert from DMPolytopeType to CeedElemTopology
24440d80af1SJames Wright 
24540d80af1SJames Wright   @param[in]  cell_type  DMPolytopeType for the cell
24640d80af1SJames Wright 
24740d80af1SJames Wright   @return CeedElemTopology, or 0 if no equivelent CeedElemTopology was found
24840d80af1SJames Wright **/
24940d80af1SJames Wright static inline CeedElemTopology PolytopeTypePetscToCeed(DMPolytopeType cell_type) {
25040d80af1SJames Wright   switch (cell_type) {
25140d80af1SJames Wright     case DM_POLYTOPE_TRIANGLE:
25240d80af1SJames Wright       return CEED_TOPOLOGY_TRIANGLE;
25340d80af1SJames Wright     case DM_POLYTOPE_QUADRILATERAL:
25440d80af1SJames Wright       return CEED_TOPOLOGY_QUAD;
25540d80af1SJames Wright     case DM_POLYTOPE_TETRAHEDRON:
25640d80af1SJames Wright       return CEED_TOPOLOGY_TET;
25740d80af1SJames Wright     case DM_POLYTOPE_HEXAHEDRON:
25840d80af1SJames Wright       return CEED_TOPOLOGY_HEX;
25940d80af1SJames Wright     default:
26040d80af1SJames Wright       return 0;
26140d80af1SJames Wright   }
26240d80af1SJames Wright }
263