xref: /libCEED/examples/petsc/include/petscutils.h (revision 4d00b080eb3f95d2e04e55c0ff369c5c847bb288)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Utility functions for using PETSc with libCEED
10 #pragma once
11 
12 #include <ceed.h>
13 #include <petsc.h>
14 
15 #include "structs.h"
16 #if PETSC_VERSION_LT(3, 21, 0)
17 #define DMSetCoordinateDisc(a, b, c) DMProjectCoordinates(a, b)
18 #endif
19 
20 CeedMemType      MemTypeP2C(PetscMemType mtype);
21 PetscErrorCode   VecP2C(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed);
22 PetscErrorCode   VecC2P(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc);
23 PetscErrorCode   VecReadP2C(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed);
24 PetscErrorCode   VecReadC2P(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc);
25 PetscErrorCode   Kershaw(DM dm_orig, PetscScalar eps);
26 PetscErrorCode   SetupDMByDegree(DM dm, PetscInt p_degree, PetscInt q_extra, PetscInt num_comp_u, PetscInt topo_dim, bool enforce_bc);
27 PetscErrorCode   CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr);
28 CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type);
29 PetscErrorCode   DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field);
30 PetscErrorCode   BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe,
31                                            PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis);
32 PetscErrorCode   CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, BPData bp_data,
33                                      CeedBasis *basis);
34 PetscErrorCode   CreateDistributedDM(RunParams rp, DM *dm);
35 
36 /**
37   @brief Translate array of `PetscInt` to `CeedInt`.
38     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
39     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
40 
41   Not collective across MPI processes.
42 
43   @param[in]      num_entries  Number of array entries
44   @param[in,out]  array_petsc  Array of `PetscInt`
45   @param[out]     array_ceed   Array of `CeedInt`
46 
47   @return An error code: 0 - success, otherwise - failure
48 **/
49 static inline PetscErrorCode IntArrayCeedToPetsc(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) {
50   const CeedInt  int_c = 0;
51   const PetscInt int_p = 0;
52 
53   PetscFunctionBeginUser;
54   if (sizeof(int_c) == sizeof(int_p)) {
55     *array_petsc = (PetscInt *)*array_ceed;
56   } else {
57     *array_petsc = malloc(num_entries * sizeof(PetscInt));
58     for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i];
59     free(*array_ceed);
60   }
61   *array_ceed = NULL;
62   PetscFunctionReturn(PETSC_SUCCESS);
63 }
64 
65 /**
66   @brief Translate array of `PetscInt` to `CeedInt`.
67     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
68     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
69 
70   Not collective across MPI processes.
71 
72   @param[in]      num_entries  Number of array entries
73   @param[in,out]  array_petsc  Array of `PetscInt`
74   @param[out]     array_ceed   Array of `CeedInt`
75 
76   @return An error code: 0 - success, otherwise - failure
77 **/
78 static inline PetscErrorCode IntArrayPetscToCeed(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) {
79   const CeedInt  int_c = 0;
80   const PetscInt int_p = 0;
81 
82   PetscFunctionBeginUser;
83   if (sizeof(int_c) == sizeof(int_p)) {
84     *array_ceed = (CeedInt *)*array_petsc;
85   } else {
86     PetscCall(PetscMalloc1(num_entries, array_ceed));
87     for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i];
88     PetscCall(PetscFree(*array_petsc));
89   }
90   *array_petsc = NULL;
91   PetscFunctionReturn(PETSC_SUCCESS);
92 }
93