xref: /libCEED/examples/petsc/src/petscutils.c (revision 7ed3e4cde27d28430628eaa24b22da48dc51cc32)
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;
12965233696SJeremy L Thompson   MPI_Comm comm;
130e83e87a5Sjeremylt 
131e83e87a5Sjeremylt   PetscFunctionBeginUser;
132e83e87a5Sjeremylt 
133e83e87a5Sjeremylt   // Setup FE
13465233696SJeremy L Thompson   ierr = PetscObjectGetComm((PetscObject)dm, &comm); CHKERRQ(ierr);
13565233696SJeremy L Thompson   ierr = PetscFECreateLagrange(comm, dim, num_comp_u, PETSC_FALSE, degree, degree,
13665233696SJeremy L Thompson                                &fe); CHKERRQ(ierr);
137e83e87a5Sjeremylt   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
138e83e87a5Sjeremylt   ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr);
139*7ed3e4cdSJeremy L Thompson   {
140*7ed3e4cdSJeremy L Thompson     /* create FE field for coordinates */
141*7ed3e4cdSJeremy L Thompson     PetscFE fe_coords;
142*7ed3e4cdSJeremy L Thompson     PetscInt num_comp_coord;
143*7ed3e4cdSJeremy L Thompson     ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr);
144*7ed3e4cdSJeremy L Thompson     ierr = PetscFECreateLagrange(comm, dim, num_comp_coord, PETSC_FALSE, 1, 1,
145*7ed3e4cdSJeremy L Thompson                                  &fe_coords); CHKERRQ(ierr);
146*7ed3e4cdSJeremy L Thompson     ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr);
147*7ed3e4cdSJeremy L Thompson     ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr);
148*7ed3e4cdSJeremy L Thompson   }
149e83e87a5Sjeremylt 
150e83e87a5Sjeremylt   // Setup DM
151e83e87a5Sjeremylt   ierr = DMCreateDS(dm); CHKERRQ(ierr);
1529b072555Sjeremylt   if (enforce_bc) {
1539b072555Sjeremylt     PetscBool has_label;
1549b072555Sjeremylt     DMHasLabel(dm, "marker", &has_label);
1559b072555Sjeremylt     if (!has_label) {CreateBCLabel(dm, "marker");}
15669f5adf1Sjeremylt     DMLabel label;
15769f5adf1Sjeremylt     ierr = DMGetLabel(dm, "marker", &label); CHKERRQ(ierr);
15869f5adf1Sjeremylt     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, "marker", 1,
15969f5adf1Sjeremylt                          marker_ids, 0, 0, NULL,
16069f5adf1Sjeremylt                          (void(*)(void))bc_func, NULL, NULL, NULL);
161e83e87a5Sjeremylt     CHKERRQ(ierr);
162e83e87a5Sjeremylt   }
163e83e87a5Sjeremylt   ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
164e83e87a5Sjeremylt   CHKERRQ(ierr);
165e83e87a5Sjeremylt   ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
166e83e87a5Sjeremylt 
167e83e87a5Sjeremylt   PetscFunctionReturn(0);
168e83e87a5Sjeremylt };
169e83e87a5Sjeremylt 
170e83e87a5Sjeremylt // -----------------------------------------------------------------------------
171e83e87a5Sjeremylt // Utility function - essential BC dofs are encoded in closure indices as -(i+1)
172e83e87a5Sjeremylt // -----------------------------------------------------------------------------
173e83e87a5Sjeremylt PetscInt Involute(PetscInt i) {
174e83e87a5Sjeremylt   return i >= 0 ? i : -(i + 1);
175e83e87a5Sjeremylt };
176e83e87a5Sjeremylt 
177e83e87a5Sjeremylt // -----------------------------------------------------------------------------
178e83e87a5Sjeremylt // Get CEED restriction data from DMPlex
179e83e87a5Sjeremylt // -----------------------------------------------------------------------------
180*7ed3e4cdSJeremy L Thompson PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height,
181*7ed3e4cdSJeremy L Thompson     DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr) {
182*7ed3e4cdSJeremy L Thompson   PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets;
183e83e87a5Sjeremylt   PetscErrorCode ierr;
184e83e87a5Sjeremylt 
185e83e87a5Sjeremylt   PetscFunctionBeginUser;
186e83e87a5Sjeremylt 
187*7ed3e4cdSJeremy L Thompson   ierr = DMPlexGetLocalOffsets(dm, domain_label, value, height, 0, &num_elem,
188*7ed3e4cdSJeremy L Thompson                                &elem_size, &num_comp, &num_dof, &elem_restr_offsets);
189*7ed3e4cdSJeremy L Thompson   CHKERRQ(ierr);
190e83e87a5Sjeremylt 
191*7ed3e4cdSJeremy L Thompson   CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp,
192*7ed3e4cdSJeremy L Thompson                             1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
1939b072555Sjeremylt                             elem_restr_offsets, elem_restr);
1949b072555Sjeremylt   ierr = PetscFree(elem_restr_offsets); CHKERRQ(ierr);
195*7ed3e4cdSJeremy L Thompson 
196e83e87a5Sjeremylt   PetscFunctionReturn(0);
197e83e87a5Sjeremylt };
198e83e87a5Sjeremylt 
199e83e87a5Sjeremylt // -----------------------------------------------------------------------------
200