xref: /libCEED/examples/fluids/include/petsc-ceed-utils.h (revision 5037d55df331588b72bed903a854832c66414480)
1*5037d55dSJames Wright // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2*5037d55dSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3*5037d55dSJames Wright //
4*5037d55dSJames Wright // SPDX-License-Identifier: BSD-2-Clause
5*5037d55dSJames Wright //
6*5037d55dSJames Wright // This file is part of CEED:  http://github.com/ceed
7*5037d55dSJames Wright #pragma once
8*5037d55dSJames Wright 
9*5037d55dSJames Wright #include <ceed.h>
10*5037d55dSJames Wright #include <petscdm.h>
11*5037d55dSJames Wright 
12*5037d55dSJames Wright /**
13*5037d55dSJames Wright   @brief Copy the reference to a `Vec`.
14*5037d55dSJames Wright          Note: If `vec_copy` is non-null, it is assumed to be a valid pointer to a `Vec` and `VecDestroy()` will be called.
15*5037d55dSJames Wright 
16*5037d55dSJames Wright   Collective across MPI processes.
17*5037d55dSJames Wright 
18*5037d55dSJames Wright   @param[in]   vec       `Vec` to reference
19*5037d55dSJames Wright   @param[out]  vec_copy  Copy of reference
20*5037d55dSJames Wright 
21*5037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
22*5037d55dSJames Wright **/
23*5037d55dSJames Wright static inline PetscErrorCode VecReferenceCopy(Vec vec, Vec *vec_copy) {
24*5037d55dSJames Wright   PetscFunctionBeginUser;
25*5037d55dSJames Wright   PetscCall(PetscObjectReference((PetscObject)vec));
26*5037d55dSJames Wright   PetscCall(VecDestroy(vec_copy));
27*5037d55dSJames Wright   *vec_copy = vec;
28*5037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
29*5037d55dSJames Wright }
30*5037d55dSJames Wright 
31*5037d55dSJames Wright /**
32*5037d55dSJames Wright   @brief Copy the reference to a `DM`.
33*5037d55dSJames Wright          Note: If `dm_copy` is non-null, it is assumed to be a valid pointer to a `DM` and `DMDestroy()` will be called.
34*5037d55dSJames Wright 
35*5037d55dSJames Wright   Collective across MPI processes.
36*5037d55dSJames Wright 
37*5037d55dSJames Wright   @param[in]   dm       `DM` to reference
38*5037d55dSJames Wright   @param[out]  dm_copy  Copy of reference
39*5037d55dSJames Wright 
40*5037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
41*5037d55dSJames Wright **/
42*5037d55dSJames Wright static inline PetscErrorCode DMReferenceCopy(DM dm, DM *dm_copy) {
43*5037d55dSJames Wright   PetscFunctionBeginUser;
44*5037d55dSJames Wright   PetscCall(PetscObjectReference((PetscObject)dm));
45*5037d55dSJames Wright   PetscCall(DMDestroy(dm_copy));
46*5037d55dSJames Wright   *dm_copy = dm;
47*5037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
48*5037d55dSJames Wright }
49*5037d55dSJames Wright 
50*5037d55dSJames Wright /**
51*5037d55dSJames Wright   @brief Translate PetscMemType to CeedMemType
52*5037d55dSJames Wright 
53*5037d55dSJames Wright   @param[in]  mem_type  PetscMemType
54*5037d55dSJames Wright 
55*5037d55dSJames Wright   @return Equivalent CeedMemType
56*5037d55dSJames Wright **/
57*5037d55dSJames Wright /// @ingroup RatelInternal
58*5037d55dSJames Wright static inline CeedMemType MemTypePetscToCeed(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; }
59*5037d55dSJames Wright 
60*5037d55dSJames Wright /**
61*5037d55dSJames Wright   @brief Translate array of `PetscInt` to `CeedInt`.
62*5037d55dSJames Wright     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
63*5037d55dSJames Wright     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
64*5037d55dSJames Wright 
65*5037d55dSJames Wright   Not collective across MPI processes.
66*5037d55dSJames Wright 
67*5037d55dSJames Wright   @param[in]      num_entries  Number of array entries
68*5037d55dSJames Wright   @param[in,out]  array_petsc  Array of `PetscInt`
69*5037d55dSJames Wright   @param[out]     array_ceed   Array of `CeedInt`
70*5037d55dSJames Wright 
71*5037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
72*5037d55dSJames Wright **/
73*5037d55dSJames Wright static inline PetscErrorCode IntArrayCeedToPetsc(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) {
74*5037d55dSJames Wright   const CeedInt  int_c = 0;
75*5037d55dSJames Wright   const PetscInt int_p = 0;
76*5037d55dSJames Wright 
77*5037d55dSJames Wright   PetscFunctionBeginUser;
78*5037d55dSJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
79*5037d55dSJames Wright     *array_petsc = (PetscInt *)*array_ceed;
80*5037d55dSJames Wright   } else {
81*5037d55dSJames Wright     *array_petsc = malloc(num_entries * sizeof(PetscInt));
82*5037d55dSJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i];
83*5037d55dSJames Wright     free(*array_ceed);
84*5037d55dSJames Wright   }
85*5037d55dSJames Wright   *array_ceed = NULL;
86*5037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
87*5037d55dSJames Wright }
88*5037d55dSJames Wright 
89*5037d55dSJames Wright /**
90*5037d55dSJames Wright   @brief Translate array of `PetscInt` to `CeedInt`.
91*5037d55dSJames Wright     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
92*5037d55dSJames Wright     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
93*5037d55dSJames Wright 
94*5037d55dSJames Wright   Not collective across MPI processes.
95*5037d55dSJames Wright 
96*5037d55dSJames Wright   @param[in]      num_entries  Number of array entries
97*5037d55dSJames Wright   @param[in,out]  array_petsc  Array of `PetscInt`
98*5037d55dSJames Wright   @param[out]     array_ceed   Array of `CeedInt`
99*5037d55dSJames Wright 
100*5037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
101*5037d55dSJames Wright **/
102*5037d55dSJames Wright static inline PetscErrorCode IntArrayPetscToCeed(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) {
103*5037d55dSJames Wright   const CeedInt  int_c = 0;
104*5037d55dSJames Wright   const PetscInt int_p = 0;
105*5037d55dSJames Wright 
106*5037d55dSJames Wright   PetscFunctionBeginUser;
107*5037d55dSJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
108*5037d55dSJames Wright     *array_ceed = (CeedInt *)*array_petsc;
109*5037d55dSJames Wright   } else {
110*5037d55dSJames Wright     PetscCall(PetscMalloc1(num_entries, array_ceed));
111*5037d55dSJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i];
112*5037d55dSJames Wright     PetscCall(PetscFree(*array_petsc));
113*5037d55dSJames Wright   }
114*5037d55dSJames Wright   *array_petsc = NULL;
115*5037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
116*5037d55dSJames Wright }
117*5037d55dSJames Wright 
118*5037d55dSJames Wright /**
119*5037d55dSJames Wright   @brief Transfer array from PETSc `Vec` to `CeedVector`.
120*5037d55dSJames Wright 
121*5037d55dSJames Wright   Collective across MPI processes.
122*5037d55dSJames Wright 
123*5037d55dSJames Wright   @param[in]   X_petsc   PETSc `Vec`
124*5037d55dSJames Wright   @param[out]  mem_type  PETSc `MemType`
125*5037d55dSJames Wright   @param[out]  x_ceed    `CeedVector`
126*5037d55dSJames Wright 
127*5037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
128*5037d55dSJames Wright **/
129*5037d55dSJames Wright static inline PetscErrorCode VecPetscToCeed(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed) {
130*5037d55dSJames Wright   PetscScalar *x;
131*5037d55dSJames Wright 
132*5037d55dSJames Wright   PetscFunctionBeginUser;
133*5037d55dSJames Wright   PetscCall(VecGetArrayAndMemType(X_petsc, &x, mem_type));
134*5037d55dSJames Wright   PetscCallCeed(CeedVectorReturnCeed(x_ceed), CeedVectorSetArray(x_ceed, MemTypePetscToCeed(*mem_type), CEED_USE_POINTER, x));
135*5037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
136*5037d55dSJames Wright }
137*5037d55dSJames Wright 
138*5037d55dSJames Wright /**
139*5037d55dSJames Wright   @brief Transfer array from `CeedVector` to PETSc `Vec`.
140*5037d55dSJames Wright 
141*5037d55dSJames Wright   Collective across MPI processes.
142*5037d55dSJames Wright 
143*5037d55dSJames Wright   @param[in]   x_ceed    `CeedVector`
144*5037d55dSJames Wright   @param[in]   mem_type  PETSc `MemType`
145*5037d55dSJames Wright   @param[out]  X_petsc   PETSc `Vec`
146*5037d55dSJames Wright 
147*5037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
148*5037d55dSJames Wright **/
149*5037d55dSJames Wright static inline PetscErrorCode VecCeedToPetsc(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc) {
150*5037d55dSJames Wright   PetscScalar *x;
151*5037d55dSJames Wright 
152*5037d55dSJames Wright   PetscFunctionBeginUser;
153*5037d55dSJames Wright   PetscCallCeed(CeedVectorReturnCeed(x_ceed), CeedVectorTakeArray(x_ceed, MemTypePetscToCeed(mem_type), &x));
154*5037d55dSJames Wright   PetscCall(VecRestoreArrayAndMemType(X_petsc, &x));
155*5037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
156*5037d55dSJames Wright }
157*5037d55dSJames Wright 
158*5037d55dSJames Wright /**
159*5037d55dSJames Wright   @brief Transfer read only array from PETSc `Vec` to `CeedVector`.
160*5037d55dSJames Wright 
161*5037d55dSJames Wright   Collective across MPI processes.
162*5037d55dSJames Wright 
163*5037d55dSJames Wright   @param[in]   X_petsc   PETSc `Vec`
164*5037d55dSJames Wright   @param[out]  mem_type  PETSc `MemType`
165*5037d55dSJames Wright   @param[out]  x_ceed    `CeedVector`
166*5037d55dSJames Wright 
167*5037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
168*5037d55dSJames Wright **/
169*5037d55dSJames Wright static inline PetscErrorCode VecReadPetscToCeed(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed) {
170*5037d55dSJames Wright   PetscScalar *x;
171*5037d55dSJames Wright 
172*5037d55dSJames Wright   PetscFunctionBeginUser;
173*5037d55dSJames Wright   PetscCall(VecGetArrayReadAndMemType(X_petsc, (const PetscScalar **)&x, mem_type));
174*5037d55dSJames Wright   PetscCallCeed(CeedVectorReturnCeed(x_ceed), CeedVectorSetArray(x_ceed, MemTypePetscToCeed(*mem_type), CEED_USE_POINTER, x));
175*5037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
176*5037d55dSJames Wright }
177*5037d55dSJames Wright 
178*5037d55dSJames Wright /**
179*5037d55dSJames Wright   @brief Transfer read only array from `CeedVector` to PETSc `Vec`.
180*5037d55dSJames Wright 
181*5037d55dSJames Wright   Collective across MPI processes.
182*5037d55dSJames Wright 
183*5037d55dSJames Wright   @param[in]   x_ceed    `CeedVector`
184*5037d55dSJames Wright   @param[in]   mem_type  PETSc `MemType`
185*5037d55dSJames Wright   @param[out]  X_petsc   PETSc `Vec`
186*5037d55dSJames Wright 
187*5037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
188*5037d55dSJames Wright **/
189*5037d55dSJames Wright static inline PetscErrorCode VecReadCeedToPetsc(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc) {
190*5037d55dSJames Wright   PetscScalar *x;
191*5037d55dSJames Wright 
192*5037d55dSJames Wright   PetscFunctionBeginUser;
193*5037d55dSJames Wright   PetscCallCeed(CeedVectorReturnCeed(x_ceed), CeedVectorTakeArray(x_ceed, MemTypePetscToCeed(mem_type), &x));
194*5037d55dSJames Wright   PetscCall(VecRestoreArrayReadAndMemType(X_petsc, (const PetscScalar **)&x));
195*5037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
196*5037d55dSJames Wright }
197*5037d55dSJames Wright 
198*5037d55dSJames Wright /**
199*5037d55dSJames Wright   @brief Copy PETSc `Vec` data into `CeedVector`
200*5037d55dSJames Wright 
201*5037d55dSJames Wright   @param[in]   X_petsc PETSc `Vec`
202*5037d55dSJames Wright   @param[out]  x_ceed  `CeedVector`
203*5037d55dSJames Wright 
204*5037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
205*5037d55dSJames Wright **/
206*5037d55dSJames Wright static inline PetscErrorCode VecCopyPetscToCeed(Vec X_petsc, CeedVector x_ceed) {
207*5037d55dSJames Wright   PetscScalar *x;
208*5037d55dSJames Wright   PetscMemType mem_type;
209*5037d55dSJames Wright   PetscInt     X_size;
210*5037d55dSJames Wright   CeedSize     x_size;
211*5037d55dSJames Wright   Ceed         ceed;
212*5037d55dSJames Wright 
213*5037d55dSJames Wright   PetscFunctionBeginUser;
214*5037d55dSJames Wright   PetscCall(CeedVectorGetCeed(x_ceed, &ceed));
215*5037d55dSJames Wright   PetscCall(VecGetLocalSize(X_petsc, &X_size));
216*5037d55dSJames Wright   PetscCallCeed(ceed, CeedVectorGetLength(x_ceed, &x_size));
217*5037d55dSJames Wright   PetscCheck(X_size == x_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "X_petsc (%" PetscInt_FMT ") and x_ceed (%" CeedSize_FMT ") must be same size",
218*5037d55dSJames Wright              X_size, x_size);
219*5037d55dSJames Wright 
220*5037d55dSJames Wright   PetscCall(VecGetArrayReadAndMemType(X_petsc, (const PetscScalar **)&x, &mem_type));
221*5037d55dSJames Wright   PetscCallCeed(ceed, CeedVectorSetArray(x_ceed, MemTypePetscToCeed(mem_type), CEED_COPY_VALUES, x));
222*5037d55dSJames Wright   PetscCall(VecRestoreArrayReadAndMemType(X_petsc, (const PetscScalar **)&x));
223*5037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
224*5037d55dSJames Wright }
225*5037d55dSJames Wright 
226*5037d55dSJames Wright /**
227*5037d55dSJames Wright   @brief Return the quadrature size from the eval mode, dimension, and number of components
228*5037d55dSJames Wright 
229*5037d55dSJames Wright   @param[in]  eval_mode       The basis evaluation mode
230*5037d55dSJames Wright   @param[in]  dim             The basis dimension
231*5037d55dSJames Wright   @param[in]  num_components  The basis number of components
232*5037d55dSJames Wright 
233*5037d55dSJames Wright   @return The maximum of the two integers
234*5037d55dSJames Wright 
235*5037d55dSJames Wright **/
236*5037d55dSJames Wright /// @ingroup RatelInternal
237*5037d55dSJames Wright static inline CeedInt GetCeedQuadratureSize(CeedEvalMode eval_mode, CeedInt dim, CeedInt num_components) {
238*5037d55dSJames Wright   switch (eval_mode) {
239*5037d55dSJames Wright     case CEED_EVAL_INTERP:
240*5037d55dSJames Wright       return num_components;
241*5037d55dSJames Wright     case CEED_EVAL_GRAD:
242*5037d55dSJames Wright       return dim * num_components;
243*5037d55dSJames Wright     default:
244*5037d55dSJames Wright       return -1;
245*5037d55dSJames Wright   }
246*5037d55dSJames Wright }
247*5037d55dSJames Wright 
248*5037d55dSJames Wright /**
249*5037d55dSJames Wright   @brief Convert from DMPolytopeType to CeedElemTopology
250*5037d55dSJames Wright 
251*5037d55dSJames Wright   @param[in]  cell_type  DMPolytopeType for the cell
252*5037d55dSJames Wright 
253*5037d55dSJames Wright   @return CeedElemTopology, or 0 if no equivelent CeedElemTopology was found
254*5037d55dSJames Wright **/
255*5037d55dSJames Wright static inline CeedElemTopology PolytopeTypePetscToCeed(DMPolytopeType cell_type) {
256*5037d55dSJames Wright   switch (cell_type) {
257*5037d55dSJames Wright     case DM_POLYTOPE_TRIANGLE:
258*5037d55dSJames Wright       return CEED_TOPOLOGY_TRIANGLE;
259*5037d55dSJames Wright     case DM_POLYTOPE_QUADRILATERAL:
260*5037d55dSJames Wright       return CEED_TOPOLOGY_QUAD;
261*5037d55dSJames Wright     case DM_POLYTOPE_TETRAHEDRON:
262*5037d55dSJames Wright       return CEED_TOPOLOGY_TET;
263*5037d55dSJames Wright     case DM_POLYTOPE_HEXAHEDRON:
264*5037d55dSJames Wright       return CEED_TOPOLOGY_HEX;
265*5037d55dSJames Wright     default:
266*5037d55dSJames Wright       return 0;
267*5037d55dSJames Wright   }
268*5037d55dSJames Wright }
269