xref: /libCEED/examples/petsc/src/petscutils.c (revision 652336964d6b29ac5b7f2a27ba4c23e6e61c06af)
1e83e87a5Sjeremylt #include "../include/petscutils.h"
269f5adf1Sjeremylt #include "../include/petscmacros.h"
3e83e87a5Sjeremylt 
4e83e87a5Sjeremylt // -----------------------------------------------------------------------------
5e83e87a5Sjeremylt // Convert PETSc MemType to libCEED MemType
6e83e87a5Sjeremylt // -----------------------------------------------------------------------------
79b072555Sjeremylt CeedMemType MemTypeP2C(PetscMemType mem_type) {
89b072555Sjeremylt   return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST;
9e83e87a5Sjeremylt }
10e83e87a5Sjeremylt 
11e83e87a5Sjeremylt // -----------------------------------------------------------------------------
12e83e87a5Sjeremylt // Utility function taken from petsc/src/dm/impls/plex/examples/tutorials/ex7.c
13e83e87a5Sjeremylt // -----------------------------------------------------------------------------
14e83e87a5Sjeremylt PetscErrorCode ProjectToUnitSphere(DM dm) {
15e83e87a5Sjeremylt   Vec            coordinates;
16e83e87a5Sjeremylt   PetscScalar   *coords;
17e83e87a5Sjeremylt   PetscInt       Nv, v, dim, d;
18e83e87a5Sjeremylt   PetscErrorCode ierr;
19e83e87a5Sjeremylt 
20e83e87a5Sjeremylt   PetscFunctionBeginUser;
21e83e87a5Sjeremylt   ierr = DMGetCoordinatesLocal(dm, &coordinates); CHKERRQ(ierr);
22e83e87a5Sjeremylt   ierr = VecGetLocalSize(coordinates, &Nv); CHKERRQ(ierr);
23e83e87a5Sjeremylt   ierr = VecGetBlockSize(coordinates, &dim); CHKERRQ(ierr);
24e83e87a5Sjeremylt   Nv  /= dim;
25e83e87a5Sjeremylt   ierr = VecGetArray(coordinates, &coords); CHKERRQ(ierr);
26e83e87a5Sjeremylt   for (v = 0; v < Nv; ++v) {
27e83e87a5Sjeremylt     PetscReal r = 0.0;
28e83e87a5Sjeremylt 
29e83e87a5Sjeremylt     for (d = 0; d < dim; ++d) r += PetscSqr(PetscRealPart(coords[v*dim+d]));
30e83e87a5Sjeremylt     r = PetscSqrtReal(r);
31e83e87a5Sjeremylt     for (d = 0; d < dim; ++d) coords[v*dim+d] /= r;
32e83e87a5Sjeremylt   }
33e83e87a5Sjeremylt   ierr = VecRestoreArray(coordinates, &coords); CHKERRQ(ierr);
34e83e87a5Sjeremylt   PetscFunctionReturn(0);
35e83e87a5Sjeremylt };
36e83e87a5Sjeremylt 
37e83e87a5Sjeremylt // -----------------------------------------------------------------------------
382c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation
392c58efb6Sjeremylt // -----------------------------------------------------------------------------
402c58efb6Sjeremylt // Transition from a value of "a" for x=0, to a value of "b" for x=1.  Optionally
412c58efb6Sjeremylt // smooth -- see the commented versions at the end.
422c58efb6Sjeremylt static double step(const double a, const double b, double x) {
432c58efb6Sjeremylt   if (x <= 0) return a;
442c58efb6Sjeremylt   if (x >= 1) return b;
452c58efb6Sjeremylt   return a + (b-a) * (x);
462c58efb6Sjeremylt }
472c58efb6Sjeremylt 
482c58efb6Sjeremylt // 1D transformation at the right boundary
492c58efb6Sjeremylt static double right(const double eps, const double x) {
502c58efb6Sjeremylt   return (x <= 0.5) ? (2-eps) * x : 1 + eps*(x-1);
512c58efb6Sjeremylt }
522c58efb6Sjeremylt 
532c58efb6Sjeremylt // 1D transformation at the left boundary
542c58efb6Sjeremylt static double left(const double eps, const double x) {
552c58efb6Sjeremylt   return 1-right(eps,1-x);
562c58efb6Sjeremylt }
572c58efb6Sjeremylt 
582c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation
592c58efb6Sjeremylt // The eps parameters are in (0, 1]
602c58efb6Sjeremylt // Uniform mesh is recovered for eps=1
619b072555Sjeremylt PetscErrorCode Kershaw(DM dm_orig, PetscScalar eps) {
622c58efb6Sjeremylt   PetscErrorCode ierr;
632c58efb6Sjeremylt   Vec coord;
642c58efb6Sjeremylt   PetscInt ncoord;
652c58efb6Sjeremylt   PetscScalar *c;
662c58efb6Sjeremylt 
672c58efb6Sjeremylt   PetscFunctionBeginUser;
689b072555Sjeremylt   ierr = DMGetCoordinatesLocal(dm_orig, &coord); CHKERRQ(ierr);
692c58efb6Sjeremylt   ierr = VecGetLocalSize(coord, &ncoord); CHKERRQ(ierr);
702c58efb6Sjeremylt   ierr = VecGetArray(coord, &c); CHKERRQ(ierr);
712c58efb6Sjeremylt 
722c58efb6Sjeremylt   for (PetscInt i = 0; i < ncoord; i += 3) {
732c58efb6Sjeremylt     PetscScalar x = c[i], y = c[i+1], z = c[i+2];
742c58efb6Sjeremylt     PetscInt layer = x*6;
752c58efb6Sjeremylt     PetscScalar lambda = (x-layer/6.0)*6;
762c58efb6Sjeremylt     c[i] = x;
772c58efb6Sjeremylt 
782c58efb6Sjeremylt     switch (layer) {
792c58efb6Sjeremylt     case 0:
802c58efb6Sjeremylt       c[i+1] = left(eps, y);
812c58efb6Sjeremylt       c[i+2] = left(eps, z);
822c58efb6Sjeremylt       break;
832c58efb6Sjeremylt     case 1:
842c58efb6Sjeremylt     case 4:
852c58efb6Sjeremylt       c[i+1] = step(left(eps, y), right(eps, y), lambda);
862c58efb6Sjeremylt       c[i+2] = step(left(eps, z), right(eps, z), lambda);
872c58efb6Sjeremylt       break;
882c58efb6Sjeremylt     case 2:
892c58efb6Sjeremylt       c[i+1] = step(right(eps, y), left(eps, y), lambda/2);
902c58efb6Sjeremylt       c[i+2] = step(right(eps, z), left(eps, z), lambda/2);
912c58efb6Sjeremylt       break;
922c58efb6Sjeremylt     case 3:
932c58efb6Sjeremylt       c[i+1] = step(right(eps, y), left(eps, y), (1+lambda)/2);
942c58efb6Sjeremylt       c[i+2] = step(right(eps, z), left(eps, z), (1+lambda)/2);
952c58efb6Sjeremylt       break;
962c58efb6Sjeremylt     default:
972c58efb6Sjeremylt       c[i+1] = right(eps, y);
982c58efb6Sjeremylt       c[i+2] = right(eps, z);
992c58efb6Sjeremylt     }
1002c58efb6Sjeremylt   }
1012c58efb6Sjeremylt   ierr = VecRestoreArray(coord, &c); CHKERRQ(ierr);
1022c58efb6Sjeremylt   PetscFunctionReturn(0);
1032c58efb6Sjeremylt }
1042c58efb6Sjeremylt 
1052c58efb6Sjeremylt // -----------------------------------------------------------------------------
106e83e87a5Sjeremylt // Create BC label
107e83e87a5Sjeremylt // -----------------------------------------------------------------------------
108e83e87a5Sjeremylt static PetscErrorCode CreateBCLabel(DM dm, const char name[]) {
109e83e87a5Sjeremylt   int ierr;
110e83e87a5Sjeremylt   DMLabel label;
111e83e87a5Sjeremylt 
112e83e87a5Sjeremylt   PetscFunctionBeginUser;
113e83e87a5Sjeremylt 
114e83e87a5Sjeremylt   ierr = DMCreateLabel(dm, name); CHKERRQ(ierr);
115e83e87a5Sjeremylt   ierr = DMGetLabel(dm, name, &label); CHKERRQ(ierr);
116e83e87a5Sjeremylt   ierr = DMPlexMarkBoundaryFaces(dm, 1, label); CHKERRQ(ierr);
117e83e87a5Sjeremylt   ierr = DMPlexLabelComplete(dm, label); CHKERRQ(ierr);
118e83e87a5Sjeremylt 
119e83e87a5Sjeremylt   PetscFunctionReturn(0);
120e83e87a5Sjeremylt };
121e83e87a5Sjeremylt 
122e83e87a5Sjeremylt // -----------------------------------------------------------------------------
123e83e87a5Sjeremylt // This function sets up a DM for a given degree
124e83e87a5Sjeremylt // -----------------------------------------------------------------------------
1259b072555Sjeremylt PetscErrorCode SetupDMByDegree(DM dm, PetscInt degree, PetscInt num_comp_u,
1269b072555Sjeremylt                                PetscInt dim, bool enforce_bc, BCFunction bc_func) {
127e83e87a5Sjeremylt   PetscInt ierr, marker_ids[1] = {1};
128e83e87a5Sjeremylt   PetscFE fe;
129*65233696SJeremy L Thompson   MPI_Comm comm;
130e83e87a5Sjeremylt 
131e83e87a5Sjeremylt   PetscFunctionBeginUser;
132e83e87a5Sjeremylt 
133e83e87a5Sjeremylt   // Setup FE
134*65233696SJeremy L Thompson   ierr = PetscObjectGetComm((PetscObject)dm, &comm); CHKERRQ(ierr);
135*65233696SJeremy L Thompson   ierr = PetscFECreateLagrange(comm, dim, num_comp_u, PETSC_FALSE, degree, degree,
136*65233696SJeremy L Thompson                                &fe); CHKERRQ(ierr);
137e83e87a5Sjeremylt   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
138e83e87a5Sjeremylt   ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr);
139e83e87a5Sjeremylt 
140e83e87a5Sjeremylt   // Setup DM
141e83e87a5Sjeremylt   ierr = DMCreateDS(dm); CHKERRQ(ierr);
1429b072555Sjeremylt   if (enforce_bc) {
1439b072555Sjeremylt     PetscBool has_label;
1449b072555Sjeremylt     DMHasLabel(dm, "marker", &has_label);
1459b072555Sjeremylt     if (!has_label) {CreateBCLabel(dm, "marker");}
14669f5adf1Sjeremylt     DMLabel label;
14769f5adf1Sjeremylt     ierr = DMGetLabel(dm, "marker", &label); CHKERRQ(ierr);
14869f5adf1Sjeremylt     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, "marker", 1,
14969f5adf1Sjeremylt                          marker_ids, 0, 0, NULL,
15069f5adf1Sjeremylt                          (void(*)(void))bc_func, NULL, NULL, NULL);
151e83e87a5Sjeremylt     CHKERRQ(ierr);
152e83e87a5Sjeremylt   }
153e83e87a5Sjeremylt   ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
154e83e87a5Sjeremylt   CHKERRQ(ierr);
155e83e87a5Sjeremylt   ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
156e83e87a5Sjeremylt 
157e83e87a5Sjeremylt   PetscFunctionReturn(0);
158e83e87a5Sjeremylt };
159e83e87a5Sjeremylt 
160e83e87a5Sjeremylt // -----------------------------------------------------------------------------
161e83e87a5Sjeremylt // Utility function - essential BC dofs are encoded in closure indices as -(i+1)
162e83e87a5Sjeremylt // -----------------------------------------------------------------------------
163e83e87a5Sjeremylt PetscInt Involute(PetscInt i) {
164e83e87a5Sjeremylt   return i >= 0 ? i : -(i + 1);
165e83e87a5Sjeremylt };
166e83e87a5Sjeremylt 
167e83e87a5Sjeremylt // -----------------------------------------------------------------------------
168e83e87a5Sjeremylt // Get CEED restriction data from DMPlex
169e83e87a5Sjeremylt // -----------------------------------------------------------------------------
170e83e87a5Sjeremylt PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
1719b072555Sjeremylt     CeedInt topo_dim, CeedInt height, DMLabel domain_label, CeedInt value,
1729b072555Sjeremylt     CeedElemRestriction *elem_restr) {
173e83e87a5Sjeremylt   PetscSection section;
1749b072555Sjeremylt   PetscInt p, num_elem, num_dof, *elem_restr_offsets, e_offset, num_fields, dim,
1759b072555Sjeremylt            depth;
1769b072555Sjeremylt   DMLabel depth_label;
1779b072555Sjeremylt   IS depth_is, iter_is;
1789b072555Sjeremylt   Vec U_loc;
1799b072555Sjeremylt   const PetscInt *iter_indices;
180e83e87a5Sjeremylt   PetscErrorCode ierr;
181e83e87a5Sjeremylt 
182e83e87a5Sjeremylt   PetscFunctionBeginUser;
183e83e87a5Sjeremylt 
184e83e87a5Sjeremylt   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
185e83e87a5Sjeremylt   dim -= height;
186e83e87a5Sjeremylt   ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
1879b072555Sjeremylt   ierr = PetscSectionGetNumFields(section, &num_fields); CHKERRQ(ierr);
1889b072555Sjeremylt   PetscInt num_comp[num_fields], field_off[num_fields+1];
1899b072555Sjeremylt   field_off[0] = 0;
1909b072555Sjeremylt   for (PetscInt f = 0; f < num_fields; f++) {
1919b072555Sjeremylt     ierr = PetscSectionGetFieldComponents(section, f, &num_comp[f]); CHKERRQ(ierr);
1929b072555Sjeremylt     field_off[f+1] = field_off[f] + num_comp[f];
193e83e87a5Sjeremylt   }
194e83e87a5Sjeremylt 
195e83e87a5Sjeremylt   ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
1969b072555Sjeremylt   ierr = DMPlexGetDepthLabel(dm, &depth_label); CHKERRQ(ierr);
1979b072555Sjeremylt   ierr = DMLabelGetStratumIS(depth_label, depth - height, &depth_is);
1989b072555Sjeremylt   CHKERRQ(ierr);
1999b072555Sjeremylt   if (domain_label) {
2009b072555Sjeremylt     IS domain_is;
2019b072555Sjeremylt     ierr = DMLabelGetStratumIS(domain_label, value, &domain_is); CHKERRQ(ierr);
2029b072555Sjeremylt     if (domain_is) { // domain_is is non-empty
2039b072555Sjeremylt       ierr = ISIntersect(depth_is, domain_is, &iter_is); CHKERRQ(ierr);
2049b072555Sjeremylt       ierr = ISDestroy(&domain_is); CHKERRQ(ierr);
2059b072555Sjeremylt     } else { // domain_is is NULL (empty)
2069b072555Sjeremylt       iter_is = NULL;
207e83e87a5Sjeremylt     }
2089b072555Sjeremylt     ierr = ISDestroy(&depth_is); CHKERRQ(ierr);
209e83e87a5Sjeremylt   } else {
2109b072555Sjeremylt     iter_is = depth_is;
211e83e87a5Sjeremylt   }
2129b072555Sjeremylt   if (iter_is) {
2139b072555Sjeremylt     ierr = ISGetLocalSize(iter_is, &num_elem); CHKERRQ(ierr);
2149b072555Sjeremylt     ierr = ISGetIndices(iter_is, &iter_indices); CHKERRQ(ierr);
215e83e87a5Sjeremylt   } else {
2169b072555Sjeremylt     num_elem = 0;
2179b072555Sjeremylt     iter_indices = NULL;
218e83e87a5Sjeremylt   }
2199b072555Sjeremylt   ierr = PetscMalloc1(num_elem*PetscPowInt(P, topo_dim), &elem_restr_offsets);
2209b072555Sjeremylt   CHKERRQ(ierr);
2219b072555Sjeremylt   for (p = 0, e_offset = 0; p < num_elem; p++) {
2229b072555Sjeremylt     PetscInt c = iter_indices[p];
2239b072555Sjeremylt     PetscInt num_indices, *indices, num_nodes;
224e83e87a5Sjeremylt     ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE,
2259b072555Sjeremylt                                    &num_indices, &indices, NULL, NULL);
226e83e87a5Sjeremylt     CHKERRQ(ierr);
227e83e87a5Sjeremylt     bool flip = false;
228e83e87a5Sjeremylt     if (height > 0) {
2299b072555Sjeremylt       PetscInt num_cells, num_faces, start = -1;
230e83e87a5Sjeremylt       const PetscInt *orients, *faces, *cells;
231e83e87a5Sjeremylt       ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr);
2329b072555Sjeremylt       ierr = DMPlexGetSupportSize(dm, c, &num_cells); CHKERRQ(ierr);
2339b072555Sjeremylt       if (num_cells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
234e83e87a5Sjeremylt                                      "Expected one cell in support of exterior face, but got %D cells",
2359b072555Sjeremylt                                      num_cells);
236e83e87a5Sjeremylt       ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr);
2379b072555Sjeremylt       ierr = DMPlexGetConeSize(dm, cells[0], &num_faces); CHKERRQ(ierr);
2389b072555Sjeremylt       for (PetscInt i=0; i<num_faces; i++) {if (faces[i] == c) start = i;}
239e83e87a5Sjeremylt       if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT,
240e83e87a5Sjeremylt                                 "Could not find face %D in cone of its support",
241e83e87a5Sjeremylt                                 c);
242e83e87a5Sjeremylt       ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr);
243e83e87a5Sjeremylt       if (orients[start] < 0) flip = true;
244e83e87a5Sjeremylt     }
2459b072555Sjeremylt     if (num_indices % field_off[num_fields]) SETERRQ1(PETSC_COMM_SELF,
246e83e87a5Sjeremylt           PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D",
247e83e87a5Sjeremylt           c);
2489b072555Sjeremylt     num_nodes = num_indices / field_off[num_fields];
2499b072555Sjeremylt     for (PetscInt i = 0; i < num_nodes; i++) {
250e83e87a5Sjeremylt       PetscInt ii = i;
251e83e87a5Sjeremylt       if (flip) {
2529b072555Sjeremylt         if (P == num_nodes) ii = num_nodes - 1 - i;
2539b072555Sjeremylt         else if (P*P == num_nodes) {
254e83e87a5Sjeremylt           PetscInt row = i / P, col = i % P;
255e83e87a5Sjeremylt           ii = row + col * P;
256e83e87a5Sjeremylt         } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP,
257e83e87a5Sjeremylt                           "No support for flipping point with %D nodes != P (%D) or P^2",
2589b072555Sjeremylt                           num_nodes, P);
259e83e87a5Sjeremylt       }
260e83e87a5Sjeremylt       // Check that indices are blocked by node and thus can be coalesced as a single field with
2619b072555Sjeremylt       // field_off[num_fields] = sum(num_comp) components.
2629b072555Sjeremylt       for (PetscInt f = 0; f < num_fields; f++) {
2639b072555Sjeremylt         for (PetscInt j = 0; j < num_comp[f]; j++) {
2649b072555Sjeremylt           if (Involute(indices[field_off[f]*num_nodes + ii*num_comp[f] + j])
2659b072555Sjeremylt               != Involute(indices[ii*num_comp[0]]) + field_off[f] + j)
266e83e87a5Sjeremylt             SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
267e83e87a5Sjeremylt                      "Cell %D closure indices not interlaced for node %D field %D component %D",
268e83e87a5Sjeremylt                      c, ii, f, j);
269e83e87a5Sjeremylt         }
270e83e87a5Sjeremylt       }
271e83e87a5Sjeremylt       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
2729b072555Sjeremylt       PetscInt loc = Involute(indices[ii*num_comp[0]]);
2739b072555Sjeremylt       elem_restr_offsets[e_offset++] = loc;
274e83e87a5Sjeremylt     }
275e83e87a5Sjeremylt     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE,
2769b072555Sjeremylt                                        &num_indices, &indices, NULL, NULL);
277e83e87a5Sjeremylt     CHKERRQ(ierr);
278e83e87a5Sjeremylt   }
2799b072555Sjeremylt   if (e_offset != num_elem*PetscPowInt(P, topo_dim))
280e83e87a5Sjeremylt     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
2819b072555Sjeremylt              "ElemRestriction of size (%D,%D) initialized %D nodes", num_elem,
2829b072555Sjeremylt              PetscPowInt(P, topo_dim),e_offset);
2839b072555Sjeremylt   if (iter_is) {
2849b072555Sjeremylt     ierr = ISRestoreIndices(iter_is, &iter_indices); CHKERRQ(ierr);
285e83e87a5Sjeremylt   }
2869b072555Sjeremylt   ierr = ISDestroy(&iter_is); CHKERRQ(ierr);
287e83e87a5Sjeremylt 
2889b072555Sjeremylt   ierr = DMGetLocalVector(dm, &U_loc); CHKERRQ(ierr);
2899b072555Sjeremylt   ierr = VecGetLocalSize(U_loc, &num_dof); CHKERRQ(ierr);
2909b072555Sjeremylt   ierr = DMRestoreLocalVector(dm, &U_loc); CHKERRQ(ierr);
2919b072555Sjeremylt   CeedElemRestrictionCreate(ceed, num_elem, PetscPowInt(P, topo_dim),
2929b072555Sjeremylt                             field_off[num_fields], 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
2939b072555Sjeremylt                             elem_restr_offsets, elem_restr);
2949b072555Sjeremylt   ierr = PetscFree(elem_restr_offsets); CHKERRQ(ierr);
295e83e87a5Sjeremylt   PetscFunctionReturn(0);
296e83e87a5Sjeremylt };
297e83e87a5Sjeremylt 
298e83e87a5Sjeremylt // -----------------------------------------------------------------------------
299