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