xref: /libCEED/examples/fluids/include/petsc-ceed-utils.h (revision 9ba83ac0e4b1fca39d6fa6737a318a9f0cbc172d)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
25037d55dSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
35037d55dSJames Wright //
45037d55dSJames Wright // SPDX-License-Identifier: BSD-2-Clause
55037d55dSJames Wright //
65037d55dSJames Wright // This file is part of CEED:  http://github.com/ceed
75037d55dSJames Wright #pragma once
85037d55dSJames Wright 
95037d55dSJames Wright #include <ceed.h>
105037d55dSJames Wright #include <petscdm.h>
115037d55dSJames Wright 
125037d55dSJames Wright /**
135037d55dSJames Wright   @brief Copy the reference to a `Vec`.
145037d55dSJames Wright          Note: If `vec_copy` is non-null, it is assumed to be a valid pointer to a `Vec` and `VecDestroy()` will be called.
155037d55dSJames Wright 
165037d55dSJames Wright   Collective across MPI processes.
175037d55dSJames Wright 
185037d55dSJames Wright   @param[in]   vec       `Vec` to reference
195037d55dSJames Wright   @param[out]  vec_copy  Copy of reference
205037d55dSJames Wright 
215037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
225037d55dSJames Wright **/
235037d55dSJames Wright static inline PetscErrorCode VecReferenceCopy(Vec vec, Vec *vec_copy) {
245037d55dSJames Wright   PetscFunctionBeginUser;
255037d55dSJames Wright   PetscCall(PetscObjectReference((PetscObject)vec));
265037d55dSJames Wright   PetscCall(VecDestroy(vec_copy));
275037d55dSJames Wright   *vec_copy = vec;
285037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
295037d55dSJames Wright }
305037d55dSJames Wright 
315037d55dSJames Wright /**
325037d55dSJames Wright   @brief Copy the reference to a `DM`.
335037d55dSJames Wright          Note: If `dm_copy` is non-null, it is assumed to be a valid pointer to a `DM` and `DMDestroy()` will be called.
345037d55dSJames Wright 
355037d55dSJames Wright   Collective across MPI processes.
365037d55dSJames Wright 
375037d55dSJames Wright   @param[in]   dm       `DM` to reference
385037d55dSJames Wright   @param[out]  dm_copy  Copy of reference
395037d55dSJames Wright 
405037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
415037d55dSJames Wright **/
425037d55dSJames Wright static inline PetscErrorCode DMReferenceCopy(DM dm, DM *dm_copy) {
435037d55dSJames Wright   PetscFunctionBeginUser;
445037d55dSJames Wright   PetscCall(PetscObjectReference((PetscObject)dm));
455037d55dSJames Wright   PetscCall(DMDestroy(dm_copy));
465037d55dSJames Wright   *dm_copy = dm;
475037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
485037d55dSJames Wright }
495037d55dSJames Wright 
505037d55dSJames Wright /**
515037d55dSJames Wright   @brief Translate PetscMemType to CeedMemType
525037d55dSJames Wright 
535037d55dSJames Wright   @param[in]  mem_type  PetscMemType
545037d55dSJames Wright 
555037d55dSJames Wright   @return Equivalent CeedMemType
565037d55dSJames Wright **/
575037d55dSJames Wright /// @ingroup RatelInternal
585037d55dSJames Wright static inline CeedMemType MemTypePetscToCeed(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; }
595037d55dSJames Wright 
605037d55dSJames Wright /**
615037d55dSJames Wright   @brief Translate array of `PetscInt` to `CeedInt`.
625037d55dSJames Wright     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
635037d55dSJames Wright     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
645037d55dSJames Wright 
655037d55dSJames Wright   Not collective across MPI processes.
665037d55dSJames Wright 
675037d55dSJames Wright   @param[in]      num_entries  Number of array entries
685037d55dSJames Wright   @param[in,out]  array_petsc  Array of `PetscInt`
695037d55dSJames Wright   @param[out]     array_ceed   Array of `CeedInt`
705037d55dSJames Wright 
715037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
725037d55dSJames Wright **/
735037d55dSJames Wright static inline PetscErrorCode IntArrayCeedToPetsc(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) {
745037d55dSJames Wright   const CeedInt  int_c = 0;
755037d55dSJames Wright   const PetscInt int_p = 0;
765037d55dSJames Wright 
775037d55dSJames Wright   PetscFunctionBeginUser;
785037d55dSJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
795037d55dSJames Wright     *array_petsc = (PetscInt *)*array_ceed;
805037d55dSJames Wright   } else {
815037d55dSJames Wright     *array_petsc = malloc(num_entries * sizeof(PetscInt));
825037d55dSJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i];
835037d55dSJames Wright     free(*array_ceed);
845037d55dSJames Wright   }
855037d55dSJames Wright   *array_ceed = NULL;
865037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
875037d55dSJames Wright }
885037d55dSJames Wright 
895037d55dSJames Wright /**
905037d55dSJames Wright   @brief Translate array of `PetscInt` to `CeedInt`.
915037d55dSJames Wright     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
925037d55dSJames Wright     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
935037d55dSJames Wright 
945037d55dSJames Wright   Not collective across MPI processes.
955037d55dSJames Wright 
965037d55dSJames Wright   @param[in]      num_entries  Number of array entries
975037d55dSJames Wright   @param[in,out]  array_petsc  Array of `PetscInt`
985037d55dSJames Wright   @param[out]     array_ceed   Array of `CeedInt`
995037d55dSJames Wright 
1005037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
1015037d55dSJames Wright **/
1025037d55dSJames Wright static inline PetscErrorCode IntArrayPetscToCeed(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) {
1035037d55dSJames Wright   const CeedInt  int_c = 0;
1045037d55dSJames Wright   const PetscInt int_p = 0;
1055037d55dSJames Wright 
1065037d55dSJames Wright   PetscFunctionBeginUser;
1075037d55dSJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
1085037d55dSJames Wright     *array_ceed = (CeedInt *)*array_petsc;
1095037d55dSJames Wright   } else {
1105037d55dSJames Wright     PetscCall(PetscMalloc1(num_entries, array_ceed));
1115037d55dSJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i];
1125037d55dSJames Wright     PetscCall(PetscFree(*array_petsc));
1135037d55dSJames Wright   }
1145037d55dSJames Wright   *array_petsc = NULL;
1155037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1165037d55dSJames Wright }
1175037d55dSJames Wright 
1185037d55dSJames Wright /**
1195037d55dSJames Wright   @brief Transfer array from PETSc `Vec` to `CeedVector`.
1205037d55dSJames Wright 
1215037d55dSJames Wright   Collective across MPI processes.
1225037d55dSJames Wright 
1235037d55dSJames Wright   @param[in]   X_petsc   PETSc `Vec`
1245037d55dSJames Wright   @param[out]  mem_type  PETSc `MemType`
1255037d55dSJames Wright   @param[out]  x_ceed    `CeedVector`
1265037d55dSJames Wright 
1275037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
1285037d55dSJames Wright **/
1295037d55dSJames Wright static inline PetscErrorCode VecPetscToCeed(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed) {
1305037d55dSJames Wright   PetscScalar *x;
1315037d55dSJames Wright 
1325037d55dSJames Wright   PetscFunctionBeginUser;
1335037d55dSJames Wright   PetscCall(VecGetArrayAndMemType(X_petsc, &x, mem_type));
1345037d55dSJames Wright   PetscCallCeed(CeedVectorReturnCeed(x_ceed), CeedVectorSetArray(x_ceed, MemTypePetscToCeed(*mem_type), CEED_USE_POINTER, x));
1355037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1365037d55dSJames Wright }
1375037d55dSJames Wright 
1385037d55dSJames Wright /**
1395037d55dSJames Wright   @brief Transfer array from `CeedVector` to PETSc `Vec`.
1405037d55dSJames Wright 
1415037d55dSJames Wright   Collective across MPI processes.
1425037d55dSJames Wright 
1435037d55dSJames Wright   @param[in]   x_ceed    `CeedVector`
1445037d55dSJames Wright   @param[in]   mem_type  PETSc `MemType`
1455037d55dSJames Wright   @param[out]  X_petsc   PETSc `Vec`
1465037d55dSJames Wright 
1475037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
1485037d55dSJames Wright **/
1495037d55dSJames Wright static inline PetscErrorCode VecCeedToPetsc(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc) {
1505037d55dSJames Wright   PetscScalar *x;
1515037d55dSJames Wright 
1525037d55dSJames Wright   PetscFunctionBeginUser;
1535037d55dSJames Wright   PetscCallCeed(CeedVectorReturnCeed(x_ceed), CeedVectorTakeArray(x_ceed, MemTypePetscToCeed(mem_type), &x));
1545037d55dSJames Wright   PetscCall(VecRestoreArrayAndMemType(X_petsc, &x));
1555037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1565037d55dSJames Wright }
1575037d55dSJames Wright 
1585037d55dSJames Wright /**
1595037d55dSJames Wright   @brief Transfer read only array from PETSc `Vec` to `CeedVector`.
1605037d55dSJames Wright 
1615037d55dSJames Wright   Collective across MPI processes.
1625037d55dSJames Wright 
1635037d55dSJames Wright   @param[in]   X_petsc   PETSc `Vec`
1645037d55dSJames Wright   @param[out]  mem_type  PETSc `MemType`
1655037d55dSJames Wright   @param[out]  x_ceed    `CeedVector`
1665037d55dSJames Wright 
1675037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
1685037d55dSJames Wright **/
1695037d55dSJames Wright static inline PetscErrorCode VecReadPetscToCeed(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed) {
1705037d55dSJames Wright   PetscScalar *x;
1715037d55dSJames Wright 
1725037d55dSJames Wright   PetscFunctionBeginUser;
1735037d55dSJames Wright   PetscCall(VecGetArrayReadAndMemType(X_petsc, (const PetscScalar **)&x, mem_type));
1745037d55dSJames Wright   PetscCallCeed(CeedVectorReturnCeed(x_ceed), CeedVectorSetArray(x_ceed, MemTypePetscToCeed(*mem_type), CEED_USE_POINTER, x));
1755037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1765037d55dSJames Wright }
1775037d55dSJames Wright 
1785037d55dSJames Wright /**
1795037d55dSJames Wright   @brief Transfer read only array from `CeedVector` to PETSc `Vec`.
1805037d55dSJames Wright 
1815037d55dSJames Wright   Collective across MPI processes.
1825037d55dSJames Wright 
1835037d55dSJames Wright   @param[in]   x_ceed    `CeedVector`
1845037d55dSJames Wright   @param[in]   mem_type  PETSc `MemType`
1855037d55dSJames Wright   @param[out]  X_petsc   PETSc `Vec`
1865037d55dSJames Wright 
1875037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
1885037d55dSJames Wright **/
1895037d55dSJames Wright static inline PetscErrorCode VecReadCeedToPetsc(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc) {
1905037d55dSJames Wright   PetscScalar *x;
1915037d55dSJames Wright 
1925037d55dSJames Wright   PetscFunctionBeginUser;
1935037d55dSJames Wright   PetscCallCeed(CeedVectorReturnCeed(x_ceed), CeedVectorTakeArray(x_ceed, MemTypePetscToCeed(mem_type), &x));
1945037d55dSJames Wright   PetscCall(VecRestoreArrayReadAndMemType(X_petsc, (const PetscScalar **)&x));
1955037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1965037d55dSJames Wright }
1975037d55dSJames Wright 
1985037d55dSJames Wright /**
1995037d55dSJames Wright   @brief Copy PETSc `Vec` data into `CeedVector`
2005037d55dSJames Wright 
2015037d55dSJames Wright   @param[in]   X_petsc PETSc `Vec`
2025037d55dSJames Wright   @param[out]  x_ceed  `CeedVector`
2035037d55dSJames Wright 
2045037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
2055037d55dSJames Wright **/
2065037d55dSJames Wright static inline PetscErrorCode VecCopyPetscToCeed(Vec X_petsc, CeedVector x_ceed) {
2075037d55dSJames Wright   PetscScalar *x;
2085037d55dSJames Wright   PetscMemType mem_type;
2095037d55dSJames Wright   PetscInt     X_size;
2105037d55dSJames Wright   CeedSize     x_size;
2115037d55dSJames Wright   Ceed         ceed;
2125037d55dSJames Wright 
2135037d55dSJames Wright   PetscFunctionBeginUser;
2145037d55dSJames Wright   PetscCall(CeedVectorGetCeed(x_ceed, &ceed));
2155037d55dSJames Wright   PetscCall(VecGetLocalSize(X_petsc, &X_size));
2165037d55dSJames Wright   PetscCallCeed(ceed, CeedVectorGetLength(x_ceed, &x_size));
2175037d55dSJames 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",
2185037d55dSJames Wright              X_size, x_size);
2195037d55dSJames Wright 
2205037d55dSJames Wright   PetscCall(VecGetArrayReadAndMemType(X_petsc, (const PetscScalar **)&x, &mem_type));
2215037d55dSJames Wright   PetscCallCeed(ceed, CeedVectorSetArray(x_ceed, MemTypePetscToCeed(mem_type), CEED_COPY_VALUES, x));
2225037d55dSJames Wright   PetscCall(VecRestoreArrayReadAndMemType(X_petsc, (const PetscScalar **)&x));
2239bc66399SJeremy L Thompson   PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PetscObjectComm((PetscObject)X_petsc), PETSC_ERR_LIB, "Destroying Ceed object failed");
2245037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2255037d55dSJames Wright }
2265037d55dSJames Wright 
2275037d55dSJames Wright /**
2285037d55dSJames Wright   @brief Return the quadrature size from the eval mode, dimension, and number of components
2295037d55dSJames Wright 
2305037d55dSJames Wright   @param[in]  eval_mode       The basis evaluation mode
2315037d55dSJames Wright   @param[in]  dim             The basis dimension
2325037d55dSJames Wright   @param[in]  num_components  The basis number of components
2335037d55dSJames Wright 
2345037d55dSJames Wright   @return The maximum of the two integers
2355037d55dSJames Wright 
2365037d55dSJames Wright **/
2375037d55dSJames Wright /// @ingroup RatelInternal
2385037d55dSJames Wright static inline CeedInt GetCeedQuadratureSize(CeedEvalMode eval_mode, CeedInt dim, CeedInt num_components) {
2395037d55dSJames Wright   switch (eval_mode) {
2405037d55dSJames Wright     case CEED_EVAL_INTERP:
2415037d55dSJames Wright       return num_components;
2425037d55dSJames Wright     case CEED_EVAL_GRAD:
2435037d55dSJames Wright       return dim * num_components;
2445037d55dSJames Wright     default:
2455037d55dSJames Wright       return -1;
2465037d55dSJames Wright   }
2475037d55dSJames Wright }
2485037d55dSJames Wright 
2495037d55dSJames Wright /**
2505037d55dSJames Wright   @brief Convert from DMPolytopeType to CeedElemTopology
2515037d55dSJames Wright 
2525037d55dSJames Wright   @param[in]  cell_type  DMPolytopeType for the cell
2535037d55dSJames Wright 
2545037d55dSJames Wright   @return CeedElemTopology, or 0 if no equivelent CeedElemTopology was found
2555037d55dSJames Wright **/
2565037d55dSJames Wright static inline CeedElemTopology PolytopeTypePetscToCeed(DMPolytopeType cell_type) {
2575037d55dSJames Wright   switch (cell_type) {
2585037d55dSJames Wright     case DM_POLYTOPE_TRIANGLE:
2595037d55dSJames Wright       return CEED_TOPOLOGY_TRIANGLE;
2605037d55dSJames Wright     case DM_POLYTOPE_QUADRILATERAL:
2615037d55dSJames Wright       return CEED_TOPOLOGY_QUAD;
2625037d55dSJames Wright     case DM_POLYTOPE_TETRAHEDRON:
2635037d55dSJames Wright       return CEED_TOPOLOGY_TET;
2645037d55dSJames Wright     case DM_POLYTOPE_HEXAHEDRON:
2655037d55dSJames Wright       return CEED_TOPOLOGY_HEX;
2665037d55dSJames Wright     default:
2675037d55dSJames Wright       return 0;
2685037d55dSJames Wright   }
2695037d55dSJames Wright }
270