xref: /libCEED/examples/petsc/include/petscutils.h (revision 4d00b080eb3f95d2e04e55c0ff369c5c847bb288)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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"
1649a40c8aSKenneth E. Jansen #if PETSC_VERSION_LT(3, 21, 0)
1749a40c8aSKenneth E. Jansen #define DMSetCoordinateDisc(a, b, c) DMProjectCoordinates(a, b)
1849a40c8aSKenneth E. Jansen #endif
19e83e87a5Sjeremylt 
20e83e87a5Sjeremylt CeedMemType      MemTypeP2C(PetscMemType mtype);
21b6972d74SZach Atkins PetscErrorCode   VecP2C(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed);
22b6972d74SZach Atkins PetscErrorCode   VecC2P(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc);
23b6972d74SZach Atkins PetscErrorCode   VecReadP2C(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed);
24b6972d74SZach Atkins PetscErrorCode   VecReadC2P(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc);
259b072555Sjeremylt PetscErrorCode   Kershaw(DM dm_orig, PetscScalar eps);
262b730f8bSJeremy L Thompson PetscErrorCode   SetupDMByDegree(DM dm, PetscInt p_degree, PetscInt q_extra, PetscInt num_comp_u, PetscInt topo_dim, bool enforce_bc);
272b730f8bSJeremy L Thompson PetscErrorCode   CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr);
28de1229c5Srezgarshakeri CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type);
292b730f8bSJeremy L Thompson PetscErrorCode   DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field);
302b730f8bSJeremy L Thompson PetscErrorCode   BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe,
312b730f8bSJeremy L Thompson                                            PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis);
322b730f8bSJeremy L Thompson PetscErrorCode   CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, BPData bp_data,
33f755c37aSrezgarshakeri                                      CeedBasis *basis);
34de1229c5Srezgarshakeri PetscErrorCode   CreateDistributedDM(RunParams rp, DM *dm);
35*4d00b080SJeremy L Thompson 
36*4d00b080SJeremy L Thompson /**
37*4d00b080SJeremy L Thompson   @brief Translate array of `PetscInt` to `CeedInt`.
38*4d00b080SJeremy L Thompson     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
39*4d00b080SJeremy L Thompson     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
40*4d00b080SJeremy L Thompson 
41*4d00b080SJeremy L Thompson   Not collective across MPI processes.
42*4d00b080SJeremy L Thompson 
43*4d00b080SJeremy L Thompson   @param[in]      num_entries  Number of array entries
44*4d00b080SJeremy L Thompson   @param[in,out]  array_petsc  Array of `PetscInt`
45*4d00b080SJeremy L Thompson   @param[out]     array_ceed   Array of `CeedInt`
46*4d00b080SJeremy L Thompson 
47*4d00b080SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
48*4d00b080SJeremy L Thompson **/
49*4d00b080SJeremy L Thompson static inline PetscErrorCode IntArrayCeedToPetsc(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) {
50*4d00b080SJeremy L Thompson   const CeedInt  int_c = 0;
51*4d00b080SJeremy L Thompson   const PetscInt int_p = 0;
52*4d00b080SJeremy L Thompson 
53*4d00b080SJeremy L Thompson   PetscFunctionBeginUser;
54*4d00b080SJeremy L Thompson   if (sizeof(int_c) == sizeof(int_p)) {
55*4d00b080SJeremy L Thompson     *array_petsc = (PetscInt *)*array_ceed;
56*4d00b080SJeremy L Thompson   } else {
57*4d00b080SJeremy L Thompson     *array_petsc = malloc(num_entries * sizeof(PetscInt));
58*4d00b080SJeremy L Thompson     for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i];
59*4d00b080SJeremy L Thompson     free(*array_ceed);
60*4d00b080SJeremy L Thompson   }
61*4d00b080SJeremy L Thompson   *array_ceed = NULL;
62*4d00b080SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
63*4d00b080SJeremy L Thompson }
64*4d00b080SJeremy L Thompson 
65*4d00b080SJeremy L Thompson /**
66*4d00b080SJeremy L Thompson   @brief Translate array of `PetscInt` to `CeedInt`.
67*4d00b080SJeremy L Thompson     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
68*4d00b080SJeremy L Thompson     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
69*4d00b080SJeremy L Thompson 
70*4d00b080SJeremy L Thompson   Not collective across MPI processes.
71*4d00b080SJeremy L Thompson 
72*4d00b080SJeremy L Thompson   @param[in]      num_entries  Number of array entries
73*4d00b080SJeremy L Thompson   @param[in,out]  array_petsc  Array of `PetscInt`
74*4d00b080SJeremy L Thompson   @param[out]     array_ceed   Array of `CeedInt`
75*4d00b080SJeremy L Thompson 
76*4d00b080SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
77*4d00b080SJeremy L Thompson **/
78*4d00b080SJeremy L Thompson static inline PetscErrorCode IntArrayPetscToCeed(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) {
79*4d00b080SJeremy L Thompson   const CeedInt  int_c = 0;
80*4d00b080SJeremy L Thompson   const PetscInt int_p = 0;
81*4d00b080SJeremy L Thompson 
82*4d00b080SJeremy L Thompson   PetscFunctionBeginUser;
83*4d00b080SJeremy L Thompson   if (sizeof(int_c) == sizeof(int_p)) {
84*4d00b080SJeremy L Thompson     *array_ceed = (CeedInt *)*array_petsc;
85*4d00b080SJeremy L Thompson   } else {
86*4d00b080SJeremy L Thompson     PetscCall(PetscMalloc1(num_entries, array_ceed));
87*4d00b080SJeremy L Thompson     for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i];
88*4d00b080SJeremy L Thompson     PetscCall(PetscFree(*array_petsc));
89*4d00b080SJeremy L Thompson   }
90*4d00b080SJeremy L Thompson   *array_petsc = NULL;
91*4d00b080SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
92*4d00b080SJeremy L Thompson }
93