xref: /libCEED/examples/petsc/include/petscutils.h (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
298285ab4SZach Atkins // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
398285ab4SZach Atkins //
498285ab4SZach Atkins // SPDX-License-Identifier: BSD-2-Clause
598285ab4SZach Atkins //
698285ab4SZach Atkins // This file is part of CEED:  http://github.com/ceed
798285ab4SZach Atkins 
898285ab4SZach Atkins /// @file
998285ab4SZach Atkins /// Utility functions for using PETSc with libCEED
1024a65d3dSJeremy L Thompson #pragma once
11e83e87a5Sjeremylt 
12e83e87a5Sjeremylt #include <ceed.h>
13e83e87a5Sjeremylt #include <petsc.h>
142b730f8bSJeremy L Thompson 
15de1229c5Srezgarshakeri #include "structs.h"
16e83e87a5Sjeremylt 
17e83e87a5Sjeremylt CeedMemType      MemTypeP2C(PetscMemType mtype);
18b6972d74SZach Atkins PetscErrorCode   VecP2C(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed);
19b6972d74SZach Atkins PetscErrorCode   VecC2P(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc);
20b6972d74SZach Atkins PetscErrorCode   VecReadP2C(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed);
21b6972d74SZach Atkins PetscErrorCode   VecReadC2P(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc);
229b072555Sjeremylt PetscErrorCode   Kershaw(DM dm_orig, PetscScalar eps);
232b730f8bSJeremy L Thompson PetscErrorCode   SetupDMByDegree(DM dm, PetscInt p_degree, PetscInt q_extra, PetscInt num_comp_u, PetscInt topo_dim, bool enforce_bc);
242b730f8bSJeremy L Thompson PetscErrorCode   CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr);
25de1229c5Srezgarshakeri CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type);
262b730f8bSJeremy L Thompson PetscErrorCode   DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field);
272b730f8bSJeremy L Thompson PetscErrorCode   BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe,
282b730f8bSJeremy L Thompson                                            PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis);
292b730f8bSJeremy L Thompson PetscErrorCode   CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, BPData bp_data,
30f755c37aSrezgarshakeri                                      CeedBasis *basis);
31de1229c5Srezgarshakeri PetscErrorCode   CreateDistributedDM(RunParams rp, DM *dm);
324d00b080SJeremy L Thompson 
334d00b080SJeremy L Thompson /**
344d00b080SJeremy L Thompson   @brief Translate array of `PetscInt` to `CeedInt`.
354d00b080SJeremy L Thompson     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
364d00b080SJeremy L Thompson     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
374d00b080SJeremy L Thompson 
384d00b080SJeremy L Thompson   Not collective across MPI processes.
394d00b080SJeremy L Thompson 
404d00b080SJeremy L Thompson   @param[in]      num_entries  Number of array entries
414d00b080SJeremy L Thompson   @param[in,out]  array_petsc  Array of `PetscInt`
424d00b080SJeremy L Thompson   @param[out]     array_ceed   Array of `CeedInt`
434d00b080SJeremy L Thompson 
444d00b080SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
454d00b080SJeremy L Thompson **/
IntArrayCeedToPetsc(PetscInt num_entries,CeedInt ** array_ceed,PetscInt ** array_petsc)464d00b080SJeremy L Thompson static inline PetscErrorCode IntArrayCeedToPetsc(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) {
474d00b080SJeremy L Thompson   const CeedInt  int_c = 0;
484d00b080SJeremy L Thompson   const PetscInt int_p = 0;
494d00b080SJeremy L Thompson 
504d00b080SJeremy L Thompson   PetscFunctionBeginUser;
514d00b080SJeremy L Thompson   if (sizeof(int_c) == sizeof(int_p)) {
524d00b080SJeremy L Thompson     *array_petsc = (PetscInt *)*array_ceed;
534d00b080SJeremy L Thompson   } else {
544d00b080SJeremy L Thompson     *array_petsc = malloc(num_entries * sizeof(PetscInt));
554d00b080SJeremy L Thompson     for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i];
564d00b080SJeremy L Thompson     free(*array_ceed);
574d00b080SJeremy L Thompson   }
584d00b080SJeremy L Thompson   *array_ceed = NULL;
594d00b080SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
604d00b080SJeremy L Thompson }
614d00b080SJeremy L Thompson 
624d00b080SJeremy L Thompson /**
634d00b080SJeremy L Thompson   @brief Translate array of `PetscInt` to `CeedInt`.
644d00b080SJeremy L Thompson     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
654d00b080SJeremy L Thompson     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
664d00b080SJeremy L Thompson 
674d00b080SJeremy L Thompson   Not collective across MPI processes.
684d00b080SJeremy L Thompson 
694d00b080SJeremy L Thompson   @param[in]      num_entries  Number of array entries
704d00b080SJeremy L Thompson   @param[in,out]  array_petsc  Array of `PetscInt`
714d00b080SJeremy L Thompson   @param[out]     array_ceed   Array of `CeedInt`
724d00b080SJeremy L Thompson 
734d00b080SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
744d00b080SJeremy L Thompson **/
IntArrayPetscToCeed(PetscInt num_entries,PetscInt ** array_petsc,CeedInt ** array_ceed)754d00b080SJeremy L Thompson static inline PetscErrorCode IntArrayPetscToCeed(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) {
764d00b080SJeremy L Thompson   const CeedInt  int_c = 0;
774d00b080SJeremy L Thompson   const PetscInt int_p = 0;
784d00b080SJeremy L Thompson 
794d00b080SJeremy L Thompson   PetscFunctionBeginUser;
804d00b080SJeremy L Thompson   if (sizeof(int_c) == sizeof(int_p)) {
814d00b080SJeremy L Thompson     *array_ceed = (CeedInt *)*array_petsc;
824d00b080SJeremy L Thompson   } else {
834d00b080SJeremy L Thompson     PetscCall(PetscMalloc1(num_entries, array_ceed));
844d00b080SJeremy L Thompson     for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i];
854d00b080SJeremy L Thompson     PetscCall(PetscFree(*array_petsc));
864d00b080SJeremy L Thompson   }
874d00b080SJeremy L Thompson   *array_petsc = NULL;
884d00b080SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
894d00b080SJeremy L Thompson }
90