xref: /libCEED/examples/petsc/src/petscutils.c (revision f755c37a38f44d41e63f3bc842c03351ea78dc06)
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 // Apply 3D Kershaw mesh transformation
12 // -----------------------------------------------------------------------------
13 // Transition from a value of "a" for x=0, to a value of "b" for x=1.  Optionally
14 // smooth -- see the commented versions at the end.
15 static double step(const double a, const double b, double x) {
16   if (x <= 0) return a;
17   if (x >= 1) return b;
18   return a + (b-a) * (x);
19 }
20 
21 // 1D transformation at the right boundary
22 static double right(const double eps, const double x) {
23   return (x <= 0.5) ? (2-eps) * x : 1 + eps*(x-1);
24 }
25 
26 // 1D transformation at the left boundary
27 static double left(const double eps, const double x) {
28   return 1-right(eps,1-x);
29 }
30 
31 // Apply 3D Kershaw mesh transformation
32 // The eps parameters are in (0, 1]
33 // Uniform mesh is recovered for eps=1
34 PetscErrorCode Kershaw(DM dm_orig, PetscScalar eps) {
35   PetscErrorCode ierr;
36   Vec coord;
37   PetscInt ncoord;
38   PetscScalar *c;
39 
40   PetscFunctionBeginUser;
41   ierr = DMGetCoordinatesLocal(dm_orig, &coord); CHKERRQ(ierr);
42   ierr = VecGetLocalSize(coord, &ncoord); CHKERRQ(ierr);
43   ierr = VecGetArray(coord, &c); CHKERRQ(ierr);
44 
45   for (PetscInt i = 0; i < ncoord; i += 3) {
46     PetscScalar x = c[i], y = c[i+1], z = c[i+2];
47     PetscInt layer = x*6;
48     PetscScalar lambda = (x-layer/6.0)*6;
49     c[i] = x;
50 
51     switch (layer) {
52     case 0:
53       c[i+1] = left(eps, y);
54       c[i+2] = left(eps, z);
55       break;
56     case 1:
57     case 4:
58       c[i+1] = step(left(eps, y), right(eps, y), lambda);
59       c[i+2] = step(left(eps, z), right(eps, z), lambda);
60       break;
61     case 2:
62       c[i+1] = step(right(eps, y), left(eps, y), lambda/2);
63       c[i+2] = step(right(eps, z), left(eps, z), lambda/2);
64       break;
65     case 3:
66       c[i+1] = step(right(eps, y), left(eps, y), (1+lambda)/2);
67       c[i+2] = step(right(eps, z), left(eps, z), (1+lambda)/2);
68       break;
69     default:
70       c[i+1] = right(eps, y);
71       c[i+2] = right(eps, z);
72     }
73   }
74   ierr = VecRestoreArray(coord, &c); CHKERRQ(ierr);
75   PetscFunctionReturn(0);
76 }
77 
78 // -----------------------------------------------------------------------------
79 // Create BC label
80 // -----------------------------------------------------------------------------
81 static PetscErrorCode CreateBCLabel(DM dm, const char name[]) {
82 
83   DMLabel label;
84 
85   PetscFunctionBeginUser;
86 
87   PetscCall(DMCreateLabel(dm, name));
88   PetscCall(DMGetLabel(dm, name, &label));
89   PetscCall(DMPlexMarkBoundaryFaces(dm, PETSC_DETERMINE, label));
90   PetscCall(DMPlexLabelComplete(dm, label));
91 
92   PetscFunctionReturn(0);
93 };
94 
95 // -----------------------------------------------------------------------------
96 // This function sets up a DM for a given degree
97 // -----------------------------------------------------------------------------
98 PetscErrorCode SetupDMByDegree(DM dm, PetscInt p_degree, PetscInt q_extra,
99                                PetscInt num_comp_u,
100                                PetscInt dim, bool enforce_bc, BCFunction bc_func) {
101   PetscInt ierr, marker_ids[1] = {1};
102   PetscInt q_degree = p_degree + q_extra;
103   PetscFE fe;
104   MPI_Comm comm;
105   PetscBool      is_simplex = PETSC_TRUE;
106 
107   PetscFunctionBeginUser;
108 
109   // Check if simplex or tensor-product mesh
110   ierr = DMPlexIsSimplex(dm, &is_simplex); CHKERRQ(ierr);
111   // Setup FE
112   ierr = PetscObjectGetComm((PetscObject)dm, &comm); CHKERRQ(ierr);
113   ierr = PetscFECreateLagrange(comm, dim, num_comp_u, is_simplex, p_degree,
114                                q_degree, &fe); CHKERRQ(ierr);
115   ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr);
116   ierr = DMCreateDS(dm); CHKERRQ(ierr);
117 
118   {
119     // create FE field for coordinates
120     PetscFE fe_coords;
121     PetscInt num_comp_coord;
122     ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr);
123     ierr = PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, 1, q_degree,
124                                  &fe_coords); CHKERRQ(ierr);
125     ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr);
126     ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr);
127   }
128 
129   // Setup DM
130   if (enforce_bc) {
131     PetscBool has_label;
132     DMHasLabel(dm, "marker", &has_label);
133     if (!has_label) {CreateBCLabel(dm, "marker");}
134     DMLabel label;
135     ierr = DMGetLabel(dm, "marker", &label); CHKERRQ(ierr);
136     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1,
137                          marker_ids, 0, 0, NULL, (void(*)(void))bc_func,
138                          NULL, NULL, NULL); CHKERRQ(ierr);
139     PetscCall(DMSetOptionsPrefix(dm, "final_"));
140     PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
141   }
142 
143   if (!is_simplex) {
144     DM dm_coord;
145     ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr);
146     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
147     CHKERRQ(ierr);
148     ierr = DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL);
149     CHKERRQ(ierr);
150   }
151   ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
152 
153   PetscFunctionReturn(0);
154 };
155 
156 // -----------------------------------------------------------------------------
157 // Get CEED restriction data from DMPlex
158 // -----------------------------------------------------------------------------
159 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height,
160     DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr) {
161   PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets;
162   PetscErrorCode ierr;
163 
164   PetscFunctionBeginUser;
165 
166   ierr = DMPlexGetLocalOffsets(dm, domain_label, value, height, 0, &num_elem,
167                                &elem_size, &num_comp, &num_dof, &elem_restr_offsets);
168   CHKERRQ(ierr);
169 
170   CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp,
171                             1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
172                             elem_restr_offsets, elem_restr);
173   ierr = PetscFree(elem_restr_offsets); CHKERRQ(ierr);
174 
175   PetscFunctionReturn(0);
176 };
177 
178 // -----------------------------------------------------------------------------
179 // Utility function - convert from DMPolytopeType to CeedElemTopology
180 // -----------------------------------------------------------------------------
181 CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) {
182   switch (cell_type) {
183   case DM_POLYTOPE_TRIANGLE:      return CEED_TOPOLOGY_TRIANGLE;
184   case DM_POLYTOPE_QUADRILATERAL: return CEED_TOPOLOGY_QUAD;
185   case DM_POLYTOPE_TETRAHEDRON:   return CEED_TOPOLOGY_TET;
186   case DM_POLYTOPE_HEXAHEDRON:    return CEED_TOPOLOGY_HEX;
187   default:                        return 0;
188   }
189 }
190 
191 // -----------------------------------------------------------------------------
192 // Convert DM field to DS field
193 // -----------------------------------------------------------------------------
194 PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field,
195                                 PetscInt *ds_field) {
196   PetscDS         ds;
197   IS              field_is;
198   const PetscInt *fields;
199   PetscInt        num_fields;
200 
201   PetscFunctionBeginUser;
202 
203   // Translate dm_field to ds_field
204   PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds));
205   PetscCall(ISGetIndices(field_is, &fields));
206   PetscCall(ISGetSize(field_is, &num_fields));
207   for (PetscInt i = 0; i < num_fields; i++) {
208     if (dm_field == fields[i]) {
209       *ds_field = i;
210       break;
211     }
212   }
213   PetscCall(ISRestoreIndices(field_is, &fields));
214 
215   if (*ds_field == -1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP,
216                                  "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);
217 
218   PetscFunctionReturn(0);
219 }
220 
221 // -----------------------------------------------------------------------------
222 // Create libCEED Basis from PetscTabulation
223 // -----------------------------------------------------------------------------
224 PetscErrorCode BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label,
225     PetscInt label_value, PetscInt height, PetscInt face,
226     PetscFE fe, PetscTabulation basis_tabulation, PetscQuadrature quadrature,
227     CeedBasis *basis) {
228 
229   PetscInt           first_point;
230   PetscInt           ids[1] = {label_value};
231   DMLabel            depth_label;
232   DMPolytopeType     cell_type;
233   CeedElemTopology   elem_topo;
234   PetscScalar       *q_points, *interp, *grad;
235   const PetscScalar *q_weights;
236   PetscDualSpace     dual_space;
237   PetscInt           num_dual_basis_vectors;
238   PetscInt           dim, num_comp, P, Q;
239 
240   PetscFunctionBeginUser;
241 
242   // General basis information
243   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
244   PetscCall(PetscFEGetNumComponents(fe, &num_comp));
245   PetscCall(PetscFEGetDualSpace(fe, &dual_space));
246   PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
247   P = num_dual_basis_vectors / num_comp;
248 
249   // Use depth label if no domain label present
250   if (!domain_label) {
251     PetscInt depth;
252 
253     PetscCall(DMPlexGetDepth(dm, &depth));
254     PetscCall(DMPlexGetDepthLabel(dm, &depth_label));
255     ids[0] = depth - height;
256   }
257 
258   // Get cell interp, grad, and quadrature data
259   PetscCall(DMGetFirstLabeledPoint(dm, dm,
260                                    domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL));
261   PetscCall(DMPlexGetCellType(dm, first_point, &cell_type));
262   elem_topo = ElemTopologyP2C(cell_type);
263   if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP,
264                             "DMPlex topology not supported");
265   {
266     size_t             q_points_size;
267     const PetscScalar *q_points_petsc;
268     PetscInt           q_dim;
269 
270     PetscCall(PetscQuadratureGetData(quadrature, &q_dim, NULL, &Q, &q_points_petsc,
271                                      &q_weights));
272     q_points_size = Q * dim * sizeof(CeedScalar);
273     PetscCall(PetscCalloc(q_points_size, &q_points));
274     for (PetscInt q = 0; q < Q; q++) {
275       for (PetscInt d = 0; d < q_dim;
276            d++) q_points[q * dim + d] = q_points_petsc[q * q_dim + d];
277     }
278   }
279 
280   // Convert to libCEED orientation
281   {
282     PetscBool       is_simplex  = PETSC_FALSE;
283     IS              permutation = NULL;
284     const PetscInt *permutation_indices;
285 
286     PetscCall(DMPlexIsSimplex(dm, &is_simplex));
287     if (!is_simplex) {
288       PetscSection section;
289 
290       // -- Get permutation
291       PetscCall(DMGetLocalSection(dm, &section));
292       PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, dim,
293                 num_comp * P, &permutation));
294       PetscCall(ISGetIndices(permutation, &permutation_indices));
295     }
296 
297     // -- Copy interp, grad matrices
298     PetscCall(PetscCalloc(P * Q * sizeof(CeedScalar), &interp));
299     PetscCall(PetscCalloc(P * Q * dim * sizeof(CeedScalar), &grad));
300     const CeedInt c = 0;
301     for (CeedInt q = 0; q < Q; q++) {
302       for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) {
303         CeedInt p_petsc = is_simplex ? (p_ceed * num_comp) : permutation_indices[p_ceed
304                           * num_comp];
305 
306         interp[q * P + p_ceed] = basis_tabulation->T[0][((face * Q + q) * P * num_comp +
307                                  p_petsc) * num_comp + c];
308         for (CeedInt d = 0; d < dim; d++) {
309           grad[(d * Q + q) * P + p_ceed] = basis_tabulation->T[1][(((
310                                              face * Q + q) * P * num_comp + p_petsc) * num_comp + c) * dim + d];
311         }
312       }
313     }
314 
315     // -- Cleanup
316     if (permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices));
317     PetscCall(ISDestroy(&permutation));
318   }
319 
320   // Finally, create libCEED basis
321   CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points,
322                     q_weights, basis);
323   PetscCall(PetscFree(q_points));
324   PetscCall(PetscFree(interp));
325   PetscCall(PetscFree(grad));
326 
327   PetscFunctionReturn(0);
328 }
329 
330 // -----------------------------------------------------------------------------
331 // Get CEED Basis from DMPlex
332 // -----------------------------------------------------------------------------
333 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label,
334                                    CeedInt label_value, CeedInt height,
335                                    CeedInt dm_field, CeedBasis *basis) {
336   PetscDS         ds;
337   PetscFE         fe;
338   PetscQuadrature quadrature;
339   PetscBool       is_simplex = PETSC_TRUE;
340   PetscInt        ds_field   = -1;
341 
342   PetscFunctionBeginUser;
343 
344   // Get element information
345   PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds));
346   PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field));
347   PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
348   PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
349   PetscCall(PetscFEGetQuadrature(fe, &quadrature));
350 
351   // Check if simplex or tensor-product mesh
352   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
353 
354   // Build libCEED basis
355   if (is_simplex) {
356     PetscTabulation basis_tabulation;
357     PetscInt        num_derivatives = 1, face = 0;
358 
359     PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation));
360     PetscCall(BasisCreateFromTabulation(ceed, dm, domain_label, label_value, height,
361                                         face, fe, basis_tabulation, quadrature, basis));
362   } else {
363     PetscDualSpace dual_space;
364     PetscInt       num_dual_basis_vectors;
365     PetscInt       dim, num_comp, P, Q;
366 
367     PetscCall(PetscFEGetSpatialDimension(fe, &dim));
368     PetscCall(PetscFEGetNumComponents(fe, &num_comp));
369     PetscCall(PetscFEGetDualSpace(fe, &dual_space));
370     PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
371     P = num_dual_basis_vectors / num_comp;
372     PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL));
373 
374     CeedInt P_1d = (CeedInt)round(pow(P, 1.0 / dim));
375     CeedInt Q_1d = (CeedInt)round(pow(Q, 1.0 / dim));
376 
377     CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, CEED_GAUSS,
378                                     basis);
379   }
380 
381   PetscFunctionReturn(0);
382 }
383 
384 // -----------------------------------------------------------------------------
385 // Utilities
386 // -----------------------------------------------------------------------------
387 
388 // Utility function, compute three factors of an integer
389 static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
390   for (PetscInt d=0, size_left=size; d<3; d++) {
391     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1./(3 - d)));
392     while (try * (size_left / try) != size_left) try++;
393     m[reverse ? 2-d : d] = try;
394     size_left /= try;
395   }
396 }
397 
398 static int Max3(const PetscInt a[3]) {
399   return PetscMax(a[0], PetscMax(a[1], a[2]));
400 }
401 
402 static int Min3(const PetscInt a[3]) {
403   return PetscMin(a[0], PetscMin(a[1], a[2]));
404 }
405 
406 // -----------------------------------------------------------------------------
407 // Create distribute dm
408 // -----------------------------------------------------------------------------
409 PetscErrorCode CreateDistributedDM(RunParams rp, DM *dm) {
410   PetscErrorCode   ierr;
411 
412   PetscFunctionBeginUser;
413   // Setup DM
414   if (rp->read_mesh) {
415     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, NULL, PETSC_TRUE,
416                                 dm);
417     CHKERRQ(ierr);
418   } else {
419     if (rp->user_l_nodes) {
420       // Find a nicely composite number of elements no less than global nodes
421       PetscMPIInt size;
422       ierr = MPI_Comm_size(rp->comm, &size); CHKERRQ(ierr);
423       for (PetscInt g_elem =
424              PetscMax(1, size * rp->local_nodes / PetscPowInt(rp->degree, rp->dim));
425            ;
426            g_elem++) {
427         Split3(g_elem, rp->mesh_elem, true);
428         if (Max3(rp->mesh_elem) / Min3(rp->mesh_elem) <= 2) break;
429       }
430     }
431     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, rp->simplex,
432                                rp->mesh_elem,
433                                NULL, NULL, NULL, PETSC_TRUE, dm); CHKERRQ(ierr);
434   }
435 
436   ierr = DMSetFromOptions(*dm); CHKERRQ(ierr);
437   ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr);
438 
439   PetscFunctionReturn(0);
440 }
441 
442 // -----------------------------------------------------------------------------
443