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