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 std::vector<moab::EntityHandle> hsets; 34 35 PetscFunctionBegin; 36 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 37 dmmoab = (DM_Moab*)(dm)->data; 38 39 if (!ldegrees) { 40 ierr = PetscMalloc1(nlevels,&pdegrees);CHKERRQ(ierr); 41 for (i=0; i<nlevels; i++) pdegrees[i]=2; /* default = Degree 2 refinement */ 42 } 43 else pdegrees=ldegrees; 44 45 /* initialize set level refinement data for hierarchy */ 46 dmmoab->nhlevels=nlevels; 47 48 /* Instantiate the nested refinement class */ 49 dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset); 50 51 ierr = PetscMalloc1(nlevels+1,&dmmoab->hsets);CHKERRQ(ierr); 52 53 /* generate the mesh hierarchy */ 54 merr = dmmoab->hierarchy->generate_mesh_hierarchy(pdegrees, nlevels, hsets);MBERRNM(merr); 55 56 PetscInfo2(NULL, "Exchanging ghost cells (dim %d) with %d rings\n",dmmoab->dim,dmmoab->nghostrings); 57 for (i=0; i<=nlevels; i++) { 58 /* resolve the shared entities by exchanging information to adjacent processors */ 59 merr = dmmoab->pcomm->exchange_ghost_cells(dmmoab->dim,0,dmmoab->nghostrings,dmmoab->dim,true,false,&hsets[i]);MBERRV(dmmoab->mbiface,merr); 60 61 dmmoab->hsets[i]=hsets[i]; 62 } 63 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 dm1,DM dm2,Mat* interpl,Vec* vec) 159 { 160 DM_Moab *dmb1, *dmb2; 161 PetscErrorCode ierr; 162 moab::ErrorCode merr; 163 PetscInt dim; 164 PetscReal factor; 165 PetscBool eonbnd; 166 PetscInt innz, *nnz, ionz, *onz; 167 PetscInt nlsiz1, nlsiz2, nlghsiz1, nlghsiz2, ngsiz1, ngsiz2; 168 std::vector<int> bndrows; 169 std::vector<PetscBool> dbdry; 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 nlsiz1 = dmb1->nloc*dmb1->numFields; 177 nlsiz2 = dmb2->nloc*dmb2->numFields; 178 ngsiz1 = dmb1->n*dmb1->numFields; 179 ngsiz2 = dmb2->n*dmb2->numFields; 180 nlghsiz1 = (dmb1->nloc+dmb1->nghost)*dmb1->numFields; 181 nlghsiz2 = (dmb2->nloc+dmb2->nghost)*dmb2->numFields; 182 183 PetscInfo4(dm1,"Creating interpolation matrix %D X %D to apply transformation between levels %D -> %D.\n",dmb1->nloc,dmb2->nloc,dmb1->hlevel,dmb2->hlevel); 184 185 ierr = MatCreate(PetscObjectComm((PetscObject)dm1), interpl);CHKERRQ(ierr); 186 ierr = MatSetType(*interpl, dm1->mattype);CHKERRQ(ierr); 187 ierr = MatSetSizes(*interpl, nlsiz2, nlsiz1, ngsiz2, ngsiz1);CHKERRQ(ierr); 188 189 /* TODO: This is a hack for the rectangular system - decipher NNZ pattern better */ 190 ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr); 191 192 /* allocate the nnz, onz arrays based on block size and local nodes */ 193 ierr = PetscCalloc2(nlghsiz2,&nnz,nlghsiz2,&onz);CHKERRQ(ierr); 194 195 /* Loop through the local elements and compute the relation between the current parent and the refined_level. */ 196 for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) { 197 198 const moab::EntityHandle ehandle = *iter; 199 std::vector<moab::EntityHandle> children; 200 std::vector<moab::EntityHandle> connp, connc; 201 202 /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 203 merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr); 204 205 /* Get connectivity and coordinates of the parent vertices */ 206 merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr); 207 for (unsigned ic=0; ic < children.size(); ic++) { 208 std::vector<moab::EntityHandle> tconnc; 209 /* Get coordinates of the parent vertices in canonical order */ 210 merr = dmb2->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr); 211 for (unsigned tc=0; tc<tconnc.size(); tc++) { 212 if (std::find(connc.begin(), connc.end(), tconnc[tc]) == connc.end()) 213 connc.push_back(tconnc[tc]); 214 } 215 } 216 217 std::vector<int> dofsp(connp.size()), dofsc(connc.size()); 218 /* TODO: specific to scalar system - use GetDofs */ 219 //ierr = DMMoabGetFieldDofs(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr); 220 ierr = DMMoabGetFieldDofsLocal(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr); 221 222 for (unsigned tp=0;tp<connp.size(); tp++) { 223 // ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr); 224 if (dmb1->vowned->find(connp[tp]) != dmb1->vowned->end()) { 225 for (unsigned tc=0;tc<connc.size(); tc++) nnz[dofsc[tc]]++; 226 } 227 else if (dmb1->vghost->find(connp[tp]) != dmb1->vghost->end()) { 228 for (unsigned tc=0;tc<connc.size(); tc++) onz[dofsc[tc]]++; 229 } 230 else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid entity in parent level %D\n", connc[tp]); 231 } 232 for(int tc = 0; tc < connc.size(); tc++) { 233 if (dmb2->vowned->find(connc[tc]) != dmb2->vowned->end()) nnz[dofsc[tc]]++; 234 else if (dmb2->vghost->find(connc[tc]) != dmb2->vghost->end()) onz[dofsc[tc]]++; 235 else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid entity in child level %D\n", connc[tc]); 236 } 237 } 238 239 /* 240 int i=0; 241 std::vector<moab::EntityHandle> adjs; 242 for(moab::Range::iterator iter = dmb1->vowned->begin(); iter!= dmb1->vowned->end(); iter++, i++) { 243 merr = dmb1->hierarchy->get_adjacencies(*iter, 0, adjs);MBERRNM(merr); 244 nnz[i] -= adjs.size(); 245 adjs.clear(); 246 } 247 for(moab::Range::iterator iter = dmb1->vghost->begin(); iter!= dmb1->vghost->end(); iter++, i++) { 248 //merr = dmb1->hierarchy->get_adjacencies(&(*iter), 1, 0, false, adjs, moab::Interface::UNION);MBERRNM(merr); 249 merr = dmb1->hierarchy->get_adjacencies(*iter, 0, adjs);MBERRNM(merr); 250 onz[i] -= adjs.size(); 251 adjs.clear(); 252 } 253 */ 254 255 ionz=onz[0]; 256 innz=nnz[0]; 257 for (int tc=0; tc < nlsiz2; tc++) { 258 // check for maximum allowed sparsity = fully dense 259 nnz[tc] = std::min(nlsiz1,nnz[tc]); 260 onz[tc] = std::min(nlsiz1,onz[tc]); 261 262 innz = (innz < nnz[tc] ? nnz[tc] : innz); 263 ionz = (ionz < onz[tc] ? onz[tc] : ionz); 264 //PetscPrintf(PETSC_COMM_SELF, "[%D] NNZ = %D, ONZ = %D\n", tc, nnz[tc], onz[tc]); 265 } 266 //PetscPrintf(PETSC_COMM_SELF, "Final: INNZ = %D, IONZ = %D\n", innz, ionz); 267 268 ierr = MatSeqAIJSetPreallocation(*interpl,innz,nnz);CHKERRQ(ierr); 269 ierr = MatMPIAIJSetPreallocation(*interpl,innz,nnz,ionz,onz);CHKERRQ(ierr); 270 271 /* clean up temporary memory */ 272 ierr = PetscFree2(nnz,onz);CHKERRQ(ierr); 273 274 /* set up internal matrix data-structures */ 275 ierr = MatSetUp(*interpl);CHKERRQ(ierr); 276 ierr = MatZeroEntries(*interpl);CHKERRQ(ierr); 277 278 ierr = DMGetDimension(dm1, &dim);CHKERRQ(ierr); 279 280 factor = std::pow(2.0 /*degree_P_for_refinement*/,(dmb2->hlevel-dmb1->hlevel)*dmb1->dim*1.0); 281 282 ierr = MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);CHKERRQ(ierr); 283 284 /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */ 285 for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) { 286 287 const moab::EntityHandle ehandle = *iter; 288 std::vector<moab::EntityHandle> children; 289 std::vector<moab::EntityHandle> connp, connc; 290 291 /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 292 merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr); 293 294 /* Get connectivity and coordinates of the parent vertices */ 295 merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr); 296 for (unsigned ic=0; ic < children.size(); ic++) { 297 std::vector<moab::EntityHandle> tconnc; 298 /* Get coordinates of the parent vertices in canonical order */ 299 merr = dmb1->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr); 300 for (unsigned tc=0; tc<tconnc.size(); tc++) { 301 connc.push_back(tconnc[tc]); 302 } 303 } 304 305 std::vector<double> pcoords(connp.size()*3), ccoords(connc.size()*3), values_phi(connc.size()); 306 /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 307 merr = dmb1->hierarchy->get_coordinates(&connp[0], connp.size(), dmb1->hlevel, &pcoords[0]);MBERRNM(merr); 308 merr = dmb2->hierarchy->get_coordinates(&connc[0], connc.size(), dmb2->hlevel, &ccoords[0]);MBERRNM(merr); 309 310 std::vector<int> dofsp(connp.size()), dofsc(connc.size()); 311 /* TODO: specific to scalar system - use GetDofs */ 312 ierr = DMMoabGetFieldDofsLocal(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr); 313 ierr = DMMoabGetFieldDofsLocal(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr); 314 315 /* Compute the interpolation weights by determining distance of 1-ring 316 neighbor vertices from current vertex */ 317 for (unsigned i=0;i<connp.size(); i++) { 318 double normsum=0.0; 319 for (unsigned j=0;j<connc.size(); j++) { 320 values_phi[j] = 0.0; 321 for (unsigned k=0;k<3; k++) 322 values_phi[j] += std::pow(pcoords[i*3+k]-ccoords[k+j*3], dim); 323 if (values_phi[j] < 1e-12) { 324 values_phi[j] = 1e12; 325 } 326 else { 327 //values_phi[j] = std::pow(values_phi[j], -1.0/dim); 328 values_phi[j] = std::pow(values_phi[j], -1.0); 329 normsum += values_phi[j]; 330 } 331 } 332 for (unsigned j=0;j<connc.size(); j++) { 333 if (values_phi[j] > 1e11) 334 values_phi[j] = factor*0.5/connc.size(); 335 else 336 values_phi[j] = factor*values_phi[j]*0.5/(connc.size()*normsum); 337 } 338 ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr); 339 } 340 341 /* check if element is on the boundary */ 342 //ierr = DMMoabIsEntityOnBoundary(dm1,ehandle,&eonbnd);CHKERRQ(ierr); 343 dbdry.resize(connc.size()); 344 ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry.data());CHKERRQ(ierr); 345 eonbnd=PETSC_FALSE; 346 for (unsigned i=0; i< connc.size(); ++i) 347 if (dbdry[i]) eonbnd=PETSC_TRUE; 348 349 values_phi.clear(); 350 values_phi.resize(connp.size()); 351 /* apply dirichlet boundary conditions */ 352 if (eonbnd) { 353 354 ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 355 ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 356 /* get the list of nodes on boundary so that we can enforce dirichlet conditions strongly */ 357 //ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry);CHKERRQ(ierr); 358 for (unsigned i=0; i < connc.size(); i++) { 359 if (dbdry[i]) { /* dirichlet node */ 360 /* think about strongly imposing dirichlet */ 361 //bndrows.push_back(dofsc[i]); 362 363 ierr = MatSetValues(*interpl, 1, &dofsc[i], connp.size(), &dofsp[0], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 364 //values_phi[0]=1.0; 365 //ierr = MatSetValues(*interpl, 1, &dofsc[i], 1, &dofsc[i], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 366 } 367 } 368 369 ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 370 ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 371 } 372 373 //get interpolation weights 374 //ierr = Compute_Quad4_Basis(pcoords, 1, coord, values_phi);CHKERRQ(ierr); 375 // for (int j=0;j<dofs_per_element; j++) 376 // std::cout<<"values "<<values_phi[j]<<std::endl; 377 378 //get row and column indices, zero weights are ignored 379 /* 380 int nz_ind = 0; 381 idx = dmb2->vowned->index(vhandle); 382 for (int j=0;j<dofs_per_element; j++){ 383 idy[nz_ind] = dmb1->vowned->index(connectivity[j]); 384 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]); 385 //values_phi[nz_ind] = values_phi[j]; 386 nz_ind = nz_ind+1; 387 } 388 */ 389 390 //ierr = MatSetValues(*interpl, nz_ind, idy, 1, &idx, values_phi, INSERT_VALUES);CHKERRQ(ierr); 391 //ierr = MatSetValues(*interpl, connp.size(), dofsp, connc.size(), dofsc, &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 392 } 393 394 //PetscPrintf(PETSC_COMM_WORLD, "[Boundary vertices = %D] :: A few: %D %D %D %D \n", bndrows.size(), bndrows[0], bndrows[1], bndrows[2], bndrows[3]); 395 //ierr = MatZeroRows(*interpl, bndrows.size(), &bndrows[0], 1.0, 0, 0);CHKERRQ(ierr); 396 397 ierr = MatAssemblyBegin(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 398 ierr = MatAssemblyEnd(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 399 PetscFunctionReturn(0); 400 } 401 402 #undef __FUNCT__ 403 #define __FUNCT__ "DMCreateInjection_Moab" 404 /*@ 405 DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy 406 by succesively refining a coarse mesh, already defined in the DM object 407 provided by the user. 408 409 Collective on MPI_Comm 410 411 Input Parameter: 412 . dmb - The DMMoab object 413 414 Output Parameter: 415 . nlevels - The number of levels of refinement needed to generate the hierarchy 416 + ldegrees - The degree of refinement at each level in the hierarchy 417 418 Level: beginner 419 420 .keywords: DMMoab, create, refinement 421 @*/ 422 PetscErrorCode DMCreateInjection_Moab(DM dm1,DM dm2,VecScatter* ctx) 423 { 424 //DM_Moab *dmmoab; 425 426 PetscFunctionBegin; 427 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 428 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 429 //dmmoab = (DM_Moab*)(dm1)->data; 430 431 PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n"); 432 PetscFunctionReturn(0); 433 } 434 435 436 #undef __FUNCT__ 437 #define __FUNCT__ "DM_UMR_Moab_Private" 438 PetscErrorCode DM_UMR_Moab_Private(DM dm,MPI_Comm comm,PetscBool refine,DM *dmref) 439 { 440 PetscErrorCode ierr; 441 PetscInt i,dim; 442 DM dm2; 443 moab::ErrorCode merr; 444 DM_Moab *dmb = (DM_Moab*)dm->data,*dd2; 445 446 PetscFunctionBegin; 447 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 448 PetscValidPointer(dmref,4); 449 450 if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) ) { 451 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); 452 if (dmb->hlevel-1 < 0 && !refine) PetscInfo1(NULL,"Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n",dmb->hlevel-1); 453 *dmref = PETSC_NULL; 454 PetscFunctionReturn(0); 455 } 456 457 ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr); 458 dd2 = (DM_Moab*)dm2->data; 459 460 dd2->mbiface = dmb->mbiface; 461 dd2->pcomm = dmb->pcomm; 462 dd2->icreatedinstance = PETSC_FALSE; 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 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 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