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