xref: /libCEED/examples/solids/include/utils.h (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
33d8e8822SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
53d8e8822SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
724a65d3dSJeremy L Thompson #pragma once
85754ecacSJeremy L Thompson 
95754ecacSJeremy L Thompson #include <ceed.h>
1049aac155SJeremy L Thompson #include <petscsys.h>
115754ecacSJeremy L Thompson 
125754ecacSJeremy L Thompson // Translate PetscMemType to CeedMemType
MemTypeP2C(PetscMemType mem_type)132b730f8bSJeremy L Thompson static inline CeedMemType MemTypeP2C(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; }
144d00b080SJeremy L Thompson 
154d00b080SJeremy L Thompson /**
164d00b080SJeremy L Thompson   @brief Translate array of `PetscInt` to `CeedInt`.
174d00b080SJeremy L Thompson     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
184d00b080SJeremy L Thompson     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
194d00b080SJeremy L Thompson 
204d00b080SJeremy L Thompson   Not collective across MPI processes.
214d00b080SJeremy L Thompson 
224d00b080SJeremy L Thompson   @param[in]      num_entries  Number of array entries
234d00b080SJeremy L Thompson   @param[in,out]  array_petsc  Array of `PetscInt`
244d00b080SJeremy L Thompson   @param[out]     array_ceed   Array of `CeedInt`
254d00b080SJeremy L Thompson 
264d00b080SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
274d00b080SJeremy L Thompson **/
IntArrayCeedToPetsc(PetscInt num_entries,CeedInt ** array_ceed,PetscInt ** array_petsc)284d00b080SJeremy L Thompson static inline PetscErrorCode IntArrayCeedToPetsc(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) {
294d00b080SJeremy L Thompson   const CeedInt  int_c = 0;
304d00b080SJeremy L Thompson   const PetscInt int_p = 0;
314d00b080SJeremy L Thompson 
324d00b080SJeremy L Thompson   PetscFunctionBeginUser;
334d00b080SJeremy L Thompson   if (sizeof(int_c) == sizeof(int_p)) {
344d00b080SJeremy L Thompson     *array_petsc = (PetscInt *)*array_ceed;
354d00b080SJeremy L Thompson   } else {
364d00b080SJeremy L Thompson     *array_petsc = malloc(num_entries * sizeof(PetscInt));
374d00b080SJeremy L Thompson     for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i];
384d00b080SJeremy L Thompson     free(*array_ceed);
394d00b080SJeremy L Thompson   }
404d00b080SJeremy L Thompson   *array_ceed = NULL;
414d00b080SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
424d00b080SJeremy L Thompson }
434d00b080SJeremy L Thompson 
444d00b080SJeremy L Thompson /**
454d00b080SJeremy L Thompson   @brief Translate array of `PetscInt` to `CeedInt`.
464d00b080SJeremy L Thompson     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
474d00b080SJeremy L Thompson     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
484d00b080SJeremy L Thompson 
494d00b080SJeremy L Thompson   Not collective across MPI processes.
504d00b080SJeremy L Thompson 
514d00b080SJeremy L Thompson   @param[in]      num_entries  Number of array entries
524d00b080SJeremy L Thompson   @param[in,out]  array_petsc  Array of `PetscInt`
534d00b080SJeremy L Thompson   @param[out]     array_ceed   Array of `CeedInt`
544d00b080SJeremy L Thompson 
554d00b080SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
564d00b080SJeremy L Thompson **/
IntArrayPetscToCeed(PetscInt num_entries,PetscInt ** array_petsc,CeedInt ** array_ceed)574d00b080SJeremy L Thompson static inline PetscErrorCode IntArrayPetscToCeed(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) {
584d00b080SJeremy L Thompson   const CeedInt  int_c = 0;
594d00b080SJeremy L Thompson   const PetscInt int_p = 0;
604d00b080SJeremy L Thompson 
614d00b080SJeremy L Thompson   PetscFunctionBeginUser;
624d00b080SJeremy L Thompson   if (sizeof(int_c) == sizeof(int_p)) {
634d00b080SJeremy L Thompson     *array_ceed = (CeedInt *)*array_petsc;
644d00b080SJeremy L Thompson   } else {
654d00b080SJeremy L Thompson     PetscCall(PetscMalloc1(num_entries, array_ceed));
664d00b080SJeremy L Thompson     for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i];
674d00b080SJeremy L Thompson     PetscCall(PetscFree(*array_petsc));
684d00b080SJeremy L Thompson   }
694d00b080SJeremy L Thompson   *array_petsc = NULL;
704d00b080SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
714d00b080SJeremy L Thompson }
72