1 #include <petsc/private/dmmbimpl.h> 2 #include <petscdmmoab.h> 3 #include <MBTagConventions.hpp> 4 #include <moab/NestedRefine.hpp> 5 6 // A helper function to convert Real vector to Scalar vector (required by MatSetValues) 7 static inline std::vector<PetscScalar> VecReal_to_VecScalar(const std::vector<PetscReal> &v) 8 { 9 std::vector<PetscScalar> res(v.size()); 10 for (size_t i = 0; i < res.size(); i++) res[i] = v[i]; 11 return res; 12 } 13 14 /*@C 15 DMMoabGenerateHierarchy - Generate a multi-level uniform refinement hierarchy 16 by succesively refining a coarse mesh, already defined in the `DM` object 17 provided by the user. 18 19 Collective 20 21 Input Parameter: 22 . dm - The `DMMOAB` object 23 24 Output Parameters: 25 + nlevels - The number of levels of refinement needed to generate the hierarchy 26 - ldegrees - The degree of refinement at each level in the hierarchy 27 28 Level: beginner 29 30 .seealso: `DMMoabCreate()` 31 @*/ 32 PetscErrorCode DMMoabGenerateHierarchy(DM dm, PetscInt nlevels, PetscInt *ldegrees) 33 { 34 DM_Moab *dmmoab; 35 moab::ErrorCode merr; 36 PetscInt *pdegrees, ilevel; 37 std::vector<moab::EntityHandle> hsets; 38 39 PetscFunctionBegin; 40 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 41 dmmoab = (DM_Moab *)(dm)->data; 42 43 if (!ldegrees) { 44 PetscCall(PetscMalloc1(nlevels, &pdegrees)); 45 for (ilevel = 0; ilevel < nlevels; ilevel++) pdegrees[ilevel] = 2; /* default = Degree 2 refinement */ 46 } else pdegrees = ldegrees; 47 48 /* initialize set level refinement data for hierarchy */ 49 dmmoab->nhlevels = nlevels; 50 51 /* Instantiate the nested refinement class */ 52 #ifdef MOAB_HAVE_MPI 53 dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core *>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset); 54 #else 55 dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core *>(dmmoab->mbiface), NULL, dmmoab->fileset); 56 #endif 57 58 PetscCall(PetscMalloc1(nlevels + 1, &dmmoab->hsets)); 59 60 /* generate the mesh hierarchy */ 61 merr = dmmoab->hierarchy->generate_mesh_hierarchy(nlevels, pdegrees, hsets, false); 62 MBERRNM(merr); 63 64 #ifdef MOAB_HAVE_MPI 65 if (dmmoab->pcomm->size() > 1) { 66 merr = dmmoab->hierarchy->exchange_ghosts(hsets, dmmoab->nghostrings); 67 MBERRNM(merr); 68 } 69 #endif 70 71 /* copy the mesh sets for nested refinement hierarchy */ 72 dmmoab->hsets[0] = hsets[0]; 73 for (ilevel = 1; ilevel <= nlevels; ilevel++) { 74 dmmoab->hsets[ilevel] = hsets[ilevel]; 75 76 #ifdef MOAB_HAVE_MPI 77 merr = dmmoab->pcomm->assign_global_ids(hsets[ilevel], dmmoab->dim, 0, false, true, false); 78 MBERRNM(merr); 79 #endif 80 81 /* Update material and other geometric tags from parent to child sets */ 82 merr = dmmoab->hierarchy->update_special_tags(ilevel, hsets[ilevel]); 83 MBERRNM(merr); 84 } 85 86 hsets.clear(); 87 if (!ldegrees) PetscCall(PetscFree(pdegrees)); 88 PetscFunctionReturn(PETSC_SUCCESS); 89 } 90 91 // PetscClangLinter pragma ignore: -fdoc-* 92 /* 93 DMRefineHierarchy_Moab - Generate a multi-level `DM` hierarchy 94 by succesively refining a coarse mesh. 95 96 Collective 97 98 Input Parameter: 99 . dm - The `DMMOAB` object 100 101 Output Parameters: 102 + nlevels - The number of levels of refinement needed to generate the hierarchy 103 - dmf - The DM objects after successive refinement of the hierarchy 104 105 Level: beginner 106 */ 107 PETSC_EXTERN PetscErrorCode DMRefineHierarchy_Moab(DM dm, PetscInt nlevels, DM dmf[]) 108 { 109 PetscInt i; 110 111 PetscFunctionBegin; 112 113 PetscCall(DMRefine(dm, PetscObjectComm((PetscObject)dm), &dmf[0])); 114 for (i = 1; i < nlevels; i++) PetscCall(DMRefine(dmf[i - 1], PetscObjectComm((PetscObject)dm), &dmf[i])); 115 PetscFunctionReturn(PETSC_SUCCESS); 116 } 117 118 // PetscClangLinter pragma ignore: -fdoc-* 119 /* 120 DMCoarsenHierarchy_Moab - Generate a multi-level `DM` hierarchy 121 by succesively coarsening a refined mesh. 122 123 Collective 124 125 Input Parameter: 126 . dm - The `DMMOAB` object 127 128 Output Parameters: 129 + nlevels - The number of levels of refinement needed to generate the hierarchy 130 - dmc - The `DM` objects after successive coarsening of the hierarchy 131 132 Level: beginner 133 */ 134 PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy_Moab(DM dm, PetscInt nlevels, DM dmc[]) 135 { 136 PetscInt i; 137 138 PetscFunctionBegin; 139 140 PetscCall(DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &dmc[0])); 141 for (i = 1; i < nlevels; i++) PetscCall(DMCoarsen(dmc[i - 1], PetscObjectComm((PetscObject)dm), &dmc[i])); 142 PetscFunctionReturn(PETSC_SUCCESS); 143 } 144 145 PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscBool); 146 147 // PetscClangLinter pragma ignore: -fdoc-* 148 /* 149 DMCreateInterpolation_Moab - Generate the interpolation operators to transform 150 operators (matrices, vectors) from parent level to child level as defined by 151 the `DM` inputs provided by the user. 152 153 Collective 154 155 Input Parameters: 156 + dmp - The `DMMOAB` object 157 - dmc - the second, finer `DMMOAB` object 158 159 Output Parameters: 160 + interpl - The interpolation operator for transferring data between the levels 161 - vec - The scaling vector (optional) 162 163 Level: developer 164 */ 165 PETSC_EXTERN PetscErrorCode DMCreateInterpolation_Moab(DM dmp, DM dmc, Mat *interpl, Vec *vec) 166 { 167 DM_Moab *dmbp, *dmbc; 168 moab::ErrorCode merr; 169 PetscInt dim; 170 PetscReal factor; 171 PetscInt innz, *nnz, ionz, *onz; 172 PetscInt nlsizp, nlsizc, nlghsizp, ngsizp, ngsizc; 173 const PetscBool use_consistent_bases = PETSC_TRUE; 174 175 PetscFunctionBegin; 176 PetscValidHeaderSpecific(dmp, DM_CLASSID, 1); 177 PetscValidHeaderSpecific(dmc, DM_CLASSID, 2); 178 dmbp = (DM_Moab *)(dmp)->data; 179 dmbc = (DM_Moab *)(dmc)->data; 180 nlsizp = dmbp->nloc; // *dmb1->numFields; 181 nlsizc = dmbc->nloc; // *dmb2->numFields; 182 ngsizp = dmbp->n; // *dmb1->numFields; 183 ngsizc = dmbc->n; // *dmb2->numFields; 184 nlghsizp = (dmbp->nloc + dmbp->nghost); // *dmb1->numFields; 185 186 // Columns = Parent DoFs ; Rows = Child DoFs 187 // Interpolation matrix: \sum_{i=1}^P Owned(Child) * (Owned(Parent) + Ghosted(Parent)) 188 // Size: nlsizc * nlghsizp 189 PetscCall(PetscInfo(NULL, "Creating interpolation matrix %" PetscInt_FMT " X %" PetscInt_FMT " to apply transformation between levels %" PetscInt_FMT " -> %" PetscInt_FMT ".\n", ngsizc, nlghsizp, dmbp->hlevel, dmbc->hlevel)); 190 191 PetscCall(DMGetDimension(dmp, &dim)); 192 193 /* allocate the nnz, onz arrays based on block size and local nodes */ 194 PetscCall(PetscCalloc2(nlsizc, &nnz, nlsizc, &onz)); 195 196 /* Loop through the local elements and compute the relation between the current parent and the refined_level. */ 197 for (moab::Range::iterator iter = dmbc->vowned->begin(); iter != dmbc->vowned->end(); iter++) { 198 const moab::EntityHandle vhandle = *iter; 199 /* define local variables */ 200 moab::EntityHandle parent; 201 std::vector<moab::EntityHandle> adjs; 202 moab::Range found; 203 204 /* store the vertex DoF number */ 205 const int ldof = dmbc->lidmap[vhandle - dmbc->seqstart]; 206 207 /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects 208 to the current vertex. We can then decipher if a vertex is ghosted or not and compute the 209 non-zero pattern accordingly. */ 210 merr = dmbc->hierarchy->get_adjacencies(vhandle, dmbc->dim, adjs); 211 MBERRNM(merr); 212 213 /* loop over vertices and update the number of connectivity */ 214 for (unsigned jter = 0; jter < adjs.size(); jter++) { 215 const moab::EntityHandle jhandle = adjs[jter]; 216 217 /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 218 merr = dmbc->hierarchy->child_to_parent(jhandle, dmbc->hlevel, dmbp->hlevel, &parent); 219 MBERRNM(merr); 220 221 /* Get connectivity information in canonical ordering for the local element */ 222 std::vector<moab::EntityHandle> connp; 223 merr = dmbp->hierarchy->get_connectivity(parent, dmbp->hlevel, connp); 224 MBERRNM(merr); 225 226 for (unsigned ic = 0; ic < connp.size(); ++ic) { 227 /* loop over each element connected to the adjacent vertex and update as needed */ 228 /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */ 229 if (found.find(connp[ic]) != found.end()) continue; /* make sure we don't double count shared vertices */ 230 if (dmbp->vghost->find(connp[ic]) != dmbp->vghost->end()) onz[ldof]++; /* update out-of-proc onz */ 231 else nnz[ldof]++; /* else local vertex */ 232 found.insert(connp[ic]); 233 } 234 } 235 } 236 237 for (int i = 0; i < nlsizc; i++) nnz[i] += 1; /* self count the node */ 238 239 ionz = onz[0]; 240 innz = nnz[0]; 241 for (int tc = 0; tc < nlsizc; tc++) { 242 // check for maximum allowed sparsity = fully dense 243 nnz[tc] = std::min(nlsizp, nnz[tc]); 244 onz[tc] = std::min(ngsizp - nlsizp, onz[tc]); 245 246 PetscCall(PetscInfo(NULL, " %d: NNZ = %d, ONZ = %d\n", tc, nnz[tc], onz[tc])); 247 248 innz = (innz < nnz[tc] ? nnz[tc] : innz); 249 ionz = (ionz < onz[tc] ? onz[tc] : ionz); 250 } 251 252 /* create interpolation matrix */ 253 PetscCall(MatCreate(PetscObjectComm((PetscObject)dmc), interpl)); 254 PetscCall(MatSetSizes(*interpl, nlsizc, nlsizp, ngsizc, ngsizp)); 255 PetscCall(MatSetType(*interpl, MATAIJ)); 256 PetscCall(MatSetFromOptions(*interpl)); 257 258 PetscCall(MatSeqAIJSetPreallocation(*interpl, innz, nnz)); 259 PetscCall(MatMPIAIJSetPreallocation(*interpl, innz, nnz, ionz, onz)); 260 261 /* clean up temporary memory */ 262 PetscCall(PetscFree2(nnz, onz)); 263 264 /* set up internal matrix data-structures */ 265 PetscCall(MatSetUp(*interpl)); 266 267 /* Define variables for assembly */ 268 std::vector<moab::EntityHandle> children; 269 std::vector<moab::EntityHandle> connp, connc; 270 std::vector<PetscReal> pcoords, ccoords, values_phi; 271 std::vector<PetscScalar> values_phi_scalar; 272 273 if (use_consistent_bases) { 274 const moab::EntityHandle ehandle = dmbp->elocal->front(); 275 276 merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children); 277 MBERRNM(merr); 278 279 /* Get connectivity and coordinates of the parent vertices */ 280 merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp); 281 MBERRNM(merr); 282 merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc); 283 MBERRNM(merr); 284 285 std::vector<PetscReal> natparam(3 * connc.size(), 0.0); 286 pcoords.resize(connp.size() * 3); 287 ccoords.resize(connc.size() * 3); 288 values_phi.resize(connp.size() * connc.size()); 289 /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 290 merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]); 291 MBERRNM(merr); 292 merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]); 293 MBERRNM(merr); 294 295 /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */ 296 for (unsigned tc = 0; tc < connc.size(); tc++) { 297 const PetscInt offset = tc * 3; 298 299 /* Scale ccoords relative to pcoords */ 300 PetscCall(DMMoabPToRMapping(dim, connp.size(), &pcoords[0], &ccoords[offset], &natparam[offset], &values_phi[connp.size() * tc])); 301 } 302 } else { 303 factor = std::pow(2.0 /*degree_P_for_refinement*/, (dmbc->hlevel - dmbp->hlevel) * dmbp->dim * 1.0); 304 } 305 306 /* TODO: Decipher the correct non-zero pattern. There is still some issue with onz allocation */ 307 PetscCall(MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 308 309 /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */ 310 values_phi_scalar = VecReal_to_VecScalar(values_phi); 311 for (moab::Range::iterator iter = dmbp->elocal->begin(); iter != dmbp->elocal->end(); iter++) { 312 const moab::EntityHandle ehandle = *iter; 313 314 /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 315 children.clear(); 316 connc.clear(); 317 merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children); 318 MBERRNM(merr); 319 320 /* Get connectivity and coordinates of the parent vertices */ 321 merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp); 322 MBERRNM(merr); 323 merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc); 324 MBERRNM(merr); 325 326 pcoords.resize(connp.size() * 3); 327 ccoords.resize(connc.size() * 3); 328 /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 329 merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]); 330 MBERRNM(merr); 331 merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]); 332 MBERRNM(merr); 333 334 std::vector<int> dofsp(connp.size()), dofsc(connc.size()); 335 /* TODO: specific to scalar system - use GetDofs */ 336 PetscCall(DMMoabGetDofsBlocked(dmp, connp.size(), &connp[0], &dofsp[0])); 337 PetscCall(DMMoabGetDofsBlocked(dmc, connc.size(), &connc[0], &dofsc[0])); 338 339 /* Compute the actual interpolation weights when projecting solution/residual between levels */ 340 if (use_consistent_bases) { 341 /* Use the cached values of natural parametric coordinates and basis pre-evaluated. 342 We are making an assumption here that UMR used in GMG to generate the hierarchy uses 343 the same template for all elements; This will fail for mixed element meshes (TRI/QUAD). 344 345 TODO: Fix the above assumption by caching data for families (especially for Tets and mixed meshes) 346 */ 347 348 /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */ 349 for (unsigned tc = 0; tc < connc.size(); tc++) { 350 /* TODO: Check if we should be using INSERT_VALUES instead */ 351 PetscCall(MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi_scalar[connp.size() * tc], ADD_VALUES)); 352 } 353 } else { 354 /* Compute the interpolation weights by determining distance of 1-ring 355 neighbor vertices from current vertex 356 357 This should be used only when FEM basis is not used for the discretization. 358 Else, the consistent interface to compute the basis function for interpolation 359 between the levels should be evaluated correctly to preserve convergence of GMG. 360 Shephard's basis will be terrible for any unsmooth problems. 361 */ 362 values_phi.resize(connp.size()); 363 for (unsigned tc = 0; tc < connc.size(); tc++) { 364 PetscReal normsum = 0.0; 365 std::vector<PetscScalar> values_phi_scalar2; 366 for (unsigned tp = 0; tp < connp.size(); tp++) { 367 values_phi[tp] = 0.0; 368 for (unsigned k = 0; k < 3; k++) values_phi[tp] += std::pow(pcoords[tp * 3 + k] - ccoords[k + tc * 3], dim); 369 if (values_phi[tp] < 1e-12) { 370 values_phi[tp] = 1e12; 371 } else { 372 //values_phi[tp] = std::pow(values_phi[tp], -1.0/dim); 373 values_phi[tp] = std::pow(values_phi[tp], -1.0); 374 normsum += values_phi[tp]; 375 } 376 } 377 for (unsigned tp = 0; tp < connp.size(); tp++) { 378 if (values_phi[tp] > 1e11) values_phi[tp] = factor * 0.5 / connp.size(); 379 else values_phi[tp] = factor * values_phi[tp] * 0.5 / (connp.size() * normsum); 380 } 381 values_phi_scalar2 = VecReal_to_VecScalar(values_phi); 382 PetscCall(MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi_scalar2[0], ADD_VALUES)); 383 } 384 } 385 } 386 if (vec) *vec = NULL; // TODO: <-- is it safe/appropriate? 387 PetscCall(MatAssemblyBegin(*interpl, MAT_FINAL_ASSEMBLY)); 388 PetscCall(MatAssemblyEnd(*interpl, MAT_FINAL_ASSEMBLY)); 389 PetscFunctionReturn(PETSC_SUCCESS); 390 } 391 392 // PetscClangLinter pragma ignore: -fdoc-* 393 /* 394 DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy 395 by succesively refining a coarse mesh, already defined in the `DM` object 396 provided by the user. 397 398 Collective 399 400 Input Parameter: 401 . dmb - The `DMMOAB` object 402 403 Output Parameters: 404 + nlevels - The number of levels of refinement needed to generate the hierarchy 405 - ldegrees - The degree of refinement at each level in the hierarchy 406 407 Level: beginner 408 */ 409 PETSC_EXTERN PetscErrorCode DMCreateInjection_Moab(DM dm1, DM dm2, VecScatter *ctx) 410 { 411 //DM_Moab *dmmoab; 412 413 PetscFunctionBegin; 414 PetscValidHeaderSpecific(dm1, DM_CLASSID, 1); 415 PetscValidHeaderSpecific(dm2, DM_CLASSID, 2); 416 //dmmoab = (DM_Moab*)(dm1)->data; 417 418 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n")); 419 PetscFunctionReturn(PETSC_SUCCESS); 420 } 421 422 static PetscErrorCode DMMoab_UMR_Private(DM dm, MPI_Comm comm, PetscBool refine, DM *dmref) 423 { 424 PetscInt i, dim; 425 DM dm2; 426 moab::ErrorCode merr; 427 DM_Moab *dmb = (DM_Moab *)dm->data, *dd2; 428 429 PetscFunctionBegin; 430 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 431 PetscAssertPointer(dmref, 4); 432 433 if ((dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine)) { 434 if (dmb->hlevel + 1 > dmb->nhlevels && refine) { 435 PetscCall(PetscInfo(NULL, "Invalid multigrid refinement hierarchy level specified (%" PetscInt_FMT "). MOAB UMR max levels = %" PetscInt_FMT ". Creating a NULL object.\n", dmb->hlevel + 1, dmb->nhlevels)); 436 } 437 if (dmb->hlevel - 1 < 0 && !refine) PetscCall(PetscInfo(NULL, "Invalid multigrid coarsen hierarchy level specified (%" PetscInt_FMT "). Creating a NULL object.\n", dmb->hlevel - 1)); 438 *dmref = NULL; 439 PetscFunctionReturn(PETSC_SUCCESS); 440 } 441 442 PetscCall(DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2)); 443 dd2 = (DM_Moab *)dm2->data; 444 445 dd2->mbiface = dmb->mbiface; 446 #ifdef MOAB_HAVE_MPI 447 dd2->pcomm = dmb->pcomm; 448 #endif 449 dd2->icreatedinstance = PETSC_FALSE; 450 dd2->nghostrings = dmb->nghostrings; 451 452 /* set the new level based on refinement/coarsening */ 453 if (refine) { 454 dd2->hlevel = dmb->hlevel + 1; 455 } else { 456 dd2->hlevel = dmb->hlevel - 1; 457 } 458 459 /* Copy the multilevel hierarchy pointers in MOAB */ 460 dd2->hierarchy = dmb->hierarchy; 461 dd2->nhlevels = dmb->nhlevels; 462 PetscCall(PetscMalloc1(dd2->nhlevels + 1, &dd2->hsets)); 463 for (i = 0; i <= dd2->nhlevels; i++) dd2->hsets[i] = dmb->hsets[i]; 464 dd2->fileset = dd2->hsets[dd2->hlevel]; 465 466 /* do the remaining initializations for DMMoab */ 467 dd2->bs = dmb->bs; 468 dd2->numFields = dmb->numFields; 469 dd2->rw_dbglevel = dmb->rw_dbglevel; 470 dd2->partition_by_rank = dmb->partition_by_rank; 471 PetscCall(PetscStrncpy(dd2->extra_read_options, dmb->extra_read_options, sizeof(dd2->extra_read_options))); 472 PetscCall(PetscStrncpy(dd2->extra_write_options, dmb->extra_write_options, sizeof(dd2->extra_write_options))); 473 dd2->read_mode = dmb->read_mode; 474 dd2->write_mode = dmb->write_mode; 475 476 /* set global ID tag handle */ 477 PetscCall(DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag)); 478 479 merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag); 480 MBERRNM(merr); 481 482 PetscCall(DMSetOptionsPrefix(dm2, ((PetscObject)dm)->prefix)); 483 PetscCall(DMGetDimension(dm, &dim)); 484 PetscCall(DMSetDimension(dm2, dim)); 485 486 /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 487 dm2->ops->creatematrix = dm->ops->creatematrix; 488 489 /* copy fill information if given */ 490 PetscCall(DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill)); 491 492 /* copy vector type information */ 493 PetscCall(DMSetMatType(dm2, dm->mattype)); 494 PetscCall(DMSetVecType(dm2, dm->vectype)); 495 dd2->numFields = dmb->numFields; 496 if (dmb->numFields) PetscCall(DMMoabSetFieldNames(dm2, dmb->numFields, dmb->fieldNames)); 497 498 PetscCall(DMSetFromOptions(dm2)); 499 500 /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */ 501 PetscCall(DMSetUp(dm2)); 502 503 *dmref = dm2; 504 PetscFunctionReturn(PETSC_SUCCESS); 505 } 506 507 // PetscClangLinter pragma ignore: -fdoc-* 508 /* 509 DMRefine_Moab - Generate a multi-level uniform refinement hierarchy 510 by succesively refining a coarse mesh, already defined in the `DM` object 511 provided by the user. 512 513 Collective 514 515 Input Parameters: 516 + dm - The `DMMOAB` object 517 - comm - the communicator to contain the new DM object (or `MPI_COMM_NULL`) 518 519 Output Parameter: 520 . dmf - the refined `DM`, or `NULL` 521 522 Level: developer 523 524 Note: 525 If no refinement was done, the return value is `NULL` 526 */ 527 PETSC_EXTERN PetscErrorCode DMRefine_Moab(DM dm, MPI_Comm comm, DM *dmf) 528 { 529 PetscFunctionBegin; 530 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 531 532 PetscCall(DMMoab_UMR_Private(dm, comm, PETSC_TRUE, dmf)); 533 PetscFunctionReturn(PETSC_SUCCESS); 534 } 535 536 // PetscClangLinter pragma ignore: -fdoc-* 537 /* 538 DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy 539 by succesively refining a coarse mesh, already defined in the `DM` object 540 provided by the user. 541 542 Collective 543 544 Input Parameters: 545 + dm - The `DMMOAB` object 546 - comm - the communicator to contain the new `DM` object (or `MPI_COMM_NULL`) 547 548 Output Parameter: 549 . dmc - the coarsened `DM`, or `NULL` 550 551 Level: developer 552 553 Note: 554 If no coarsening was done, the return value is `NULL` 555 */ 556 PETSC_EXTERN PetscErrorCode DMCoarsen_Moab(DM dm, MPI_Comm comm, DM *dmc) 557 { 558 PetscFunctionBegin; 559 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 560 PetscCall(DMMoab_UMR_Private(dm, comm, PETSC_FALSE, dmc)); 561 PetscFunctionReturn(PETSC_SUCCESS); 562 } 563