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