xref: /libCEED/examples/petsc/src/petscutils.c (revision b761d2cab99a31c30bb32e57e285fc6533e68118)
1 #include "../include/petscutils.h"
2 
3 // -----------------------------------------------------------------------------
4 // Convert PETSc MemType to libCEED MemType
5 // -----------------------------------------------------------------------------
6 CeedMemType MemTypeP2C(PetscMemType mem_type) {
7   return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST;
8 }
9 
10 // -----------------------------------------------------------------------------
11 // Utility function taken from petsc/src/dm/impls/plex/examples/tutorials/ex7.c
12 // -----------------------------------------------------------------------------
13 PetscErrorCode ProjectToUnitSphere(DM dm) {
14   Vec            coordinates;
15   PetscScalar   *coords;
16   PetscInt       Nv, v, dim, d;
17   PetscErrorCode ierr;
18 
19   PetscFunctionBeginUser;
20   ierr = DMGetCoordinatesLocal(dm, &coordinates); CHKERRQ(ierr);
21   ierr = VecGetLocalSize(coordinates, &Nv); CHKERRQ(ierr);
22   ierr = VecGetBlockSize(coordinates, &dim); CHKERRQ(ierr);
23   Nv  /= dim;
24   ierr = VecGetArray(coordinates, &coords); CHKERRQ(ierr);
25   for (v = 0; v < Nv; ++v) {
26     PetscReal r = 0.0;
27 
28     for (d = 0; d < dim; ++d) r += PetscSqr(PetscRealPart(coords[v*dim+d]));
29     r = PetscSqrtReal(r);
30     for (d = 0; d < dim; ++d) coords[v*dim+d] /= r;
31   }
32   ierr = VecRestoreArray(coordinates, &coords); CHKERRQ(ierr);
33   PetscFunctionReturn(0);
34 };
35 
36 // -----------------------------------------------------------------------------
37 // Apply 3D Kershaw mesh transformation
38 // -----------------------------------------------------------------------------
39 // Transition from a value of "a" for x=0, to a value of "b" for x=1.  Optionally
40 // smooth -- see the commented versions at the end.
41 static double step(const double a, const double b, double x) {
42   if (x <= 0) return a;
43   if (x >= 1) return b;
44   return a + (b-a) * (x);
45 }
46 
47 // 1D transformation at the right boundary
48 static double right(const double eps, const double x) {
49   return (x <= 0.5) ? (2-eps) * x : 1 + eps*(x-1);
50 }
51 
52 // 1D transformation at the left boundary
53 static double left(const double eps, const double x) {
54   return 1-right(eps,1-x);
55 }
56 
57 // Apply 3D Kershaw mesh transformation
58 // The eps parameters are in (0, 1]
59 // Uniform mesh is recovered for eps=1
60 PetscErrorCode Kershaw(DM dm_orig, PetscScalar eps) {
61   PetscErrorCode ierr;
62   Vec coord;
63   PetscInt ncoord;
64   PetscScalar *c;
65 
66   PetscFunctionBeginUser;
67   ierr = DMGetCoordinatesLocal(dm_orig, &coord); CHKERRQ(ierr);
68   ierr = VecGetLocalSize(coord, &ncoord); CHKERRQ(ierr);
69   ierr = VecGetArray(coord, &c); CHKERRQ(ierr);
70 
71   for (PetscInt i = 0; i < ncoord; i += 3) {
72     PetscScalar x = c[i], y = c[i+1], z = c[i+2];
73     PetscInt layer = x*6;
74     PetscScalar lambda = (x-layer/6.0)*6;
75     c[i] = x;
76 
77     switch (layer) {
78     case 0:
79       c[i+1] = left(eps, y);
80       c[i+2] = left(eps, z);
81       break;
82     case 1:
83     case 4:
84       c[i+1] = step(left(eps, y), right(eps, y), lambda);
85       c[i+2] = step(left(eps, z), right(eps, z), lambda);
86       break;
87     case 2:
88       c[i+1] = step(right(eps, y), left(eps, y), lambda/2);
89       c[i+2] = step(right(eps, z), left(eps, z), lambda/2);
90       break;
91     case 3:
92       c[i+1] = step(right(eps, y), left(eps, y), (1+lambda)/2);
93       c[i+2] = step(right(eps, z), left(eps, z), (1+lambda)/2);
94       break;
95     default:
96       c[i+1] = right(eps, y);
97       c[i+2] = right(eps, z);
98     }
99   }
100   ierr = VecRestoreArray(coord, &c); CHKERRQ(ierr);
101   PetscFunctionReturn(0);
102 }
103 
104 // -----------------------------------------------------------------------------
105 // PETSc FE Boilerplate
106 // -----------------------------------------------------------------------------
107 PetscErrorCode PetscFECreateByDegree(DM dm, PetscInt dim, PetscInt Nc,
108                                      PetscBool isSimplex, const char prefix[],
109                                      PetscInt order, PetscFE *fem) {
110   PetscQuadrature q, fq;
111   DM              K;
112   PetscSpace      P;
113   PetscDualSpace  Q;
114   PetscInt        quadPointsPerEdge;
115   PetscBool       tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
116   PetscErrorCode  ierr;
117 
118   PetscFunctionBeginUser;
119   /* Create space */
120   ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) dm), &P); CHKERRQ(ierr);
121   ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix); CHKERRQ(ierr);
122   ierr = PetscSpacePolynomialSetTensor(P, tensor); CHKERRQ(ierr);
123   ierr = PetscSpaceSetFromOptions(P); CHKERRQ(ierr);
124   ierr = PetscSpaceSetNumComponents(P, Nc); CHKERRQ(ierr);
125   ierr = PetscSpaceSetNumVariables(P, dim); CHKERRQ(ierr);
126   ierr = PetscSpaceSetDegree(P, order, order); CHKERRQ(ierr);
127   ierr = PetscSpaceSetUp(P); CHKERRQ(ierr);
128   ierr = PetscSpacePolynomialGetTensor(P, &tensor); CHKERRQ(ierr);
129   /* Create dual space */
130   ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) dm), &Q);
131   CHKERRQ(ierr);
132   ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE); CHKERRQ(ierr);
133   ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix); CHKERRQ(ierr);
134   ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K); CHKERRQ(ierr);
135   ierr = PetscDualSpaceSetDM(Q, K); CHKERRQ(ierr);
136   ierr = DMDestroy(&K); CHKERRQ(ierr);
137   ierr = PetscDualSpaceSetNumComponents(Q, Nc); CHKERRQ(ierr);
138   ierr = PetscDualSpaceSetOrder(Q, order); CHKERRQ(ierr);
139   ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor); CHKERRQ(ierr);
140   ierr = PetscDualSpaceSetFromOptions(Q); CHKERRQ(ierr);
141   ierr = PetscDualSpaceSetUp(Q); CHKERRQ(ierr);
142   /* Create element */
143   ierr = PetscFECreate(PetscObjectComm((PetscObject) dm), fem); CHKERRQ(ierr);
144   ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix); CHKERRQ(ierr);
145   ierr = PetscFESetFromOptions(*fem); CHKERRQ(ierr);
146   ierr = PetscFESetBasisSpace(*fem, P); CHKERRQ(ierr);
147   ierr = PetscFESetDualSpace(*fem, Q); CHKERRQ(ierr);
148   ierr = PetscFESetNumComponents(*fem, Nc); CHKERRQ(ierr);
149   ierr = PetscFESetUp(*fem); CHKERRQ(ierr);
150   ierr = PetscSpaceDestroy(&P); CHKERRQ(ierr);
151   ierr = PetscDualSpaceDestroy(&Q); CHKERRQ(ierr);
152   /* Create quadrature */
153   quadPointsPerEdge = PetscMax(order + 1,1);
154   if (isSimplex) {
155     ierr = PetscDTStroudConicalQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0,
156                                           &q); CHKERRQ(ierr);
157     ierr = PetscDTStroudConicalQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0,
158                                           &fq); CHKERRQ(ierr);
159   } else {
160     ierr = PetscDTGaussTensorQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0,
161                                         &q); CHKERRQ(ierr);
162     ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0,
163                                         &fq); CHKERRQ(ierr);
164   }
165   ierr = PetscFESetQuadrature(*fem, q); CHKERRQ(ierr);
166   ierr = PetscFESetFaceQuadrature(*fem, fq); CHKERRQ(ierr);
167   ierr = PetscQuadratureDestroy(&q); CHKERRQ(ierr);
168   ierr = PetscQuadratureDestroy(&fq); CHKERRQ(ierr);
169 
170   PetscFunctionReturn(0);
171 };
172 
173 // -----------------------------------------------------------------------------
174 // Create BC label
175 // -----------------------------------------------------------------------------
176 static PetscErrorCode CreateBCLabel(DM dm, const char name[]) {
177   int ierr;
178   DMLabel label;
179 
180   PetscFunctionBeginUser;
181 
182   ierr = DMCreateLabel(dm, name); CHKERRQ(ierr);
183   ierr = DMGetLabel(dm, name, &label); CHKERRQ(ierr);
184   ierr = DMPlexMarkBoundaryFaces(dm, 1, label); CHKERRQ(ierr);
185   ierr = DMPlexLabelComplete(dm, label); CHKERRQ(ierr);
186 
187   PetscFunctionReturn(0);
188 };
189 
190 // -----------------------------------------------------------------------------
191 // This function sets up a DM for a given degree
192 // -----------------------------------------------------------------------------
193 PetscErrorCode SetupDMByDegree(DM dm, PetscInt degree, PetscInt num_comp_u,
194                                PetscInt dim, bool enforce_bc, BCFunction bc_func) {
195   PetscInt ierr, marker_ids[1] = {1};
196   PetscFE fe;
197 
198   PetscFunctionBeginUser;
199 
200   // Setup FE
201   ierr = PetscFECreateByDegree(dm, dim, num_comp_u, PETSC_FALSE, NULL, degree,
202                                &fe);
203   CHKERRQ(ierr);
204   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
205   ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr);
206 
207   // Setup DM
208   ierr = DMCreateDS(dm); CHKERRQ(ierr);
209   if (enforce_bc) {
210     PetscBool has_label;
211     DMHasLabel(dm, "marker", &has_label);
212     if (!has_label) {CreateBCLabel(dm, "marker");}
213     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL,
214                          (void(*)(void))bc_func, NULL, 1, marker_ids, NULL);
215     CHKERRQ(ierr);
216   }
217   ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
218   CHKERRQ(ierr);
219   ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
220 
221   PetscFunctionReturn(0);
222 };
223 
224 // -----------------------------------------------------------------------------
225 // Utility function - essential BC dofs are encoded in closure indices as -(i+1)
226 // -----------------------------------------------------------------------------
227 PetscInt Involute(PetscInt i) {
228   return i >= 0 ? i : -(i + 1);
229 };
230 
231 // -----------------------------------------------------------------------------
232 // Get CEED restriction data from DMPlex
233 // -----------------------------------------------------------------------------
234 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
235     CeedInt topo_dim, CeedInt height, DMLabel domain_label, CeedInt value,
236     CeedElemRestriction *elem_restr) {
237   PetscSection section;
238   PetscInt p, num_elem, num_dof, *elem_restr_offsets, e_offset, num_fields, dim,
239            depth;
240   DMLabel depth_label;
241   IS depth_is, iter_is;
242   Vec U_loc;
243   const PetscInt *iter_indices;
244   PetscErrorCode ierr;
245 
246   PetscFunctionBeginUser;
247 
248   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
249   dim -= height;
250   ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
251   ierr = PetscSectionGetNumFields(section, &num_fields); CHKERRQ(ierr);
252   PetscInt num_comp[num_fields], field_off[num_fields+1];
253   field_off[0] = 0;
254   for (PetscInt f = 0; f < num_fields; f++) {
255     ierr = PetscSectionGetFieldComponents(section, f, &num_comp[f]); CHKERRQ(ierr);
256     field_off[f+1] = field_off[f] + num_comp[f];
257   }
258 
259   ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
260   ierr = DMPlexGetDepthLabel(dm, &depth_label); CHKERRQ(ierr);
261   ierr = DMLabelGetStratumIS(depth_label, depth - height, &depth_is);
262   CHKERRQ(ierr);
263   if (domain_label) {
264     IS domain_is;
265     ierr = DMLabelGetStratumIS(domain_label, value, &domain_is); CHKERRQ(ierr);
266     if (domain_is) { // domain_is is non-empty
267       ierr = ISIntersect(depth_is, domain_is, &iter_is); CHKERRQ(ierr);
268       ierr = ISDestroy(&domain_is); CHKERRQ(ierr);
269     } else { // domain_is is NULL (empty)
270       iter_is = NULL;
271     }
272     ierr = ISDestroy(&depth_is); CHKERRQ(ierr);
273   } else {
274     iter_is = depth_is;
275   }
276   if (iter_is) {
277     ierr = ISGetLocalSize(iter_is, &num_elem); CHKERRQ(ierr);
278     ierr = ISGetIndices(iter_is, &iter_indices); CHKERRQ(ierr);
279   } else {
280     num_elem = 0;
281     iter_indices = NULL;
282   }
283   ierr = PetscMalloc1(num_elem*PetscPowInt(P, topo_dim), &elem_restr_offsets);
284   CHKERRQ(ierr);
285   for (p = 0, e_offset = 0; p < num_elem; p++) {
286     PetscInt c = iter_indices[p];
287     PetscInt num_indices, *indices, num_nodes;
288     ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE,
289                                    &num_indices, &indices, NULL, NULL);
290     CHKERRQ(ierr);
291     bool flip = false;
292     if (height > 0) {
293       PetscInt num_cells, num_faces, start = -1;
294       const PetscInt *orients, *faces, *cells;
295       ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr);
296       ierr = DMPlexGetSupportSize(dm, c, &num_cells); CHKERRQ(ierr);
297       if (num_cells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
298                                      "Expected one cell in support of exterior face, but got %D cells",
299                                      num_cells);
300       ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr);
301       ierr = DMPlexGetConeSize(dm, cells[0], &num_faces); CHKERRQ(ierr);
302       for (PetscInt i=0; i<num_faces; i++) {if (faces[i] == c) start = i;}
303       if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT,
304                                 "Could not find face %D in cone of its support",
305                                 c);
306       ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr);
307       if (orients[start] < 0) flip = true;
308     }
309     if (num_indices % field_off[num_fields]) SETERRQ1(PETSC_COMM_SELF,
310           PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D",
311           c);
312     num_nodes = num_indices / field_off[num_fields];
313     for (PetscInt i = 0; i < num_nodes; i++) {
314       PetscInt ii = i;
315       if (flip) {
316         if (P == num_nodes) ii = num_nodes - 1 - i;
317         else if (P*P == num_nodes) {
318           PetscInt row = i / P, col = i % P;
319           ii = row + col * P;
320         } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP,
321                           "No support for flipping point with %D nodes != P (%D) or P^2",
322                           num_nodes, P);
323       }
324       // Check that indices are blocked by node and thus can be coalesced as a single field with
325       // field_off[num_fields] = sum(num_comp) components.
326       for (PetscInt f = 0; f < num_fields; f++) {
327         for (PetscInt j = 0; j < num_comp[f]; j++) {
328           if (Involute(indices[field_off[f]*num_nodes + ii*num_comp[f] + j])
329               != Involute(indices[ii*num_comp[0]]) + field_off[f] + j)
330             SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
331                      "Cell %D closure indices not interlaced for node %D field %D component %D",
332                      c, ii, f, j);
333         }
334       }
335       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
336       PetscInt loc = Involute(indices[ii*num_comp[0]]);
337       elem_restr_offsets[e_offset++] = loc;
338     }
339     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE,
340                                        &num_indices, &indices, NULL, NULL);
341     CHKERRQ(ierr);
342   }
343   if (e_offset != num_elem*PetscPowInt(P, topo_dim))
344     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
345              "ElemRestriction of size (%D,%D) initialized %D nodes", num_elem,
346              PetscPowInt(P, topo_dim),e_offset);
347   if (iter_is) {
348     ierr = ISRestoreIndices(iter_is, &iter_indices); CHKERRQ(ierr);
349   }
350   ierr = ISDestroy(&iter_is); CHKERRQ(ierr);
351 
352   ierr = DMGetLocalVector(dm, &U_loc); CHKERRQ(ierr);
353   ierr = VecGetLocalSize(U_loc, &num_dof); CHKERRQ(ierr);
354   ierr = DMRestoreLocalVector(dm, &U_loc); CHKERRQ(ierr);
355   CeedElemRestrictionCreate(ceed, num_elem, PetscPowInt(P, topo_dim),
356                             field_off[num_fields], 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
357                             elem_restr_offsets, elem_restr);
358   ierr = PetscFree(elem_restr_offsets); CHKERRQ(ierr);
359   PetscFunctionReturn(0);
360 };
361 
362 // -----------------------------------------------------------------------------
363