xref: /libCEED/examples/petsc/src/petscutils.c (revision 51ad7d5b081d35e6048ae671f50ebaa73ec46d56)
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   int ierr;
83   DMLabel label;
84 
85   PetscFunctionBeginUser;
86 
87   ierr = DMCreateLabel(dm, name); CHKERRQ(ierr);
88   ierr = DMGetLabel(dm, name, &label); CHKERRQ(ierr);
89   ierr = DMPlexMarkBoundaryFaces(dm, 1, label); CHKERRQ(ierr);
90 
91   PetscFunctionReturn(0);
92 };
93 
94 // -----------------------------------------------------------------------------
95 // This function sets up a DM for a given degree
96 // -----------------------------------------------------------------------------
97 PetscErrorCode SetupDMByDegree(DM dm, PetscInt p_degree, PetscInt q_extra,
98                                PetscInt num_comp_u,
99                                PetscInt dim, bool enforce_bc, BCFunction bc_func) {
100   PetscInt ierr, marker_ids[1] = {1};
101   PetscInt q_degree = p_degree + q_extra;
102   PetscFE fe;
103   MPI_Comm comm;
104   PetscBool      is_simplex = PETSC_TRUE;
105 
106   PetscFunctionBeginUser;
107 
108   // Check if simplex or tensor-product mesh
109   ierr = DMPlexIsSimplex(dm, &is_simplex); CHKERRQ(ierr);
110   // Setup FE
111   ierr = PetscObjectGetComm((PetscObject)dm, &comm); CHKERRQ(ierr);
112   ierr = PetscFECreateLagrange(comm, dim, num_comp_u, is_simplex, p_degree,
113                                q_degree, &fe); CHKERRQ(ierr);
114   ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr);
115   ierr = DMCreateDS(dm); CHKERRQ(ierr);
116 
117   {
118     // create FE field for coordinates
119     PetscFE fe_coords;
120     PetscInt num_comp_coord;
121     ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr);
122     ierr = PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, 1, q_degree,
123                                  &fe_coords); CHKERRQ(ierr);
124     ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr);
125     ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr);
126   }
127 
128   // Setup DM
129   if (enforce_bc) {
130     PetscBool has_label;
131     DMHasLabel(dm, "marker", &has_label);
132     if (!has_label) {CreateBCLabel(dm, "marker");}
133     DMLabel label;
134     ierr = DMGetLabel(dm, "marker", &label); CHKERRQ(ierr);
135     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1,
136                          marker_ids, 0, 0, NULL, (void(*)(void))bc_func,
137                          NULL, NULL, NULL); CHKERRQ(ierr);
138   }
139 
140   if (!is_simplex) {
141     DM dm_coord;
142     ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr);
143     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
144     CHKERRQ(ierr);
145     ierr = DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL);
146     CHKERRQ(ierr);
147   }
148   ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
149 
150   PetscFunctionReturn(0);
151 };
152 
153 // -----------------------------------------------------------------------------
154 // Get CEED restriction data from DMPlex
155 // -----------------------------------------------------------------------------
156 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height,
157     DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr) {
158   PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets;
159   PetscErrorCode ierr;
160 
161   PetscFunctionBeginUser;
162 
163   ierr = DMPlexGetLocalOffsets(dm, domain_label, value, height, 0, &num_elem,
164                                &elem_size, &num_comp, &num_dof, &elem_restr_offsets);
165   CHKERRQ(ierr);
166 
167   CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp,
168                             1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
169                             elem_restr_offsets, elem_restr);
170   ierr = PetscFree(elem_restr_offsets); CHKERRQ(ierr);
171 
172   PetscFunctionReturn(0);
173 };
174 
175 // -----------------------------------------------------------------------------
176 // Utility function - convert from DMPolytopeType to CeedElemTopology
177 // -----------------------------------------------------------------------------
178 CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) {
179   switch (cell_type) {
180   case DM_POLYTOPE_TRIANGLE:      return CEED_TOPOLOGY_TRIANGLE;
181   case DM_POLYTOPE_QUADRILATERAL: return CEED_TOPOLOGY_QUAD;
182   case DM_POLYTOPE_TETRAHEDRON:   return CEED_TOPOLOGY_TET;
183   case DM_POLYTOPE_HEXAHEDRON:    return CEED_TOPOLOGY_HEX;
184   default:                        return 0;
185   }
186 }
187 
188 // -----------------------------------------------------------------------------
189 // Get CEED Basis from DMPlex
190 // -----------------------------------------------------------------------------
191 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label,
192                                    CeedInt label_value, CeedInt height,
193                                    CeedInt dm_field, CeedBasis *basis) {
194   PetscErrorCode   ierr;
195   PetscDS          ds;
196   PetscFE          fe;
197   PetscQuadrature  quadrature;
198   PetscBool        is_simplex = PETSC_TRUE;
199   PetscInt         dim, ds_field = -1, num_comp, P, Q;
200 
201   PetscFunctionBeginUser;
202 
203   // Get basis information
204   {
205     IS             field_is;
206     const PetscInt *fields;
207     PetscInt       num_fields;
208 
209     ierr = DMGetRegionDS(dm, domain_label, &field_is, &ds); CHKERRQ(ierr);
210     // Translate dm_field to ds_field
211     ierr = ISGetIndices(field_is, &fields); CHKERRQ(ierr);
212     ierr = ISGetSize(field_is, &num_fields); CHKERRQ(ierr);
213     for (PetscInt i = 0; i < num_fields; i++) {
214       if (dm_field == fields[i]) {
215         ds_field = i;
216         break;
217       }
218     }
219     ierr = ISRestoreIndices(field_is, &fields); CHKERRQ(ierr);
220   }
221   if (ds_field == -1) {
222     // LCOV_EXCL_START
223     SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP,
224             "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);
225     // LCOV_EXCL_STOP
226   }
227 
228   // Get element information
229   {
230     PetscDualSpace dual_space;
231     PetscInt       num_dual_basis_vectors;
232 
233     ierr = PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe);
234     CHKERRQ(ierr);
235     ierr = PetscFEGetHeightSubspace(fe, height, &fe); CHKERRQ(ierr);
236     ierr = PetscFEGetSpatialDimension(fe, &dim); CHKERRQ(ierr);
237     ierr = PetscFEGetNumComponents(fe, &num_comp); CHKERRQ(ierr);
238     ierr = PetscFEGetDualSpace(fe, &dual_space); CHKERRQ(ierr);
239     ierr = PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors);
240     CHKERRQ(ierr);
241     P = num_dual_basis_vectors / num_comp;
242     ierr = PetscFEGetQuadrature(fe, &quadrature); CHKERRQ(ierr);
243     ierr = PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL);
244     CHKERRQ(ierr);
245   }
246 
247   // Check if simplex or tensor-product mesh
248   ierr = DMPlexIsSimplex(dm, &is_simplex); CHKERRQ(ierr);
249   // Build libCEED basis
250   if (is_simplex) {
251     PetscInt          num_derivatives = 1, first_point;
252     PetscInt          ids[1] = {label_value};
253     PetscTabulation   basis_tabulation;
254     const PetscScalar *q_points, *q_weights;
255     DMLabel           depth_label;
256     DMPolytopeType    cell_type;
257     CeedElemTopology  elem_topo;
258     PetscScalar       *interp, *grad;
259 
260     // Use depth label if no domain label present
261     if (!domain_label) {
262       PetscInt depth;
263 
264       ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
265       ierr = DMPlexGetDepthLabel(dm, &depth_label); CHKERRQ(ierr);
266       ids[0] = depth - height;
267     }
268     // Get cell interp, grad, and quadrature data
269     ierr = PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation);
270     CHKERRQ(ierr);
271     ierr = PetscQuadratureGetData(quadrature, NULL, NULL, NULL, &q_points,
272                                   &q_weights); CHKERRQ(ierr);
273     ierr = DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label,
274                                   1, ids, height, &first_point, NULL);
275     CHKERRQ(ierr);
276     ierr = DMPlexGetCellType(dm, first_point, &cell_type); CHKERRQ(ierr);
277     elem_topo = ElemTopologyP2C(cell_type);
278     if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP,
279                               "DMPlex topology not supported");
280     // Convert to libCEED orientation
281     ierr = PetscCalloc(P * Q * sizeof(PetscScalar), &interp); CHKERRQ(ierr);
282     ierr = PetscCalloc(P * Q * dim * sizeof(PetscScalar), &grad); CHKERRQ(ierr);
283     const CeedInt c = 0;
284     for (CeedInt q = 0; q < Q; q++) {
285       for (CeedInt p = 0; p < P; p++) {
286         interp[q*P + p] = basis_tabulation->T[0][(q*P + p)*num_comp*num_comp + c];
287         for (CeedInt d = 0; d < dim; d++) {
288           grad[(d*Q + q)*P + p] = basis_tabulation->T[1][((q*P + p)*num_comp*num_comp + c)
289                                   *dim + d];
290         }
291       }
292     }
293     // Finaly, create libCEED basis
294     ierr = CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad,
295                              q_points, q_weights, basis);
296     CHKERRQ(ierr);
297     ierr = PetscFree(interp); CHKERRQ(ierr);
298     ierr = PetscFree(grad); CHKERRQ(ierr);
299   } else {
300     CeedInt P_1d = (CeedInt) round(pow(P, 1.0 / dim));
301     CeedInt Q_1d = (CeedInt) round(pow(Q, 1.0 / dim));
302 
303     ierr = CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d,
304                                            CEED_GAUSS, basis);
305     CHKERRQ(ierr);
306   }
307 
308   PetscFunctionReturn(0);
309 };
310 
311 // -----------------------------------------------------------------------------
312 // Utilities
313 // -----------------------------------------------------------------------------
314 
315 // Utility function, compute three factors of an integer
316 static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
317   for (PetscInt d=0, size_left=size; d<3; d++) {
318     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1./(3 - d)));
319     while (try * (size_left / try) != size_left) try++;
320     m[reverse ? 2-d : d] = try;
321     size_left /= try;
322   }
323 }
324 
325 static int Max3(const PetscInt a[3]) {
326   return PetscMax(a[0], PetscMax(a[1], a[2]));
327 }
328 
329 static int Min3(const PetscInt a[3]) {
330   return PetscMin(a[0], PetscMin(a[1], a[2]));
331 }
332 
333 // -----------------------------------------------------------------------------
334 // Create distribute dm
335 // -----------------------------------------------------------------------------
336 PetscErrorCode CreateDistributedDM(RunParams rp, DM *dm) {
337   PetscErrorCode   ierr;
338 
339   PetscFunctionBeginUser;
340   // Setup DM
341   if (rp->read_mesh) {
342     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, NULL, PETSC_TRUE,
343                                 dm);
344     CHKERRQ(ierr);
345   } else {
346     if (rp->user_l_nodes) {
347       // Find a nicely composite number of elements no less than global nodes
348       PetscMPIInt size;
349       ierr = MPI_Comm_size(rp->comm, &size); CHKERRQ(ierr);
350       for (PetscInt g_elem =
351              PetscMax(1, size * rp->local_nodes / PetscPowInt(rp->degree, rp->dim));
352            ;
353            g_elem++) {
354         Split3(g_elem, rp->mesh_elem, true);
355         if (Max3(rp->mesh_elem) / Min3(rp->mesh_elem) <= 2) break;
356       }
357     }
358     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, rp->simplex,
359                                rp->mesh_elem,
360                                NULL, NULL, NULL, PETSC_TRUE, dm); CHKERRQ(ierr);
361   }
362 
363   ierr = DMSetFromOptions(*dm); CHKERRQ(ierr);
364   ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr);
365 
366   PetscFunctionReturn(0);
367 }
368 
369 // -----------------------------------------------------------------------------
370