xref: /libCEED/examples/petsc/src/petscutils.c (revision 751eb813a262a2147c23ac2e55685e7f0d6b8af0)
1e83e87a5Sjeremylt #include "../include/petscutils.h"
2e83e87a5Sjeremylt 
3e83e87a5Sjeremylt // -----------------------------------------------------------------------------
4e83e87a5Sjeremylt // Convert PETSc MemType to libCEED MemType
5e83e87a5Sjeremylt // -----------------------------------------------------------------------------
69b072555Sjeremylt CeedMemType MemTypeP2C(PetscMemType mem_type) {
79b072555Sjeremylt   return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST;
8e83e87a5Sjeremylt }
9e83e87a5Sjeremylt 
10e83e87a5Sjeremylt // -----------------------------------------------------------------------------
112c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation
122c58efb6Sjeremylt // -----------------------------------------------------------------------------
132c58efb6Sjeremylt // Transition from a value of "a" for x=0, to a value of "b" for x=1.  Optionally
142c58efb6Sjeremylt // smooth -- see the commented versions at the end.
152c58efb6Sjeremylt static double step(const double a, const double b, double x) {
162c58efb6Sjeremylt   if (x <= 0) return a;
172c58efb6Sjeremylt   if (x >= 1) return b;
182c58efb6Sjeremylt   return a + (b-a) * (x);
192c58efb6Sjeremylt }
202c58efb6Sjeremylt 
212c58efb6Sjeremylt // 1D transformation at the right boundary
222c58efb6Sjeremylt static double right(const double eps, const double x) {
232c58efb6Sjeremylt   return (x <= 0.5) ? (2-eps) * x : 1 + eps*(x-1);
242c58efb6Sjeremylt }
252c58efb6Sjeremylt 
262c58efb6Sjeremylt // 1D transformation at the left boundary
272c58efb6Sjeremylt static double left(const double eps, const double x) {
282c58efb6Sjeremylt   return 1-right(eps,1-x);
292c58efb6Sjeremylt }
302c58efb6Sjeremylt 
312c58efb6Sjeremylt // Apply 3D Kershaw mesh transformation
322c58efb6Sjeremylt // The eps parameters are in (0, 1]
332c58efb6Sjeremylt // Uniform mesh is recovered for eps=1
349b072555Sjeremylt PetscErrorCode Kershaw(DM dm_orig, PetscScalar eps) {
352c58efb6Sjeremylt   PetscErrorCode ierr;
362c58efb6Sjeremylt   Vec coord;
372c58efb6Sjeremylt   PetscInt ncoord;
382c58efb6Sjeremylt   PetscScalar *c;
392c58efb6Sjeremylt 
402c58efb6Sjeremylt   PetscFunctionBeginUser;
419b072555Sjeremylt   ierr = DMGetCoordinatesLocal(dm_orig, &coord); CHKERRQ(ierr);
422c58efb6Sjeremylt   ierr = VecGetLocalSize(coord, &ncoord); CHKERRQ(ierr);
432c58efb6Sjeremylt   ierr = VecGetArray(coord, &c); CHKERRQ(ierr);
442c58efb6Sjeremylt 
452c58efb6Sjeremylt   for (PetscInt i = 0; i < ncoord; i += 3) {
462c58efb6Sjeremylt     PetscScalar x = c[i], y = c[i+1], z = c[i+2];
472c58efb6Sjeremylt     PetscInt layer = x*6;
482c58efb6Sjeremylt     PetscScalar lambda = (x-layer/6.0)*6;
492c58efb6Sjeremylt     c[i] = x;
502c58efb6Sjeremylt 
512c58efb6Sjeremylt     switch (layer) {
522c58efb6Sjeremylt     case 0:
532c58efb6Sjeremylt       c[i+1] = left(eps, y);
542c58efb6Sjeremylt       c[i+2] = left(eps, z);
552c58efb6Sjeremylt       break;
562c58efb6Sjeremylt     case 1:
572c58efb6Sjeremylt     case 4:
582c58efb6Sjeremylt       c[i+1] = step(left(eps, y), right(eps, y), lambda);
592c58efb6Sjeremylt       c[i+2] = step(left(eps, z), right(eps, z), lambda);
602c58efb6Sjeremylt       break;
612c58efb6Sjeremylt     case 2:
622c58efb6Sjeremylt       c[i+1] = step(right(eps, y), left(eps, y), lambda/2);
632c58efb6Sjeremylt       c[i+2] = step(right(eps, z), left(eps, z), lambda/2);
642c58efb6Sjeremylt       break;
652c58efb6Sjeremylt     case 3:
662c58efb6Sjeremylt       c[i+1] = step(right(eps, y), left(eps, y), (1+lambda)/2);
672c58efb6Sjeremylt       c[i+2] = step(right(eps, z), left(eps, z), (1+lambda)/2);
682c58efb6Sjeremylt       break;
692c58efb6Sjeremylt     default:
702c58efb6Sjeremylt       c[i+1] = right(eps, y);
712c58efb6Sjeremylt       c[i+2] = right(eps, z);
722c58efb6Sjeremylt     }
732c58efb6Sjeremylt   }
742c58efb6Sjeremylt   ierr = VecRestoreArray(coord, &c); CHKERRQ(ierr);
752c58efb6Sjeremylt   PetscFunctionReturn(0);
762c58efb6Sjeremylt }
772c58efb6Sjeremylt 
782c58efb6Sjeremylt // -----------------------------------------------------------------------------
79e83e87a5Sjeremylt // Create BC label
80e83e87a5Sjeremylt // -----------------------------------------------------------------------------
81e83e87a5Sjeremylt static PetscErrorCode CreateBCLabel(DM dm, const char name[]) {
82*751eb813Srezgarshakeri 
83e83e87a5Sjeremylt   DMLabel label;
84e83e87a5Sjeremylt 
85e83e87a5Sjeremylt   PetscFunctionBeginUser;
86e83e87a5Sjeremylt 
87*751eb813Srezgarshakeri   PetscCall(DMCreateLabel(dm, name));
88*751eb813Srezgarshakeri   PetscCall(DMGetLabel(dm, name, &label));
89*751eb813Srezgarshakeri   PetscCall(DMPlexMarkBoundaryFaces(dm, PETSC_DETERMINE, label));
90*751eb813Srezgarshakeri   PetscCall(DMPlexLabelComplete(dm, label));
91e83e87a5Sjeremylt 
92e83e87a5Sjeremylt   PetscFunctionReturn(0);
93e83e87a5Sjeremylt };
94e83e87a5Sjeremylt 
95e83e87a5Sjeremylt // -----------------------------------------------------------------------------
96e83e87a5Sjeremylt // This function sets up a DM for a given degree
97e83e87a5Sjeremylt // -----------------------------------------------------------------------------
98de1229c5Srezgarshakeri PetscErrorCode SetupDMByDegree(DM dm, PetscInt p_degree, PetscInt q_extra,
99de1229c5Srezgarshakeri                                PetscInt num_comp_u,
1009b072555Sjeremylt                                PetscInt dim, bool enforce_bc, BCFunction bc_func) {
101e83e87a5Sjeremylt   PetscInt ierr, marker_ids[1] = {1};
102de1229c5Srezgarshakeri   PetscInt q_degree = p_degree + q_extra;
103e83e87a5Sjeremylt   PetscFE fe;
10465233696SJeremy L Thompson   MPI_Comm comm;
105129d8736Srezgarshakeri   PetscBool      is_simplex = PETSC_TRUE;
106e83e87a5Sjeremylt 
107e83e87a5Sjeremylt   PetscFunctionBeginUser;
108e83e87a5Sjeremylt 
109129d8736Srezgarshakeri   // Check if simplex or tensor-product mesh
110129d8736Srezgarshakeri   ierr = DMPlexIsSimplex(dm, &is_simplex); CHKERRQ(ierr);
111e83e87a5Sjeremylt   // Setup FE
11265233696SJeremy L Thompson   ierr = PetscObjectGetComm((PetscObject)dm, &comm); CHKERRQ(ierr);
113de1229c5Srezgarshakeri   ierr = PetscFECreateLagrange(comm, dim, num_comp_u, is_simplex, p_degree,
114de1229c5Srezgarshakeri                                q_degree, &fe); CHKERRQ(ierr);
115e83e87a5Sjeremylt   ierr = DMAddField(dm, NULL, (PetscObject)fe); CHKERRQ(ierr);
116129d8736Srezgarshakeri   ierr = DMCreateDS(dm); CHKERRQ(ierr);
117129d8736Srezgarshakeri 
118129d8736Srezgarshakeri   {
119129d8736Srezgarshakeri     // create FE field for coordinates
1207ed3e4cdSJeremy L Thompson     PetscFE fe_coords;
1217ed3e4cdSJeremy L Thompson     PetscInt num_comp_coord;
1227ed3e4cdSJeremy L Thompson     ierr = DMGetCoordinateDim(dm, &num_comp_coord); CHKERRQ(ierr);
123de1229c5Srezgarshakeri     ierr = PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, 1, q_degree,
1247ed3e4cdSJeremy L Thompson                                  &fe_coords); CHKERRQ(ierr);
1257ed3e4cdSJeremy L Thompson     ierr = DMProjectCoordinates(dm, fe_coords); CHKERRQ(ierr);
1267ed3e4cdSJeremy L Thompson     ierr = PetscFEDestroy(&fe_coords); CHKERRQ(ierr);
1277ed3e4cdSJeremy L Thompson   }
128de1229c5Srezgarshakeri 
129e83e87a5Sjeremylt   // Setup DM
1309b072555Sjeremylt   if (enforce_bc) {
1319b072555Sjeremylt     PetscBool has_label;
1329b072555Sjeremylt     DMHasLabel(dm, "marker", &has_label);
1339b072555Sjeremylt     if (!has_label) {CreateBCLabel(dm, "marker");}
13469f5adf1Sjeremylt     DMLabel label;
13569f5adf1Sjeremylt     ierr = DMGetLabel(dm, "marker", &label); CHKERRQ(ierr);
136b8962995SJeremy L Thompson     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1,
137b8962995SJeremy L Thompson                          marker_ids, 0, 0, NULL, (void(*)(void))bc_func,
138b8962995SJeremy L Thompson                          NULL, NULL, NULL); CHKERRQ(ierr);
139*751eb813Srezgarshakeri     PetscCall(DMSetOptionsPrefix(dm, "final_"));
140*751eb813Srezgarshakeri     PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
141e83e87a5Sjeremylt   }
142129d8736Srezgarshakeri 
143129d8736Srezgarshakeri   if (!is_simplex) {
144129d8736Srezgarshakeri     DM dm_coord;
145129d8736Srezgarshakeri     ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr);
146de1229c5Srezgarshakeri     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
147de1229c5Srezgarshakeri     CHKERRQ(ierr);
148de1229c5Srezgarshakeri     ierr = DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL);
149de1229c5Srezgarshakeri     CHKERRQ(ierr);
150129d8736Srezgarshakeri   }
151e83e87a5Sjeremylt   ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
152e83e87a5Sjeremylt 
153e83e87a5Sjeremylt   PetscFunctionReturn(0);
154e83e87a5Sjeremylt };
155e83e87a5Sjeremylt 
156e83e87a5Sjeremylt // -----------------------------------------------------------------------------
157e83e87a5Sjeremylt // Get CEED restriction data from DMPlex
158e83e87a5Sjeremylt // -----------------------------------------------------------------------------
1597ed3e4cdSJeremy L Thompson PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height,
1607ed3e4cdSJeremy L Thompson     DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr) {
1617ed3e4cdSJeremy L Thompson   PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets;
162e83e87a5Sjeremylt   PetscErrorCode ierr;
163e83e87a5Sjeremylt 
164e83e87a5Sjeremylt   PetscFunctionBeginUser;
165e83e87a5Sjeremylt 
1667ed3e4cdSJeremy L Thompson   ierr = DMPlexGetLocalOffsets(dm, domain_label, value, height, 0, &num_elem,
1677ed3e4cdSJeremy L Thompson                                &elem_size, &num_comp, &num_dof, &elem_restr_offsets);
1687ed3e4cdSJeremy L Thompson   CHKERRQ(ierr);
169e83e87a5Sjeremylt 
1707ed3e4cdSJeremy L Thompson   CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp,
1717ed3e4cdSJeremy L Thompson                             1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
1729b072555Sjeremylt                             elem_restr_offsets, elem_restr);
1739b072555Sjeremylt   ierr = PetscFree(elem_restr_offsets); CHKERRQ(ierr);
1747ed3e4cdSJeremy L Thompson 
175e83e87a5Sjeremylt   PetscFunctionReturn(0);
176e83e87a5Sjeremylt };
177e83e87a5Sjeremylt 
178e83e87a5Sjeremylt // -----------------------------------------------------------------------------
179129d8736Srezgarshakeri // Utility function - convert from DMPolytopeType to CeedElemTopology
180129d8736Srezgarshakeri // -----------------------------------------------------------------------------
181de1229c5Srezgarshakeri CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) {
182129d8736Srezgarshakeri   switch (cell_type) {
183129d8736Srezgarshakeri   case DM_POLYTOPE_TRIANGLE:      return CEED_TOPOLOGY_TRIANGLE;
184129d8736Srezgarshakeri   case DM_POLYTOPE_QUADRILATERAL: return CEED_TOPOLOGY_QUAD;
185129d8736Srezgarshakeri   case DM_POLYTOPE_TETRAHEDRON:   return CEED_TOPOLOGY_TET;
186129d8736Srezgarshakeri   case DM_POLYTOPE_HEXAHEDRON:    return CEED_TOPOLOGY_HEX;
187129d8736Srezgarshakeri   default:                        return 0;
188129d8736Srezgarshakeri   }
189129d8736Srezgarshakeri }
190129d8736Srezgarshakeri 
191129d8736Srezgarshakeri // -----------------------------------------------------------------------------
192129d8736Srezgarshakeri // Get CEED Basis from DMPlex
193129d8736Srezgarshakeri // -----------------------------------------------------------------------------
194129d8736Srezgarshakeri PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label,
195129d8736Srezgarshakeri                                    CeedInt label_value, CeedInt height,
196129d8736Srezgarshakeri                                    CeedInt dm_field, CeedBasis *basis) {
197129d8736Srezgarshakeri   PetscErrorCode   ierr;
198129d8736Srezgarshakeri   PetscDS          ds;
199129d8736Srezgarshakeri   PetscFE          fe;
200129d8736Srezgarshakeri   PetscQuadrature  quadrature;
201129d8736Srezgarshakeri   PetscBool        is_simplex = PETSC_TRUE;
202129d8736Srezgarshakeri   PetscInt         dim, ds_field = -1, num_comp, P, Q;
203129d8736Srezgarshakeri 
204129d8736Srezgarshakeri   PetscFunctionBeginUser;
205129d8736Srezgarshakeri 
206129d8736Srezgarshakeri   // Get basis information
207129d8736Srezgarshakeri   {
208129d8736Srezgarshakeri     IS             field_is;
209129d8736Srezgarshakeri     const PetscInt *fields;
210129d8736Srezgarshakeri     PetscInt       num_fields;
211129d8736Srezgarshakeri 
212129d8736Srezgarshakeri     ierr = DMGetRegionDS(dm, domain_label, &field_is, &ds); CHKERRQ(ierr);
213129d8736Srezgarshakeri     // Translate dm_field to ds_field
214129d8736Srezgarshakeri     ierr = ISGetIndices(field_is, &fields); CHKERRQ(ierr);
215129d8736Srezgarshakeri     ierr = ISGetSize(field_is, &num_fields); CHKERRQ(ierr);
216129d8736Srezgarshakeri     for (PetscInt i = 0; i < num_fields; i++) {
217129d8736Srezgarshakeri       if (dm_field == fields[i]) {
218129d8736Srezgarshakeri         ds_field = i;
219129d8736Srezgarshakeri         break;
220129d8736Srezgarshakeri       }
221129d8736Srezgarshakeri     }
222129d8736Srezgarshakeri     ierr = ISRestoreIndices(field_is, &fields); CHKERRQ(ierr);
223129d8736Srezgarshakeri   }
224129d8736Srezgarshakeri   if (ds_field == -1) {
225129d8736Srezgarshakeri     // LCOV_EXCL_START
226de1229c5Srezgarshakeri     SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP,
22751ad7d5bSrezgarshakeri             "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);
228129d8736Srezgarshakeri     // LCOV_EXCL_STOP
229129d8736Srezgarshakeri   }
230129d8736Srezgarshakeri 
231129d8736Srezgarshakeri   // Get element information
232129d8736Srezgarshakeri   {
233129d8736Srezgarshakeri     PetscDualSpace dual_space;
234129d8736Srezgarshakeri     PetscInt       num_dual_basis_vectors;
235129d8736Srezgarshakeri 
236129d8736Srezgarshakeri     ierr = PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe);
237129d8736Srezgarshakeri     CHKERRQ(ierr);
238129d8736Srezgarshakeri     ierr = PetscFEGetHeightSubspace(fe, height, &fe); CHKERRQ(ierr);
239129d8736Srezgarshakeri     ierr = PetscFEGetSpatialDimension(fe, &dim); CHKERRQ(ierr);
240129d8736Srezgarshakeri     ierr = PetscFEGetNumComponents(fe, &num_comp); CHKERRQ(ierr);
241129d8736Srezgarshakeri     ierr = PetscFEGetDualSpace(fe, &dual_space); CHKERRQ(ierr);
242129d8736Srezgarshakeri     ierr = PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors);
243129d8736Srezgarshakeri     CHKERRQ(ierr);
244129d8736Srezgarshakeri     P = num_dual_basis_vectors / num_comp;
245129d8736Srezgarshakeri     ierr = PetscFEGetQuadrature(fe, &quadrature); CHKERRQ(ierr);
246129d8736Srezgarshakeri     ierr = PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL);
247129d8736Srezgarshakeri     CHKERRQ(ierr);
248129d8736Srezgarshakeri   }
249129d8736Srezgarshakeri 
250129d8736Srezgarshakeri   // Check if simplex or tensor-product mesh
251129d8736Srezgarshakeri   ierr = DMPlexIsSimplex(dm, &is_simplex); CHKERRQ(ierr);
252129d8736Srezgarshakeri   // Build libCEED basis
253129d8736Srezgarshakeri   if (is_simplex) {
254129d8736Srezgarshakeri     PetscInt          num_derivatives = 1, first_point;
255129d8736Srezgarshakeri     PetscInt          ids[1] = {label_value};
256129d8736Srezgarshakeri     PetscTabulation   basis_tabulation;
257129d8736Srezgarshakeri     const PetscScalar *q_points, *q_weights;
258129d8736Srezgarshakeri     DMLabel           depth_label;
259129d8736Srezgarshakeri     DMPolytopeType    cell_type;
260129d8736Srezgarshakeri     CeedElemTopology  elem_topo;
261129d8736Srezgarshakeri     PetscScalar       *interp, *grad;
262129d8736Srezgarshakeri 
263129d8736Srezgarshakeri     // Use depth label if no domain label present
264129d8736Srezgarshakeri     if (!domain_label) {
265129d8736Srezgarshakeri       PetscInt depth;
266129d8736Srezgarshakeri 
267129d8736Srezgarshakeri       ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
268129d8736Srezgarshakeri       ierr = DMPlexGetDepthLabel(dm, &depth_label); CHKERRQ(ierr);
269129d8736Srezgarshakeri       ids[0] = depth - height;
270129d8736Srezgarshakeri     }
271129d8736Srezgarshakeri     // Get cell interp, grad, and quadrature data
272129d8736Srezgarshakeri     ierr = PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation);
273129d8736Srezgarshakeri     CHKERRQ(ierr);
274129d8736Srezgarshakeri     ierr = PetscQuadratureGetData(quadrature, NULL, NULL, NULL, &q_points,
275129d8736Srezgarshakeri                                   &q_weights); CHKERRQ(ierr);
276129d8736Srezgarshakeri     ierr = DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label,
277129d8736Srezgarshakeri                                   1, ids, height, &first_point, NULL);
278129d8736Srezgarshakeri     CHKERRQ(ierr);
279129d8736Srezgarshakeri     ierr = DMPlexGetCellType(dm, first_point, &cell_type); CHKERRQ(ierr);
280129d8736Srezgarshakeri     elem_topo = ElemTopologyP2C(cell_type);
281129d8736Srezgarshakeri     if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP,
282129d8736Srezgarshakeri                               "DMPlex topology not supported");
283129d8736Srezgarshakeri     // Convert to libCEED orientation
284129d8736Srezgarshakeri     ierr = PetscCalloc(P * Q * sizeof(PetscScalar), &interp); CHKERRQ(ierr);
285129d8736Srezgarshakeri     ierr = PetscCalloc(P * Q * dim * sizeof(PetscScalar), &grad); CHKERRQ(ierr);
286129d8736Srezgarshakeri     const CeedInt c = 0;
287129d8736Srezgarshakeri     for (CeedInt q = 0; q < Q; q++) {
288129d8736Srezgarshakeri       for (CeedInt p = 0; p < P; p++) {
289129d8736Srezgarshakeri         interp[q*P + p] = basis_tabulation->T[0][(q*P + p)*num_comp*num_comp + c];
290129d8736Srezgarshakeri         for (CeedInt d = 0; d < dim; d++) {
291129d8736Srezgarshakeri           grad[(d*Q + q)*P + p] = basis_tabulation->T[1][((q*P + p)*num_comp*num_comp + c)
292129d8736Srezgarshakeri                                   *dim + d];
293129d8736Srezgarshakeri         }
294129d8736Srezgarshakeri       }
295129d8736Srezgarshakeri     }
296129d8736Srezgarshakeri     // Finaly, create libCEED basis
297129d8736Srezgarshakeri     ierr = CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad,
298129d8736Srezgarshakeri                              q_points, q_weights, basis);
299129d8736Srezgarshakeri     CHKERRQ(ierr);
300129d8736Srezgarshakeri     ierr = PetscFree(interp); CHKERRQ(ierr);
301129d8736Srezgarshakeri     ierr = PetscFree(grad); CHKERRQ(ierr);
302129d8736Srezgarshakeri   } else {
303129d8736Srezgarshakeri     CeedInt P_1d = (CeedInt) round(pow(P, 1.0 / dim));
304129d8736Srezgarshakeri     CeedInt Q_1d = (CeedInt) round(pow(Q, 1.0 / dim));
305129d8736Srezgarshakeri 
306129d8736Srezgarshakeri     ierr = CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d,
307129d8736Srezgarshakeri                                            CEED_GAUSS, basis);
308129d8736Srezgarshakeri     CHKERRQ(ierr);
309129d8736Srezgarshakeri   }
310129d8736Srezgarshakeri 
311129d8736Srezgarshakeri   PetscFunctionReturn(0);
312129d8736Srezgarshakeri };
313129d8736Srezgarshakeri 
314129d8736Srezgarshakeri // -----------------------------------------------------------------------------
315de1229c5Srezgarshakeri // Utilities
316de1229c5Srezgarshakeri // -----------------------------------------------------------------------------
317de1229c5Srezgarshakeri 
318de1229c5Srezgarshakeri // Utility function, compute three factors of an integer
319de1229c5Srezgarshakeri static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
320de1229c5Srezgarshakeri   for (PetscInt d=0, size_left=size; d<3; d++) {
321de1229c5Srezgarshakeri     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1./(3 - d)));
322de1229c5Srezgarshakeri     while (try * (size_left / try) != size_left) try++;
323de1229c5Srezgarshakeri     m[reverse ? 2-d : d] = try;
324de1229c5Srezgarshakeri     size_left /= try;
325de1229c5Srezgarshakeri   }
326de1229c5Srezgarshakeri }
327de1229c5Srezgarshakeri 
328de1229c5Srezgarshakeri static int Max3(const PetscInt a[3]) {
329de1229c5Srezgarshakeri   return PetscMax(a[0], PetscMax(a[1], a[2]));
330de1229c5Srezgarshakeri }
331de1229c5Srezgarshakeri 
332de1229c5Srezgarshakeri static int Min3(const PetscInt a[3]) {
333de1229c5Srezgarshakeri   return PetscMin(a[0], PetscMin(a[1], a[2]));
334de1229c5Srezgarshakeri }
335de1229c5Srezgarshakeri 
336de1229c5Srezgarshakeri // -----------------------------------------------------------------------------
337de1229c5Srezgarshakeri // Create distribute dm
338de1229c5Srezgarshakeri // -----------------------------------------------------------------------------
339de1229c5Srezgarshakeri PetscErrorCode CreateDistributedDM(RunParams rp, DM *dm) {
340de1229c5Srezgarshakeri   PetscErrorCode   ierr;
341de1229c5Srezgarshakeri 
342de1229c5Srezgarshakeri   PetscFunctionBeginUser;
343de1229c5Srezgarshakeri   // Setup DM
344de1229c5Srezgarshakeri   if (rp->read_mesh) {
345de1229c5Srezgarshakeri     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, rp->filename, NULL, PETSC_TRUE,
346de1229c5Srezgarshakeri                                 dm);
347de1229c5Srezgarshakeri     CHKERRQ(ierr);
348de1229c5Srezgarshakeri   } else {
349de1229c5Srezgarshakeri     if (rp->user_l_nodes) {
350de1229c5Srezgarshakeri       // Find a nicely composite number of elements no less than global nodes
351de1229c5Srezgarshakeri       PetscMPIInt size;
352de1229c5Srezgarshakeri       ierr = MPI_Comm_size(rp->comm, &size); CHKERRQ(ierr);
353de1229c5Srezgarshakeri       for (PetscInt g_elem =
354de1229c5Srezgarshakeri              PetscMax(1, size * rp->local_nodes / PetscPowInt(rp->degree, rp->dim));
355de1229c5Srezgarshakeri            ;
356de1229c5Srezgarshakeri            g_elem++) {
357de1229c5Srezgarshakeri         Split3(g_elem, rp->mesh_elem, true);
358de1229c5Srezgarshakeri         if (Max3(rp->mesh_elem) / Min3(rp->mesh_elem) <= 2) break;
359de1229c5Srezgarshakeri       }
360de1229c5Srezgarshakeri     }
361de1229c5Srezgarshakeri     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, rp->dim, rp->simplex,
362de1229c5Srezgarshakeri                                rp->mesh_elem,
363de1229c5Srezgarshakeri                                NULL, NULL, NULL, PETSC_TRUE, dm); CHKERRQ(ierr);
364de1229c5Srezgarshakeri   }
365de1229c5Srezgarshakeri 
366de1229c5Srezgarshakeri   ierr = DMSetFromOptions(*dm); CHKERRQ(ierr);
367de1229c5Srezgarshakeri   ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr);
368de1229c5Srezgarshakeri 
369de1229c5Srezgarshakeri   PetscFunctionReturn(0);
370de1229c5Srezgarshakeri }
371de1229c5Srezgarshakeri 
372de1229c5Srezgarshakeri // -----------------------------------------------------------------------------
373