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