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