1 #include <petsc-private/dmmbimpl.h> /*I "petscdmmoab.h" I*/ 2 3 #include <petscdmmoab.h> 4 #include <MBTagConventions.hpp> 5 #include <moab/NestedRefine.hpp> 6 7 #undef __FUNCT__ 8 #define __FUNCT__ "DMMoabGenerateHierarchy" 9 /*@ 10 DMMoabGenerateHierarchy - Generate a multi-level uniform refinement hierarchy 11 by succesively refining a coarse mesh, already defined in the DM object 12 provided by the user. 13 14 Collective on MPI_Comm 15 16 Input Parameter: 17 + dmb - The DMMoab object 18 19 Output Parameter: 20 + nlevels - The number of levels of refinement needed to generate the hierarchy 21 . ldegrees - The degree of refinement at each level in the hierarchy 22 23 Level: beginner 24 25 .keywords: DMMoab, create, refinement 26 @*/ 27 PetscErrorCode DMMoabGenerateHierarchy(DM dm,PetscInt nlevels,PetscInt *ldegrees) 28 { 29 DM_Moab *dmmoab; 30 PetscErrorCode ierr; 31 moab::ErrorCode merr; 32 PetscInt *pdegrees,i; 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->fileset); 49 50 /* generate the mesh hierarchy */ 51 merr = dmmoab->hierarchy->generate_mesh_hierarchy(pdegrees, nlevels, &dmmoab->hsets);MBERRNM(merr); 52 53 if (!ldegrees) { 54 ierr = PetscFree(pdegrees);CHKERRQ(ierr); 55 } 56 PetscFunctionReturn(0); 57 } 58 59 #undef __FUNCT__ 60 #define __FUNCT__ "DMRefineHierarchy_Moab" 61 /*@ 62 DMRefineHierarchy_Moab - Generate a multi-level uniform refinement hierarchy 63 by succesively refining a coarse mesh, already defined in the DM object 64 provided by the user. 65 66 Collective on MPI_Comm 67 68 Input Parameter: 69 . dmb - The DMMoab object 70 71 Output Parameter: 72 . nlevels - The number of levels of refinement needed to generate the hierarchy 73 + ldegrees - The degree of refinement at each level in the hierarchy 74 75 Level: beginner 76 77 .keywords: DMMoab, create, refinement 78 @*/ 79 PetscErrorCode DMRefineHierarchy_Moab(DM dm,PetscInt nlevels,DM dmf[]) 80 { 81 DM_Moab *dmmoab; 82 PetscErrorCode ierr; 83 moab::ErrorCode merr; 84 PetscInt i; 85 86 PetscFunctionBegin; 87 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 88 dmmoab = (DM_Moab*)(dm)->data; 89 90 for (i=0; i<nlevels; i++) { 91 dmf[i] = dm; 92 PetscObjectReference((PetscObject)dm); 93 } 94 95 PetscPrintf(PETSC_COMM_WORLD, "[DMRefineHierarchy_Moab] :: Placeholder\n"); 96 PetscFunctionReturn(0); 97 } 98 99 #undef __FUNCT__ 100 #define __FUNCT__ "DMCoarsenHierarchy_Moab" 101 /*@ 102 DMCoarsenHierarchy_Moab - Generate a multi-level uniform refinement hierarchy 103 by succesively refining a coarse mesh, already defined in the DM object 104 provided by the user. 105 106 Collective on MPI_Comm 107 108 Input Parameter: 109 . dmb - The DMMoab object 110 111 Output Parameter: 112 . nlevels - The number of levels of refinement needed to generate the hierarchy 113 + ldegrees - The degree of refinement at each level in the hierarchy 114 115 Level: beginner 116 117 .keywords: DMMoab, create, refinement 118 @*/ 119 PetscErrorCode DMCoarsenHierarchy_Moab(DM dm,PetscInt nlevels,DM dmc[]) 120 { 121 DM_Moab *dmmoab; 122 PetscErrorCode ierr; 123 moab::ErrorCode merr; 124 PetscInt i; 125 126 PetscFunctionBegin; 127 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 128 dmmoab = (DM_Moab*)(dm)->data; 129 130 for (i=0; i<nlevels; i++) { 131 dmc[i] = dm; 132 PetscObjectReference((PetscObject)dm); 133 } 134 135 PetscPrintf(PETSC_COMM_WORLD, "[DMCoarsenHierarchy_Moab] :: Placeholder\n"); 136 PetscFunctionReturn(0); 137 } 138 139 140 #undef __FUNCT__ 141 #define __FUNCT__ "DMCreateInterpolation_Moab" 142 /*@ 143 DMCreateInterpolation_Moab - Generate a multi-level uniform refinement hierarchy 144 by succesively refining a coarse mesh, already defined in the DM object 145 provided by the user. 146 147 Collective on MPI_Comm 148 149 Input Parameter: 150 . dmb - The DMMoab object 151 152 Output Parameter: 153 . nlevels - The number of levels of refinement needed to generate the hierarchy 154 + ldegrees - The degree of refinement at each level in the hierarchy 155 156 Level: beginner 157 158 .keywords: DMMoab, create, refinement 159 @*/ 160 PetscErrorCode DMCreateInterpolation_Moab(DM dm1,DM dm2,Mat* interpl,Vec* vec) 161 { 162 DM_Moab *dmb1, *dmb2; 163 PetscErrorCode ierr; 164 moab::ErrorCode merr; 165 Vec unitv; 166 PetscInt dim,dofs_per_element=4; 167 PetscReal unitval=1.0; 168 PetscBool eonbnd,dbdry[27]; 169 std::vector<int> bndrows; 170 171 PetscFunctionBegin; 172 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 173 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 174 dmb1 = (DM_Moab*)(dm1)->data; 175 dmb2 = (DM_Moab*)(dm2)->data; 176 177 ierr = MatCreate(PetscObjectComm((PetscObject)dm1), interpl);CHKERRQ(ierr); 178 ierr = MatSetType(*interpl, MATMPIAIJ);CHKERRQ(ierr); 179 ierr = MatSetSizes(*interpl, dmb2->nloc, dmb1->nloc, dmb2->n, dmb1->n);CHKERRQ(ierr); 180 181 /* TODO: This is a hack for the rectangular system - decipher NNZ pattern better */ 182 ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr); 183 184 ierr = MatSeqAIJSetPreallocation(*interpl,dmb2->nloc,0);CHKERRQ(ierr); 185 ierr = MatMPIAIJSetPreallocation(*interpl,dmb2->nloc,0,dmb2->nghost,0);CHKERRQ(ierr); 186 187 /* set up internal matrix data-structures */ 188 ierr = MatSetUp(*interpl);CHKERRQ(ierr); 189 190 ierr = MatZeroEntries(*interpl);CHKERRQ(ierr); 191 192 ierr = DMGetDimension(dm1, &dim);CHKERRQ(ierr); 193 194 PetscPrintf(PETSC_COMM_WORLD, "Creating a %D DIM matrix that is %D X %D and setting diagonal to 1.0", dim, dmb1->nloc, dmb2->nloc); 195 196 double factor = std::pow(2.0,(dmb2->hlevel-dmb1->hlevel)*dmb1->dim*1.0); 197 198 //Loop through the remaining vertices. These vertices appear only on the current refined_level. 199 for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) { 200 201 const moab::EntityHandle ehandle = *iter; 202 std::vector<moab::EntityHandle> children; 203 std::vector<moab::EntityHandle> connp, connc; 204 205 /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 206 merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr); 207 208 /* Get connectivity and coordinates of the parent vertices */ 209 merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr); 210 for (int ic=0; ic < children.size(); ic++) { 211 std::vector<moab::EntityHandle> tconnc; 212 /* Get coordinates of the parent vertices in canonical order */ 213 merr = dmb1->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr); 214 for (int tc=0; tc<tconnc.size(); tc++) { 215 connc.push_back(tconnc[tc]); 216 } 217 } 218 219 std::vector<double> pcoords(connp.size()*3), ccoords(connc.size()*3), values_phi(connc.size()); 220 /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 221 merr = dmb1->hierarchy->get_coordinates(&connp[0], connp.size(), dmb1->hlevel, &pcoords[0]);MBERRNM(merr); 222 merr = dmb2->hierarchy->get_coordinates(&connc[0], connc.size(), dmb2->hlevel, &ccoords[0]);MBERRNM(merr); 223 224 std::vector<int> dofsp(connp.size()), dofsc(connc.size()); 225 /* TODO: specific to scalar system - use GetDofs */ 226 ierr = DMMoabGetFieldDofs(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr); 227 ierr = DMMoabGetFieldDofs(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr); 228 229 // Compute the interoplation weights by determining distance of 1-ring neighbor vertices 230 // from current vertex 231 for (int i=0;i<connp.size(); i++) { 232 double normsum=0.0; 233 for (int j=0;j<connc.size(); j++) { 234 int offset=j; 235 values_phi[offset] = 0.0; 236 for (int k=0;k<3; k++) 237 values_phi[offset] += (pcoords[i*3+k]-ccoords[k+j*3])*(pcoords[i*3+k]-ccoords[k+j*3]); 238 if (values_phi[offset] < 1e-12) { 239 values_phi[offset] = 1e12; 240 // PetscPrintf(PETSC_COMM_WORLD, "Found coarse and fine restrictive space: %D, %D\n",dofsp[i],dofsc[j]); 241 } 242 else { 243 values_phi[offset] = 1.0/(values_phi[offset]); 244 normsum += values_phi[offset]; 245 } 246 } 247 //PetscPrintf(PETSC_COMM_WORLD, "\nRow [%D]: ", dofsp[i]); 248 for (int j=0;j<connc.size(); j++) { 249 //int offset=i*connc.size()+j; 250 int offset=j; 251 if (values_phi[offset] > 1e11) 252 values_phi[offset] = factor*0.5/connc.size(); 253 else 254 values_phi[offset] = factor*values_phi[offset]*0.5/(connc.size()*normsum); 255 //PetscPrintf(PETSC_COMM_WORLD, "%D : %g, ", dofsc[j], values_phi[offset]); 256 } 257 ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr); 258 } 259 260 /* check if element is on the boundary */ 261 //ierr = DMMoabIsEntityOnBoundary(dm1,ehandle,&eonbnd);CHKERRQ(ierr); 262 ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry);CHKERRQ(ierr); 263 eonbnd=PETSC_FALSE; 264 for (int i=0; i< connc.size(); ++i) 265 if (dbdry[i]) eonbnd=PETSC_TRUE; 266 267 values_phi.clear(); 268 values_phi.resize(connp.size()); 269 /* apply dirichlet boundary conditions */ 270 if (eonbnd) { 271 272 ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 273 ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 274 /* get the list of nodes on boundary so that we can enforce dirichlet conditions strongly */ 275 //ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry);CHKERRQ(ierr); 276 for (int i=0; i < connc.size(); i++) { 277 if (dmb2->hierarchy->is_boundary_vertex(connc[i])) { /* dirichlet node */ 278 /* think about strongly imposing dirichlet */ 279 //bndrows.push_back(dofsc[i]); 280 281 ierr = MatSetValues(*interpl, 1, &dofsc[i], connp.size(), &dofsp[0], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 282 //values_phi[0]=1.0; 283 //ierr = MatSetValues(*interpl, 1, &dofsc[i], 1, &dofsc[i], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 284 } 285 } 286 287 ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 288 ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 289 } 290 291 //get interpolation weights 292 //ierr = Compute_Quad4_Basis(pcoords, 1, coord, values_phi);CHKERRQ(ierr); 293 /*for (int j=0;j<dofs_per_element; j++) 294 std::cout<<"values "<<values_phi[j]<<std::endl;*/ 295 296 //get row and column indices, zero weights are ignored 297 /* 298 int nz_ind = 0; 299 idx = dmb2->vowned->index(vhandle); 300 for (int j=0;j<dofs_per_element; j++){ 301 idy[nz_ind] = dmb1->vowned->index(connectivity[j]); 302 PetscPrintf(PETSC_COMM_WORLD, "Finding coarse connectivity vertex %D associated with [%D, %D] - set to %D\n", connectivity[j], parent.size(), vhandle, idy[nz_ind]); 303 //values_phi[nz_ind] = values_phi[j]; 304 nz_ind = nz_ind+1; 305 } 306 */ 307 308 //PetscPrintf(PETSC_COMM_WORLD, "Setting basis for row %D: %g %g %g %g\n", idx, values_phi[0], values_phi[1], values_phi[2], values_phi[3]); 309 //PetscPrintf(PETSC_COMM_WORLD, "Setting entries for row %D: %D %D %D %D\n", idx, idy[0], idy[1], idy[2], idy[3]); 310 //ierr = MatSetValues(*interpl, nz_ind, idy, 1, &idx, values_phi, INSERT_VALUES);CHKERRQ(ierr); 311 //ierr = MatSetValues(*interpl, connp.size(), dofsp, connc.size(), dofsc, &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 312 } 313 314 //PetscPrintf(PETSC_COMM_WORLD, "[Boundary vertices = %D] :: A few: %D %D %D %D \n", bndrows.size(), bndrows[0], bndrows[1], bndrows[2], bndrows[3]); 315 //ierr = MatZeroRows(*interpl, bndrows.size(), &bndrows[0], 1.0, 0, 0);CHKERRQ(ierr); 316 317 ierr = MatAssemblyBegin(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 318 ierr = MatAssemblyEnd(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 319 //MatView(*interpl, 0); 320 321 PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInterpolation_Moab] :: Placeholder\n"); 322 PetscFunctionReturn(0); 323 } 324 325 #undef __FUNCT__ 326 #define __FUNCT__ "DMCreateInjection_Moab" 327 /*@ 328 DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy 329 by succesively refining a coarse mesh, already defined in the DM object 330 provided by the user. 331 332 Collective on MPI_Comm 333 334 Input Parameter: 335 . dmb - The DMMoab object 336 337 Output Parameter: 338 . nlevels - The number of levels of refinement needed to generate the hierarchy 339 + ldegrees - The degree of refinement at each level in the hierarchy 340 341 Level: beginner 342 343 .keywords: DMMoab, create, refinement 344 @*/ 345 PetscErrorCode DMCreateInjection_Moab(DM dm1,DM dm2,VecScatter* ctx) 346 { 347 DM_Moab *dmmoab; 348 PetscErrorCode ierr; 349 moab::ErrorCode merr; 350 351 PetscFunctionBegin; 352 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 353 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 354 dmmoab = (DM_Moab*)(dm1)->data; 355 356 PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n"); 357 PetscFunctionReturn(0); 358 } 359 360 361 #undef __FUNCT__ 362 #define __FUNCT__ "DM_UMR_Moab_Private" 363 PetscErrorCode DM_UMR_Moab_Private(DM dm,MPI_Comm comm,PetscBool refine,DM *dmref) 364 { 365 PetscErrorCode ierr; 366 PetscInt M,N,P,i,dim; 367 DM dm2; 368 moab::ErrorCode merr; 369 DM_Moab *dmb = (DM_Moab*)dm->data,*dd2; 370 371 PetscFunctionBegin; 372 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 373 PetscValidPointer(dmref,4); 374 375 //if (dmb->hlevel+1 > dmb->nhlevels && refine) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Invalid multigrid refinement hierarchy level specified (%D). MOAB UMR max levels = %D\n",dmb->hlevel+1,dmb->nhlevels); 376 //if (dmb->hlevel-1 < 0 && !refine) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Invalid multigrid coarsen hierarchy level specified (%D).\n",dmb->hlevel-1); 377 378 if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) ) { 379 dmref = PETSC_NULL; 380 PetscFunctionReturn(0); 381 } 382 383 ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr); 384 dd2 = (DM_Moab*)dm2->data; 385 386 dd2->mbiface = dmb->mbiface; 387 dd2->pcomm = dmb->pcomm; 388 dd2->icreatedinstance = PETSC_FALSE; 389 390 /* set the new level based on refinement/coarsening */ 391 if (refine) { 392 dd2->hlevel=dmb->hlevel+1; 393 } 394 else { 395 dd2->hlevel=dmb->hlevel-1; 396 } 397 398 /* Copy the multilevel hierarchy pointers in MOAB */ 399 dd2->hierarchy = dmb->hierarchy; 400 dd2->nhlevels = dmb->nhlevels; 401 ierr = PetscMalloc1(dd2->nhlevels+1,&dd2->hsets);CHKERRQ(ierr); 402 for (i=0; i<=dd2->nhlevels; i++) { 403 dd2->hsets[i]=dmb->hsets[i]; 404 } 405 dd2->fileset = dd2->hsets[dd2->hlevel]; 406 407 /* do the remaining initializations for DMMoab */ 408 dd2->bs = dmb->bs; 409 dd2->numFields = dmb->numFields; 410 dd2->rw_dbglevel = dmb->rw_dbglevel; 411 dd2->partition_by_rank = dmb->partition_by_rank; 412 ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr); 413 ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr); 414 dd2->read_mode = dmb->read_mode; 415 dd2->write_mode = dmb->write_mode; 416 417 /* set global ID tag handle */ 418 ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr); 419 420 merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr); 421 422 //ierr = DMMoabCreateMoab(PetscObjectComm((PetscObject)dm),dmb->mbiface,dmb->pcomm,&dmb->ltog_tag,NULL,&dm2);CHKERRQ(ierr); 423 ierr = DMSetOptionsPrefix(dm2,((PetscObject)dm)->prefix);CHKERRQ(ierr); 424 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 425 ierr = DMSetDimension(dm2,dim);CHKERRQ(ierr); 426 427 /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 428 dm2->ops->creatematrix = dm->ops->creatematrix; 429 430 /* copy fill information if given */ 431 ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr); 432 433 /* copy vector type information */ 434 ierr = DMSetMatType(dm2,dm->mattype);CHKERRQ(ierr); 435 ierr = DMSetVecType(dm2,dm->vectype);CHKERRQ(ierr); 436 ierr = DMMoabSetFieldNames(dm2,dmb->numFields,dmb->fieldNames);CHKERRQ(ierr); 437 438 ierr = DMSetFromOptions(dm2);CHKERRQ(ierr); 439 440 /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */ 441 ierr = DMSetUp(dm2);CHKERRQ(ierr); 442 443 *dmref = dm2; 444 PetscFunctionReturn(0); 445 } 446 447 448 #undef __FUNCT__ 449 #define __FUNCT__ "DMRefine_Moab" 450 /*@ 451 DMRefine_Moab - Generate a multi-level uniform refinement hierarchy 452 by succesively refining a coarse mesh, already defined in the DM object 453 provided by the user. 454 455 Collective on MPI_Comm 456 457 Input Parameter: 458 . dmb - The DMMoab object 459 460 Output Parameter: 461 . nlevels - The number of levels of refinement needed to generate the hierarchy 462 + ldegrees - The degree of refinement at each level in the hierarchy 463 464 Level: beginner 465 466 .keywords: DMMoab, create, refinement 467 @*/ 468 PetscErrorCode DMRefine_Moab(DM dm,MPI_Comm comm,DM* dmf) 469 { 470 DM_Moab *dmmoab; 471 PetscErrorCode ierr; 472 moab::ErrorCode merr; 473 474 PetscFunctionBegin; 475 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 476 dmmoab = (DM_Moab*)(dm)->data; 477 478 PetscPrintf(PETSC_COMM_WORLD, "[DMRefine_Moab] :: Placeholder\n"); 479 ierr = DM_UMR_Moab_Private(dm,comm,PETSC_TRUE,dmf);CHKERRQ(ierr); 480 PetscPrintf(PETSC_COMM_WORLD, "[DMRefine_Moab] :: All done\n"); 481 PetscFunctionReturn(0); 482 } 483 484 #undef __FUNCT__ 485 #define __FUNCT__ "DMCoarsen_Moab" 486 /*@ 487 DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy 488 by succesively refining a coarse mesh, already defined in the DM object 489 provided by the user. 490 491 Collective on MPI_Comm 492 493 Input Parameter: 494 . dmb - The DMMoab object 495 496 Output Parameter: 497 . nlevels - The number of levels of refinement needed to generate the hierarchy 498 + ldegrees - The degree of refinement at each level in the hierarchy 499 500 Level: beginner 501 502 .keywords: DMMoab, create, refinement 503 @*/ 504 PetscErrorCode DMCoarsen_Moab(DM dm,MPI_Comm comm,DM* dmc) 505 { 506 DM_Moab *dmmoab; 507 PetscErrorCode ierr; 508 moab::ErrorCode merr; 509 510 PetscFunctionBegin; 511 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 512 dmmoab = (DM_Moab*)(dm)->data; 513 514 PetscPrintf(PETSC_COMM_WORLD, "[DMCoarsen_Moab] :: Placeholder\n"); 515 ierr = DM_UMR_Moab_Private(dm,comm,PETSC_FALSE,dmc);CHKERRQ(ierr); 516 PetscPrintf(PETSC_COMM_WORLD, "[DMCoarsen_Moab] :: All done\n"); 517 PetscFunctionReturn(0); 518 } 519