1 #include <petsc-private/dmmbimpl.h> /*I "petscdm.h" I*/ 2 #include <petsc-private/vecimpl.h> /*I "petscdm.h" I*/ 3 4 #include <petscdmmoab.h> 5 #include <MBTagConventions.hpp> 6 #include <moab/Skinner.hpp> 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "DMDestroy_Moab" 10 PetscErrorCode DMDestroy_Moab(DM dm) 11 { 12 PetscErrorCode ierr; 13 DM_Moab *dmmoab = (DM_Moab*)dm->data; 14 PetscSection section; 15 16 PetscFunctionBegin; 17 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 18 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 19 ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 20 if (dmmoab->icreatedinstance) { 21 delete dmmoab->mbiface; 22 } 23 dmmoab->mbiface = NULL; 24 dmmoab->pcomm = NULL; 25 delete dmmoab->vlocal; 26 delete dmmoab->vowned; 27 delete dmmoab->vghost; 28 delete dmmoab->elocal; 29 delete dmmoab->eghost; 30 delete dmmoab->bndyvtx; 31 delete dmmoab->bndyfaces; 32 delete dmmoab->bndyelems; 33 34 ierr = PetscFree(dmmoab->gsindices);CHKERRQ(ierr); 35 ierr = VecScatterDestroy(&dmmoab->ltog_sendrecv);CHKERRQ(ierr); 36 ierr = ISLocalToGlobalMappingDestroy(&dmmoab->ltog_map);CHKERRQ(ierr); 37 ierr = PetscFree(dm->data);CHKERRQ(ierr); 38 PetscFunctionReturn(0); 39 } 40 41 #undef __FUNCT__ 42 #define __FUNCT__ "DMSetUp_Moab" 43 PetscErrorCode DMSetUp_Moab(DM dm) 44 { 45 PetscErrorCode ierr; 46 moab::ErrorCode merr; 47 Vec local, global; 48 IS from; 49 moab::Range::iterator iter; 50 PetscInt i,j,bs,gsiz,lsiz; 51 DM_Moab *dmmoab = (DM_Moab*)dm->data; 52 PetscInt totsize; 53 PetscSection section; 54 PetscInt gmin,lmin,lmax; 55 56 moab::Range adj; 57 58 PetscFunctionBegin; 59 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 60 /* Get the local and shared vertices and cache it */ 61 if (dmmoab->mbiface == PETSC_NULL || dmmoab->pcomm == PETSC_NULL) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ORDER, "Set the MOAB Interface and ParallelComm objects before calling SetUp."); 62 63 /* Get the entities recursively in the current part of the mesh, if user did not set the local vertices explicitly */ 64 if (dmmoab->vlocal->empty()) { 65 merr = dmmoab->mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,*dmmoab->vlocal,true);MBERRNM(merr); 66 67 /* filter based on parallel status */ 68 merr = dmmoab->pcomm->filter_pstatus(*dmmoab->vlocal,PSTATUS_NOT_OWNED,PSTATUS_NOT,-1,dmmoab->vowned);MBERRNM(merr); 69 *dmmoab->vghost = moab::subtract(*dmmoab->vlocal, *dmmoab->vowned); 70 71 dmmoab->nloc = dmmoab->vowned->size(); 72 dmmoab->nghost = dmmoab->vghost->size(); 73 ierr = MPI_Allreduce(&dmmoab->nloc, &dmmoab->n, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr); 74 75 #if 0 76 if(dmmoab->pcomm->rank() || dmmoab->pcomm->size()==1) { 77 PetscPrintf(PETSC_COMM_SELF, "Vertices: global: %D, local: %D", dmmoab->n, dmmoab->nloc+dmmoab->nghost); 78 dmmoab->vlocal->print(0); 79 PetscPrintf(PETSC_COMM_SELF, "Vertices: owned: %D", dmmoab->nloc); 80 dmmoab->vowned->print(0); 81 PetscPrintf(PETSC_COMM_SELF, "Vertices: ghost: %D", dmmoab->nghost); 82 dmmoab->vghost->print(0); 83 } 84 #endif 85 } 86 87 /* get the information about the local elements in the mesh */ 88 { 89 dmmoab->eghost->clear(); 90 91 /* first decipher the leading dimension */ 92 for (i=3;i>0;i--) { 93 dmmoab->elocal->clear(); 94 merr = dmmoab->mbiface->get_entities_by_dimension(dmmoab->fileset, i, *dmmoab->elocal, true);CHKERRQ(merr); 95 96 /* store the current mesh dimension */ 97 if (dmmoab->elocal->size()) { 98 dmmoab->dim=i; 99 break; 100 } 101 } 102 103 *dmmoab->eghost = *dmmoab->elocal; 104 merr = dmmoab->pcomm->filter_pstatus(*dmmoab->elocal,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr); 105 *dmmoab->eghost = moab::subtract(*dmmoab->eghost, *dmmoab->elocal); 106 107 dmmoab->neleloc = dmmoab->elocal->size(); 108 ierr = MPI_Allreduce(&dmmoab->neleloc, &dmmoab->nele, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr); 109 } 110 111 bs = dmmoab->bs; 112 if (!dmmoab->ltog_tag) { 113 /* Get the global ID tag. The global ID tag is applied to each 114 vertex. It acts as an global identifier which MOAB uses to 115 assemble the individual pieces of the mesh */ 116 merr = dmmoab->mbiface->tag_get_handle(GLOBAL_ID_TAG_NAME, dmmoab->ltog_tag);MBERRNM(merr); 117 } 118 119 totsize=dmmoab->vlocal->size(); 120 ierr = PetscMalloc(totsize*sizeof(PetscInt), &dmmoab->gsindices);CHKERRQ(ierr); 121 { 122 /* first get the local indices */ 123 merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,*dmmoab->vowned,&dmmoab->gsindices[0]);MBERRNM(merr); 124 /* next get the ghosted indices */ 125 if (dmmoab->nghost) { 126 merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,*dmmoab->vghost,&dmmoab->gsindices[dmmoab->nloc]);MBERRNM(merr); 127 } 128 129 /* find out the local and global minima of GLOBAL_ID */ 130 lmin=lmax=dmmoab->gsindices[0]; 131 for (i=0; i<totsize; ++i) { 132 if(lmin>dmmoab->gsindices[i]) lmin=dmmoab->gsindices[i]; 133 if(lmax<dmmoab->gsindices[i]) lmax=dmmoab->gsindices[i]; 134 } 135 136 ierr = MPI_Allreduce(&lmin, &gmin, 1, MPI_INT, MPI_MIN, ((PetscObject)dm)->comm);CHKERRQ(ierr); 137 PetscInfo3(dm, "GLOBAL_ID: Local minima - %D, Local maxima - %D, Global minima - %D.\n", lmin, lmax, gmin); 138 } 139 140 { 141 ierr = PetscSectionCreate(((PetscObject)dm)->comm, §ion);CHKERRQ(ierr); 142 ierr = PetscSectionSetNumFields(section, dmmoab->nfields);CHKERRQ(ierr); 143 ierr = PetscSectionSetChart(section, lmin, lmax+1);CHKERRQ(ierr); 144 for (j=0; j<totsize; ++j) { 145 PetscInt locgid = dmmoab->gsindices[j]; 146 for (i=0; i < dmmoab->nfields; ++i) { 147 ierr = PetscSectionSetFieldName(section, i, dmmoab->fields[i]);CHKERRQ(ierr); 148 if (bs>1) { 149 ierr = PetscSectionSetFieldDof(section, locgid, i, (locgid-gmin)*dmmoab->nfields+i);CHKERRQ(ierr); 150 ierr = PetscSectionSetFieldOffset(section, locgid, i, (locgid-gmin)*dmmoab->nfields); 151 } 152 else { 153 ierr = PetscSectionSetFieldDof(section, locgid, i, dmmoab->n*i+locgid-gmin);CHKERRQ(ierr); 154 ierr = PetscSectionSetFieldOffset(section, locgid, i, i*dmmoab->n); 155 } 156 } 157 ierr = PetscSectionSetDof(section, locgid, dmmoab->nfields);CHKERRQ(ierr); 158 } 159 ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 160 ierr = DMSetDefaultSection(dm, section);CHKERRQ(ierr); 161 } 162 163 { 164 for (i=0; i<totsize; ++i) { 165 dmmoab->gsindices[i]-=gmin; /* zero based index needed for IS */ 166 } 167 168 /* Create Global to Local Vector Scatter Context */ 169 ierr = DMCreateGlobalVector_Moab(dm, &global);CHKERRQ(ierr); 170 ierr = DMCreateLocalVector_Moab(dm, &local);CHKERRQ(ierr); 171 172 /* global to local must retrieve ghost points */ 173 ierr = ISCreateBlock(((PetscObject)dm)->comm,bs,totsize,&dmmoab->gsindices[0],PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 174 175 ierr = VecGetLocalSize(global,&gsiz);CHKERRQ(ierr); 176 ierr = VecGetLocalSize(local,&lsiz);CHKERRQ(ierr); 177 178 ierr = VecScatterCreate(local,from,global,from,&dmmoab->ltog_sendrecv);CHKERRQ(ierr); 179 ierr = ISDestroy(&from);CHKERRQ(ierr); 180 ierr = VecDestroy(&local);CHKERRQ(ierr); 181 ierr = VecDestroy(&global);CHKERRQ(ierr); 182 } 183 184 /* skin the boundary and store nodes */ 185 { 186 /* get the skin vertices of boundary faces for the current partition and then filter 187 the local, boundary faces, vertices and elements alone via PSTATUS flags; 188 this should not give us any ghosted boundary, but if user needs such a functionality 189 it would be easy to add it based on the find_skin query below */ 190 moab::Skinner skinner(dmmoab->mbiface); 191 dmmoab->bndyvtx = new moab::Range(); 192 dmmoab->bndyfaces = new moab::Range(); 193 dmmoab->bndyelems = new moab::Range(); 194 195 /* get the entities on the skin - only the faces */ 196 merr = skinner.find_skin(dmmoab->fileset, *dmmoab->elocal, false, *dmmoab->bndyfaces, NULL, false, true, false, false);MBERRNM(merr); // 'false' param indicates we want faces back, not vertices 197 198 /* filter all the non-owned and shared entities out of the list */ 199 merr = dmmoab->pcomm->filter_pstatus(*dmmoab->bndyfaces,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr); 200 merr = dmmoab->pcomm->filter_pstatus(*dmmoab->bndyfaces,PSTATUS_SHARED,PSTATUS_NOT);MBERRNM(merr); 201 202 if (dmmoab->dim == 3) { 203 // get the edges from faces and then do the same if needed 204 } 205 206 /* get all the nodes via connectivity and the parent elements via adjacency information */ 207 merr = dmmoab->mbiface->get_connectivity(*dmmoab->bndyfaces, *dmmoab->bndyvtx, false);MBERRNM(ierr); 208 merr = dmmoab->mbiface->get_adjacencies(*dmmoab->bndyfaces, dmmoab->dim, false, *dmmoab->bndyelems, moab::Interface::UNION);MBERRNM(ierr); 209 PetscInfo3(dm, "Found %D boundary vertices, %D boundary faces and %D boundary elements.\n", dmmoab->bndyvtx->size(), dmmoab->bndyvtx->size(), dmmoab->bndyelems->size()); 210 } 211 PetscFunctionReturn(0); 212 } 213 214 215 #undef __FUNCT__ 216 #define __FUNCT__ "DMCreate_Moab" 217 PETSC_EXTERN PetscErrorCode DMCreate_Moab(DM dm) 218 { 219 PetscErrorCode ierr; 220 221 PetscFunctionBegin; 222 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 223 ierr = PetscNewLog(dm,&dm->data);CHKERRQ(ierr); 224 225 ((DM_Moab*)dm->data)->bs = 1; 226 ((DM_Moab*)dm->data)->nfields = 1; 227 ((DM_Moab*)dm->data)->n = 0; 228 ((DM_Moab*)dm->data)->nloc = 0; 229 ((DM_Moab*)dm->data)->nele = 0; 230 ((DM_Moab*)dm->data)->neleloc = 0; 231 ((DM_Moab*)dm->data)->nghost = 0; 232 ((DM_Moab*)dm->data)->ltog_map = PETSC_NULL; 233 ((DM_Moab*)dm->data)->ltog_sendrecv = PETSC_NULL; 234 235 ((DM_Moab*)dm->data)->vlocal = new moab::Range(); 236 ((DM_Moab*)dm->data)->vowned = new moab::Range(); 237 ((DM_Moab*)dm->data)->vghost = new moab::Range(); 238 ((DM_Moab*)dm->data)->elocal = new moab::Range(); 239 ((DM_Moab*)dm->data)->eghost = new moab::Range(); 240 241 dm->ops->createglobalvector = DMCreateGlobalVector_Moab; 242 dm->ops->createlocalvector = DMCreateLocalVector_Moab; 243 dm->ops->creatematrix = DMCreateMatrix_Moab; 244 dm->ops->setup = DMSetUp_Moab; 245 dm->ops->destroy = DMDestroy_Moab; 246 dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Moab; 247 dm->ops->globaltolocalend = DMGlobalToLocalEnd_Moab; 248 dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Moab; 249 dm->ops->localtoglobalend = DMLocalToGlobalEnd_Moab; 250 PetscFunctionReturn(0); 251 } 252 253 #undef __FUNCT__ 254 #define __FUNCT__ "DMMoabCreate" 255 /*@ 256 DMMoabCreate - Creates a DMMoab object, which encapsulates a moab instance 257 258 Collective on MPI_Comm 259 260 Input Parameter: 261 . comm - The communicator for the DMMoab object 262 263 Output Parameter: 264 . dmb - The DMMoab object 265 266 Level: beginner 267 268 .keywords: DMMoab, create 269 @*/ 270 PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *dmb) 271 { 272 PetscErrorCode ierr; 273 274 PetscFunctionBegin; 275 PetscValidPointer(dmb,2); 276 ierr = DMCreate(comm, dmb);CHKERRQ(ierr); 277 ierr = DMSetType(*dmb, DMMOAB);CHKERRQ(ierr); 278 PetscFunctionReturn(0); 279 } 280 281 #undef __FUNCT__ 282 #define __FUNCT__ "DMMoabCreateMoab" 283 /*@ 284 DMMoabCreate - Creates a DMMoab object, optionally from an instance and other data 285 286 Collective on MPI_Comm 287 288 Input Parameter: 289 . comm - The communicator for the DMMoab object 290 . mbiface - (ptr to) the MOAB Instance; if passed in NULL, MOAB instance is created inside PETSc, and destroyed 291 along with the DMMoab 292 . pcomm - (ptr to) a ParallelComm; if NULL, creates one internally for the whole communicator 293 . ltog_tag - A tag to use to retrieve global id for an entity; if 0, will use GLOBAL_ID_TAG_NAME/tag 294 . range - If non-NULL, contains range of entities to which DOFs will be assigned 295 296 Output Parameter: 297 . dmb - The DMMoab object 298 299 Level: intermediate 300 301 .keywords: DMMoab, create 302 @*/ 303 PetscErrorCode DMMoabCreateMoab(MPI_Comm comm, moab::Interface *mbiface, moab::ParallelComm *pcomm, moab::Tag *ltog_tag, moab::Range *range, DM *dmb) 304 { 305 PetscErrorCode ierr; 306 moab::ErrorCode merr; 307 moab::EntityHandle partnset; 308 PetscInt rank, nprocs; 309 DM_Moab *dmmoab; 310 311 PetscFunctionBegin; 312 PetscValidPointer(dmb,6); 313 ierr = DMMoabCreate(comm, dmb);CHKERRQ(ierr); 314 dmmoab = (DM_Moab*)(*dmb)->data; 315 316 if (!mbiface) { 317 dmmoab->mbiface = new moab::Core(); 318 dmmoab->icreatedinstance = PETSC_TRUE; 319 } 320 else { 321 dmmoab->mbiface = mbiface; 322 dmmoab->icreatedinstance = PETSC_FALSE; 323 } 324 325 /* create a fileset to store the hierarchy of entities belonging to current DM */ 326 merr = dmmoab->mbiface->create_meshset(moab::MESHSET_ORDERED, dmmoab->fileset);MBERR("Creating file set failed", merr); 327 328 if (!pcomm) { 329 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 330 ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr); 331 332 /* Create root sets for each mesh. Then pass these 333 to the load_file functions to be populated. */ 334 merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, partnset); 335 MBERR("Creating partition set failed", merr); 336 337 /* Create the parallel communicator object with the partition handle associated with MOAB */ 338 dmmoab->pcomm = moab::ParallelComm::get_pcomm(dmmoab->mbiface, partnset, &comm); 339 } 340 else { 341 ierr = DMMoabSetParallelComm(*dmb, pcomm);CHKERRQ(ierr); 342 } 343 344 /* do the remaining initializations for DMMoab */ 345 dmmoab->bs = 1; 346 dmmoab->nfields = 1; 347 348 /* set global ID tag handle */ 349 if (!ltog_tag) { 350 merr = dmmoab->mbiface->tag_get_handle(GLOBAL_ID_TAG_NAME, dmmoab->ltog_tag);MBERRNM(merr); 351 } 352 else { 353 ierr = DMMoabSetLocalToGlobalTag(*dmb, *ltog_tag);CHKERRQ(ierr); 354 } 355 356 /* set the local range of entities (vertices) of interest */ 357 if (range) { 358 ierr = DMMoabSetLocalVertices(*dmb, range);CHKERRQ(ierr); 359 } 360 PetscFunctionReturn(0); 361 } 362 363 #undef __FUNCT__ 364 #define __FUNCT__ "DMMoabSetParallelComm" 365 /*@ 366 DMMoabSetParallelComm - Set the ParallelComm used with this DMMoab 367 368 Collective on MPI_Comm 369 370 Input Parameter: 371 . dm - The DMMoab object being set 372 . pcomm - The ParallelComm being set on the DMMoab 373 374 Level: beginner 375 376 .keywords: DMMoab, create 377 @*/ 378 PetscErrorCode DMMoabSetParallelComm(DM dm,moab::ParallelComm *pcomm) 379 { 380 DM_Moab *dmmoab = (DM_Moab*)(dm)->data; 381 382 PetscFunctionBegin; 383 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 384 PetscValidPointer(pcomm,2); 385 dmmoab->pcomm = pcomm; 386 dmmoab->mbiface = pcomm->get_moab(); 387 dmmoab->icreatedinstance = PETSC_FALSE; 388 PetscFunctionReturn(0); 389 } 390 391 392 #undef __FUNCT__ 393 #define __FUNCT__ "DMMoabGetParallelComm" 394 /*@ 395 DMMoabGetParallelComm - Get the ParallelComm used with this DMMoab 396 397 Collective on MPI_Comm 398 399 Input Parameter: 400 . dm - The DMMoab object being set 401 402 Output Parameter: 403 . pcomm - The ParallelComm for the DMMoab 404 405 Level: beginner 406 407 .keywords: DMMoab, create 408 @*/ 409 PetscErrorCode DMMoabGetParallelComm(DM dm,moab::ParallelComm **pcomm) 410 { 411 PetscFunctionBegin; 412 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 413 *pcomm = ((DM_Moab*)(dm)->data)->pcomm; 414 PetscFunctionReturn(0); 415 } 416 417 418 #undef __FUNCT__ 419 #define __FUNCT__ "DMMoabSetInterface" 420 /*@ 421 DMMoabSetInterface - Set the MOAB instance used with this DMMoab 422 423 Collective on MPI_Comm 424 425 Input Parameter: 426 . dm - The DMMoab object being set 427 . mbiface - The MOAB instance being set on this DMMoab 428 429 Level: beginner 430 431 .keywords: DMMoab, create 432 @*/ 433 PetscErrorCode DMMoabSetInterface(DM dm,moab::Interface *mbiface) 434 { 435 DM_Moab *dmmoab = (DM_Moab*)(dm)->data; 436 437 PetscFunctionBegin; 438 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 439 PetscValidPointer(mbiface,2); 440 dmmoab->pcomm = NULL; 441 dmmoab->mbiface = mbiface; 442 dmmoab->icreatedinstance = PETSC_FALSE; 443 PetscFunctionReturn(0); 444 } 445 446 447 #undef __FUNCT__ 448 #define __FUNCT__ "DMMoabGetInterface" 449 /*@ 450 DMMoabGetInterface - Get the MOAB instance used with this DMMoab 451 452 Collective on MPI_Comm 453 454 Input Parameter: 455 . dm - The DMMoab object being set 456 457 Output Parameter: 458 . mbiface - The MOAB instance set on this DMMoab 459 460 Level: beginner 461 462 .keywords: DMMoab, create 463 @*/ 464 PetscErrorCode DMMoabGetInterface(DM dm,moab::Interface **mbiface) 465 { 466 PetscErrorCode ierr; 467 static PetscBool cite = PETSC_FALSE; 468 469 PetscFunctionBegin; 470 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 471 ierr = PetscCitationsRegister("@techreport{tautges_moab:_2004,\n type = {{SAND2004-1592}},\n title = {{MOAB:} A Mesh-Oriented Database}, institution = {Sandia National Laboratories},\n author = {Tautges, T. J. and Meyers, R. and Merkley, K. and Stimpson, C. and Ernst, C.},\n year = {2004}, note = {Report}\n}\n",&cite);CHKERRQ(ierr); 472 *mbiface = ((DM_Moab*)dm->data)->mbiface; 473 PetscFunctionReturn(0); 474 } 475 476 477 #undef __FUNCT__ 478 #define __FUNCT__ "DMMoabSetLocalVertices" 479 /*@ 480 DMMoabSetLocalVertices - Set the entities having DOFs on this DMMoab 481 482 Collective on MPI_Comm 483 484 Input Parameter: 485 . dm - The DMMoab object being set 486 . range - The entities treated by this DMMoab 487 488 Level: beginner 489 490 .keywords: DMMoab, create 491 @*/ 492 PetscErrorCode DMMoabSetLocalVertices(DM dm,moab::Range *range) 493 { 494 moab::ErrorCode merr; 495 PetscErrorCode ierr; 496 DM_Moab *dmmoab = (DM_Moab*)(dm)->data; 497 498 PetscFunctionBegin; 499 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 500 dmmoab->vlocal->clear(); 501 dmmoab->vowned->clear(); 502 dmmoab->vlocal->insert(range->begin(), range->end()); 503 *dmmoab->vowned = *dmmoab->vlocal; 504 merr = dmmoab->pcomm->filter_pstatus(*dmmoab->vowned,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr); 505 *dmmoab->vghost = moab::subtract(*range, *dmmoab->vowned); 506 dmmoab->nloc=dmmoab->vowned->size(); 507 dmmoab->nghost=dmmoab->vghost->size(); 508 ierr = MPI_Allreduce(&dmmoab->nloc, &dmmoab->n, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr); 509 PetscFunctionReturn(0); 510 } 511 512 513 #undef __FUNCT__ 514 #define __FUNCT__ "DMMoabGetLocalVertices" 515 /*@ 516 DMMoabGetLocalVertices - Get the entities having DOFs on this DMMoab 517 518 Collective on MPI_Comm 519 520 Input Parameter: 521 . dm - The DMMoab object being set 522 523 Output Parameter: 524 . owned - The owned vertex entities in this DMMoab 525 . ghost - The ghosted entities (non-owned) stored locally in this partition 526 527 Level: beginner 528 529 .keywords: DMMoab, create 530 @*/ 531 PetscErrorCode DMMoabGetLocalVertices(DM dm,moab::Range *owned,moab::Range *ghost) 532 { 533 PetscFunctionBegin; 534 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 535 if (owned) *owned = *((DM_Moab*)dm->data)->vowned; 536 if (ghost) *ghost = *((DM_Moab*)dm->data)->vghost; 537 PetscFunctionReturn(0); 538 } 539 540 #undef __FUNCT__ 541 #define __FUNCT__ "DMMoabGetLocalElements" 542 /*@ 543 DMMoabGetLocalElements - Get the higher-dimensional entities that are locally owned 544 545 Collective on MPI_Comm 546 547 Input Parameter: 548 . dm - The DMMoab object being set 549 550 Output Parameter: 551 . range - The entities owned locally 552 553 Level: beginner 554 555 .keywords: DMMoab, create 556 @*/ 557 PetscErrorCode DMMoabGetLocalElements(DM dm,moab::Range *range) 558 { 559 PetscFunctionBegin; 560 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 561 if (range) *range = *((DM_Moab*)dm->data)->elocal; 562 PetscFunctionReturn(0); 563 } 564 565 566 #undef __FUNCT__ 567 #define __FUNCT__ "DMMoabSetLocalElements" 568 /*@ 569 DMMoabSetLocalElements - Set the entities having DOFs on this DMMoab 570 571 Collective on MPI_Comm 572 573 Input Parameter: 574 . dm - The DMMoab object being set 575 . range - The entities treated by this DMMoab 576 577 Level: beginner 578 579 .keywords: DMMoab, create 580 @*/ 581 PetscErrorCode DMMoabSetLocalElements(DM dm,moab::Range *range) 582 { 583 moab::ErrorCode merr; 584 PetscErrorCode ierr; 585 DM_Moab *dmmoab = (DM_Moab*)(dm)->data; 586 587 PetscFunctionBegin; 588 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 589 dmmoab->elocal->clear(); 590 dmmoab->eghost->clear(); 591 dmmoab->elocal->insert(range->begin(), range->end()); 592 PetscInfo2(dm, "Range size = %D; elocal size = %D.\n", range->size(), dmmoab->elocal->size()); 593 merr = dmmoab->pcomm->filter_pstatus(*dmmoab->elocal,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr); 594 *dmmoab->eghost = moab::subtract(*range, *dmmoab->elocal); 595 dmmoab->neleloc=dmmoab->elocal->size(); 596 ierr = MPI_Allreduce(&dmmoab->nele, &dmmoab->neleloc, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr); 597 PetscInfo2(dm, "Created %D local and %D glocal elements.\n", dmmoab->neleloc, dmmoab->nele); 598 PetscFunctionReturn(0); 599 } 600 601 602 #undef __FUNCT__ 603 #define __FUNCT__ "DMMoabSetLocalToGlobalTag" 604 /*@ 605 DMMoabSetLocalToGlobalTag - Set the tag used for local to global numbering 606 607 Collective on MPI_Comm 608 609 Input Parameter: 610 . dm - The DMMoab object being set 611 . ltogtag - The MOAB tag used for local to global ids 612 613 Level: beginner 614 615 .keywords: DMMoab, create 616 @*/ 617 PetscErrorCode DMMoabSetLocalToGlobalTag(DM dm,moab::Tag ltogtag) 618 { 619 PetscFunctionBegin; 620 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 621 ((DM_Moab*)dm->data)->ltog_tag = ltogtag; 622 PetscFunctionReturn(0); 623 } 624 625 626 #undef __FUNCT__ 627 #define __FUNCT__ "DMMoabGetLocalToGlobalTag" 628 /*@ 629 DMMoabGetLocalToGlobalTag - Get the tag used for local to global numbering 630 631 Collective on MPI_Comm 632 633 Input Parameter: 634 . dm - The DMMoab object being set 635 636 Output Parameter: 637 . ltogtag - The MOAB tag used for local to global ids 638 639 Level: beginner 640 641 .keywords: DMMoab, create 642 @*/ 643 PetscErrorCode DMMoabGetLocalToGlobalTag(DM dm,moab::Tag *ltog_tag) 644 { 645 PetscFunctionBegin; 646 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 647 *ltog_tag = ((DM_Moab*)dm->data)->ltog_tag; 648 PetscFunctionReturn(0); 649 } 650 651 652 #undef __FUNCT__ 653 #define __FUNCT__ "DMMoabSetBlockSize" 654 /*@ 655 DMMoabSetBlockSize - Set the block size used with this DMMoab 656 657 Collective on MPI_Comm 658 659 Input Parameter: 660 . dm - The DMMoab object being set 661 . bs - The block size used with this DMMoab 662 663 Level: beginner 664 665 .keywords: DMMoab, create 666 @*/ 667 PetscErrorCode DMMoabSetBlockSize(DM dm,PetscInt bs) 668 { 669 PetscFunctionBegin; 670 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 671 ((DM_Moab*)dm->data)->bs = bs; 672 PetscFunctionReturn(0); 673 } 674 675 676 #undef __FUNCT__ 677 #define __FUNCT__ "DMMoabGetBlockSize" 678 /*@ 679 DMMoabGetBlockSize - Get the block size used with this DMMoab 680 681 Collective on MPI_Comm 682 683 Input Parameter: 684 . dm - The DMMoab object being set 685 686 Output Parameter: 687 . bs - The block size used with this DMMoab 688 689 Level: beginner 690 691 .keywords: DMMoab, create 692 @*/ 693 PetscErrorCode DMMoabGetBlockSize(DM dm,PetscInt *bs) 694 { 695 PetscFunctionBegin; 696 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 697 *bs = ((DM_Moab*)dm->data)->bs; 698 PetscFunctionReturn(0); 699 } 700 701 702 #undef __FUNCT__ 703 #define __FUNCT__ "DMMoabGetSize" 704 /*@ 705 DMMoabGetSize - Get the global vertex size used with this DMMoab 706 707 Collective on MPI_Comm 708 709 Input Parameter: 710 . dm - The DMMoab object being set 711 712 Output Parameter: 713 . ng - The global size of the DMMoab instance 714 715 Level: beginner 716 717 .keywords: DMMoab, create 718 @*/ 719 PetscErrorCode DMMoabGetSize(DM dm,PetscInt *ng) 720 { 721 PetscFunctionBegin; 722 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 723 if(ng) *ng = ((DM_Moab*)dm->data)->n; 724 PetscFunctionReturn(0); 725 } 726 727 728 #undef __FUNCT__ 729 #define __FUNCT__ "DMMoabGetLocalSize" 730 /*@ 731 DMMoabGetLocalSize - Get the local and ghosted vertex size used with this DMMoab 732 733 Collective on MPI_Comm 734 735 Input Parameter: 736 . dm - The DMMoab object being set 737 738 Output Parameter: 739 . nl - The local size of the DMMoab instance 740 . ng - The ghosted size of the DMMoab instance 741 742 Level: beginner 743 744 .keywords: DMMoab, create 745 @*/ 746 PetscErrorCode DMMoabGetLocalSize(DM dm,PetscInt *nl,PetscInt *ng) 747 { 748 PetscFunctionBegin; 749 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 750 if(nl) *nl = ((DM_Moab*)dm->data)->nloc; 751 if(ng) *ng = ((DM_Moab*)dm->data)->nghost; 752 PetscFunctionReturn(0); 753 } 754 755 756 #undef __FUNCT__ 757 #define __FUNCT__ "DMMoabGetDimension" 758 /*@ 759 DMMoabGetDimension - Get the dimension of the DM Mesh 760 761 Collective on MPI_Comm 762 763 Input Parameter: 764 . dm - The DMMoab object being set 765 766 Output Parameter: 767 . dim - The dimension of DM 768 769 Level: beginner 770 771 .keywords: DMMoab, create 772 @*/ 773 PetscErrorCode DMMoabGetDimension(DM dm,PetscInt *dim) 774 { 775 PetscFunctionBegin; 776 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 777 *dim = ((DM_Moab*)dm->data)->dim; 778 PetscFunctionReturn(0); 779 } 780 781 782 783 #undef __FUNCT__ 784 #define __FUNCT__ "DMMoabSetFieldVector" 785 PetscErrorCode DMMoabSetFieldVector(DM dm, PetscInt ifield, Vec fvec) 786 { 787 DM_Moab *dmmoab; 788 moab::Tag vtag,ntag; 789 const PetscScalar *varray; 790 PetscScalar *farray; 791 moab::ErrorCode merr; 792 PetscErrorCode ierr; 793 PetscInt doff; 794 std::string tag_name; 795 moab::Range::iterator iter; 796 797 PetscFunctionBegin; 798 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 799 dmmoab = (DM_Moab*)(dm)->data; 800 801 /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 802 merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 803 moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr); 804 805 ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr); 806 807 merr = dmmoab->mbiface->tag_get_name(vtag, tag_name); 808 if (!tag_name.length() && merr !=moab::MB_SUCCESS) { 809 ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr); 810 811 for(iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++) { 812 moab::EntityHandle vtx = (*iter); 813 814 /* get field dof index */ 815 ierr = DMMoabGetFieldDof(dm, vtx, ifield, &doff); 816 817 /* use the entity handle and the Dof index to set the right value */ 818 merr = dmmoab->mbiface->tag_set_data(ntag, &vtx, 1, (const void*)&varray[doff]);MBERRNM(merr); 819 } 820 ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr); 821 } 822 else { 823 ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&farray);CHKERRQ(ierr); 824 /* we are using a MOAB Vec - directly copy the tag data to new one */ 825 merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)farray);MBERRNM(merr); 826 merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr); 827 /* make sure the parallel exchange for ghosts are done appropriately */ 828 merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr); 829 ierr = PetscFree(farray);CHKERRQ(ierr); 830 } 831 PetscFunctionReturn(0); 832 } 833 834 835 #undef __FUNCT__ 836 #define __FUNCT__ "DMMoabSetGlobalFieldVector" 837 PetscErrorCode DMMoabSetGlobalFieldVector(DM dm, Vec fvec) 838 { 839 DM_Moab *dmmoab; 840 moab::Tag vtag,ntag; 841 const PetscScalar *varray; 842 PetscScalar *farray; 843 moab::ErrorCode merr; 844 PetscErrorCode ierr; 845 PetscSection section; 846 PetscInt i,doff,ifield; 847 std::string tag_name; 848 moab::Range::iterator iter; 849 850 PetscFunctionBegin; 851 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 852 dmmoab = (DM_Moab*)(dm)->data; 853 854 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 855 856 /* get the Tag corresponding to the global vector - possible that there is no tag associated.. */ 857 ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr); 858 merr = dmmoab->mbiface->tag_get_name(vtag, tag_name); 859 if (!tag_name.length() && merr !=moab::MB_SUCCESS) { 860 /* not a MOAB vector - use VecGetSubVector to get the parts as needed */ 861 862 ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr); 863 for (ifield=0; ifield<dmmoab->nfields; ++ifield) { 864 865 /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 866 merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 867 moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr); 868 869 for(iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++) { 870 moab::EntityHandle vtx = (*iter); 871 872 /* get field dof index */ 873 ierr = DMMoabGetFieldDof(dm, vtx, ifield, &doff); 874 875 /* use the entity handle and the Dof index to set the right value */ 876 merr = dmmoab->mbiface->tag_set_data(ntag, &vtx, 1, (const void*)&varray[doff]);MBERRNM(merr); 877 } 878 } 879 ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr); 880 } 881 else { 882 ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&farray);CHKERRQ(ierr); 883 ierr = PetscMalloc(dmmoab->nloc*dmmoab->bs*sizeof(PetscScalar),&varray);CHKERRQ(ierr); 884 885 /* we are using a MOAB Vec - directly copy the tag data to new one */ 886 merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)varray);MBERRNM(merr); 887 for (ifield=0; ifield<dmmoab->nfields; ++ifield) { 888 889 /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 890 merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 891 moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr); 892 893 /* we are using a MOAB Vec - directly copy the tag data to new one */ 894 for(i=0; i < dmmoab->nloc; i++) { 895 farray[i] = varray[i*dmmoab->bs+ifield]; 896 } 897 898 merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr); 899 /* make sure the parallel exchange for ghosts are done appropriately */ 900 merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr); 901 } 902 ierr = PetscFree(farray);CHKERRQ(ierr); 903 ierr = PetscFree(varray);CHKERRQ(ierr); 904 } 905 PetscFunctionReturn(0); 906 } 907 908 909 910 #undef __FUNCT__ 911 #define __FUNCT__ "DMMoabGetVertexCoordinates" 912 PetscErrorCode DMMoabGetVertexCoordinates(DM dm,PetscInt nconn,const moab::EntityHandle *conn,PetscScalar *vpos) 913 { 914 DM_Moab *dmmoab; 915 PetscErrorCode ierr; 916 moab::ErrorCode merr; 917 918 PetscFunctionBegin; 919 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 920 PetscValidPointer(conn,3); 921 dmmoab = (DM_Moab*)(dm)->data; 922 923 if (!vpos) { 924 ierr = PetscMalloc(sizeof(PetscScalar)*nconn*3, &vpos);CHKERRQ(ierr); 925 } 926 927 /* Get connectivity information in MOAB canonical ordering */ 928 merr = dmmoab->mbiface->get_coords(conn, nconn, vpos);MBERRNM(merr); 929 PetscFunctionReturn(0); 930 } 931 932 933 #undef __FUNCT__ 934 #define __FUNCT__ "DMMoabGetElementConnectivity" 935 PetscErrorCode DMMoabGetElementConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn,const moab::EntityHandle **conn) 936 { 937 DM_Moab *dmmoab; 938 const moab::EntityHandle *connect; 939 moab::ErrorCode merr; 940 PetscInt nnodes; 941 942 PetscFunctionBegin; 943 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 944 PetscValidPointer(conn,4); 945 dmmoab = (DM_Moab*)(dm)->data; 946 947 /* Get connectivity information in MOAB canonical ordering */ 948 merr = dmmoab->mbiface->get_connectivity(ehandle, connect, nnodes);MBERRNM(merr); 949 if (conn) *conn=connect; 950 if (nconn) *nconn=nnodes; 951 PetscFunctionReturn(0); 952 } 953 954 955 #undef __FUNCT__ 956 #define __FUNCT__ "DMMoabIsEntityOnBoundary" 957 PetscErrorCode DMMoabIsEntityOnBoundary(DM dm,const moab::EntityHandle ent,PetscBool* ent_on_boundary) 958 { 959 moab::EntityType etype; 960 DM_Moab *dmmoab; 961 PetscInt edim; 962 963 PetscFunctionBegin; 964 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 965 PetscValidPointer(ent_on_boundary,3); 966 dmmoab = (DM_Moab*)(dm)->data; 967 968 /* get the entity type and handle accordingly */ 969 etype=dmmoab->mbiface->type_from_handle(ent); 970 if(etype >= moab::MBPOLYHEDRON) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Entity type on the boundary skin is invalid. EntityType = %D\n",etype); 971 972 /* get the entity dimension */ 973 edim=dmmoab->mbiface->dimension_from_handle(ent); 974 975 *ent_on_boundary=PETSC_FALSE; 976 if(etype == moab::MBVERTEX && edim == 0) { 977 moab::Range::const_iterator giter = dmmoab->bndyvtx->find(ent); 978 if (giter != dmmoab->bndyvtx->end()) *ent_on_boundary=PETSC_TRUE; 979 } 980 else { 981 if (edim == dmmoab->dim) { /* check the higher-dimensional elements first */ 982 moab::Range::const_iterator geiter = dmmoab->bndyelems->find(ent); 983 if (geiter != dmmoab->bndyelems->end()) *ent_on_boundary=PETSC_TRUE; 984 } 985 else { /* next check the lower-dimensional faces */ 986 moab::Range::const_iterator gfiter = dmmoab->bndyfaces->find(ent); 987 if (gfiter != dmmoab->bndyfaces->end()) *ent_on_boundary=PETSC_TRUE; 988 } 989 } 990 PetscFunctionReturn(0); 991 } 992 993 994 #undef __FUNCT__ 995 #define __FUNCT__ "DMMoabCheckBoundaryVertices" 996 PetscErrorCode DMMoabCheckBoundaryVertices(DM dm,PetscInt nconn,const moab::EntityHandle *cnt,PetscBool* isbdvtx) 997 { 998 DM_Moab *dmmoab; 999 PetscInt i; 1000 1001 PetscFunctionBegin; 1002 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1003 PetscValidPointer(cnt,3); 1004 PetscValidPointer(isbdvtx,4); 1005 dmmoab = (DM_Moab*)(dm)->data; 1006 1007 for (i=0; i < nconn; ++i) { 1008 moab::Range::const_iterator giter = dmmoab->bndyvtx->find(cnt[i]); 1009 if (giter != dmmoab->bndyvtx->end()) isbdvtx[i] = PETSC_TRUE; 1010 else isbdvtx[i] = PETSC_FALSE; 1011 } 1012 PetscFunctionReturn(0); 1013 } 1014 1015 1016 #undef __FUNCT__ 1017 #define __FUNCT__ "DMMoabGetBoundaryEntities" 1018 PetscErrorCode DMMoabGetBoundaryEntities(DM dm,moab::Range *bdvtx,moab::Range* bdfaces,moab::Range* bdelems) 1019 { 1020 DM_Moab *dmmoab; 1021 1022 PetscFunctionBegin; 1023 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1024 dmmoab = (DM_Moab*)(dm)->data; 1025 1026 if (bdvtx) *bdvtx = *dmmoab->bndyvtx; 1027 if (bdfaces) *bdfaces = *dmmoab->bndyfaces; 1028 if (bdelems) *bdfaces = *dmmoab->bndyelems; 1029 PetscFunctionReturn(0); 1030 } 1031 1032 1033 #undef __FUNCT__ 1034 #define __FUNCT__ "DMMoabSetFields" 1035 PetscErrorCode DMMoabSetFields(DM dm,PetscInt nfields,const char** fields) 1036 { 1037 DM_Moab *dmmoab; 1038 1039 PetscFunctionBegin; 1040 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1041 dmmoab = (DM_Moab*)(dm)->data; 1042 1043 dmmoab->fields = fields; 1044 dmmoab->nfields = nfields; 1045 PetscFunctionReturn(0); 1046 } 1047 1048 1049 #undef __FUNCT__ 1050 #define __FUNCT__ "DMMoabGetFieldDof" 1051 PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point,PetscInt field,PetscInt* dof) 1052 { 1053 PetscSection section; 1054 PetscInt gid; 1055 PetscErrorCode ierr; 1056 moab::ErrorCode merr; 1057 DM_Moab *dmmoab; 1058 1059 PetscFunctionBegin; 1060 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1061 dmmoab = (DM_Moab*)(dm)->data; 1062 1063 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1064 1065 /* first get the global ID for the point */ 1066 merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&point,1,&gid);MBERRNM(merr); 1067 1068 /* get the dof value for the field */ 1069 ierr = PetscSectionGetFieldDof(section, gid, field, dof);CHKERRQ(ierr); 1070 PetscFunctionReturn(0); 1071 } 1072 1073 1074 #undef __FUNCT__ 1075 #define __FUNCT__ "DMMoabGetFieldDofs" 1076 PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof) 1077 { 1078 PetscInt i,gid; 1079 PetscSection section; 1080 PetscErrorCode ierr; 1081 moab::ErrorCode merr; 1082 DM_Moab *dmmoab; 1083 1084 PetscFunctionBegin; 1085 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1086 PetscValidPointer(points,2); 1087 dmmoab = (DM_Moab*)(dm)->data; 1088 1089 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1090 if (!dof) { 1091 ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr); 1092 } 1093 1094 /* first get the local indices */ 1095 merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr); 1096 1097 for (i=0; i<npoints; ++i) { 1098 gid=dof[i]; 1099 ierr = PetscSectionGetFieldDof(section, gid, field, &dof[i]);CHKERRQ(ierr); 1100 } 1101 PetscFunctionReturn(0); 1102 } 1103 1104 1105 #undef __FUNCT__ 1106 #define __FUNCT__ "DMMoabGetFieldDofsLocal" 1107 PetscErrorCode DMMoabGetFieldDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof) 1108 { 1109 PetscInt i,offset; 1110 PetscErrorCode ierr; 1111 DM_Moab *dmmoab; 1112 1113 PetscFunctionBegin; 1114 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1115 PetscValidPointer(points,2); 1116 dmmoab = (DM_Moab*)(dm)->data; 1117 1118 if (!dof) { 1119 ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr); 1120 } 1121 1122 if (dmmoab->bs > 1) { 1123 for (i=0; i<npoints; ++i) 1124 dof[i] = (points[i]-1)*dmmoab->bs+field; 1125 } 1126 else { 1127 offset = field*dmmoab->n; /* assume all fields have equal distribution */ 1128 for (i=0; i<npoints; ++i) 1129 dof[i] = offset+points[i]-1; 1130 } 1131 PetscFunctionReturn(0); 1132 } 1133 1134 1135 #undef __FUNCT__ 1136 #define __FUNCT__ "DMMoabGetDofs" 1137 PetscErrorCode DMMoabGetDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 1138 { 1139 PetscInt i,f,gid; 1140 PetscSection section; 1141 PetscErrorCode ierr; 1142 moab::ErrorCode merr; 1143 DM_Moab *dmmoab; 1144 1145 PetscFunctionBegin; 1146 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1147 PetscValidPointer(points,2); 1148 dmmoab = (DM_Moab*)(dm)->data; 1149 1150 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1151 if (!dof) { 1152 ierr = PetscMalloc(sizeof(PetscScalar)*dmmoab->nfields*npoints, &dof);CHKERRQ(ierr); 1153 } 1154 1155 /* first get the local indices */ 1156 merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr); 1157 1158 for (i=0; i<npoints; ++i) { 1159 gid=dof[i]; 1160 for (f=0; f<dmmoab->nfields; ++f) { 1161 ierr = PetscSectionGetFieldDof(section, gid, f, &dof[i*dmmoab->nfields+f]);CHKERRQ(ierr); 1162 } 1163 } 1164 PetscFunctionReturn(0); 1165 } 1166 1167 1168 #undef __FUNCT__ 1169 #define __FUNCT__ "DMMoabGetDofsLocal" 1170 PetscErrorCode DMMoabGetDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 1171 { 1172 PetscInt i,f,offset; 1173 PetscErrorCode ierr; 1174 DM_Moab *dmmoab; 1175 1176 PetscFunctionBegin; 1177 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1178 PetscValidPointer(points,2); 1179 dmmoab = (DM_Moab*)(dm)->data; 1180 1181 if (!dof) { 1182 ierr = PetscMalloc(sizeof(PetscScalar)*dmmoab->nfields*npoints, &dof);CHKERRQ(ierr); 1183 } 1184 1185 if (dmmoab->bs > 1) { 1186 for (f=0; f<dmmoab->nfields; ++f) 1187 for (i=0; i<npoints; ++i) 1188 dof[i*dmmoab->nfields+f] = (points[i]-1)*dmmoab->bs+f; 1189 } 1190 else { 1191 for (f=0; f<dmmoab->nfields; ++f) { 1192 offset = f*dmmoab->n; /* assume all fields have equal distribution - say all vertex based */ 1193 for (i=0; i<npoints; ++i) 1194 dof[i*dmmoab->nfields+f] = offset+points[i]-1; 1195 } 1196 } 1197 PetscFunctionReturn(0); 1198 } 1199 1200 1201 #undef __FUNCT__ 1202 #define __FUNCT__ "DMMoabGetDofsBlocked" 1203 PetscErrorCode DMMoabGetDofsBlocked(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 1204 { 1205 PetscInt i,gid,dofindx; 1206 PetscSection section; 1207 PetscErrorCode ierr; 1208 moab::ErrorCode merr; 1209 DM_Moab *dmmoab; 1210 1211 PetscFunctionBegin; 1212 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1213 PetscValidPointer(points,2); 1214 dmmoab = (DM_Moab*)(dm)->data; 1215 1216 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1217 if (!dof) { 1218 ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr); 1219 } 1220 1221 /* first get the local indices */ 1222 merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr); 1223 1224 for (i=0; i<npoints; ++i) { 1225 gid=dof[i]; 1226 ierr = PetscSectionGetFieldDof(section, gid, 0, &dofindx);CHKERRQ(ierr); 1227 if (dmmoab->bs > 1) dof[i]=dofindx/dmmoab->bs; 1228 else dof[i]=dofindx; 1229 } 1230 PetscFunctionReturn(0); 1231 } 1232 1233 1234 #undef __FUNCT__ 1235 #define __FUNCT__ "DMMoabGetDofsBlockedLocal" 1236 PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 1237 { 1238 PetscInt i; 1239 PetscErrorCode ierr; 1240 1241 PetscFunctionBegin; 1242 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1243 PetscValidPointer(points,2); 1244 1245 if (!dof) { 1246 ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr); 1247 } 1248 1249 for (i=0; i<npoints; ++i) 1250 dof[i] = points[i]-1; 1251 PetscFunctionReturn(0); 1252 } 1253 1254 #undef __FUNCT__ 1255 #define __FUNCT__ "DMMoab_GetWriteOptions_Private" 1256 PetscErrorCode DMMoab_GetWriteOptions_Private(PetscInt fsetid, PetscInt numproc, PetscInt dim, MoabWriteMode mode, PetscInt dbglevel, const char* extra_opts, const char** write_opts) 1257 { 1258 std::ostringstream str; 1259 1260 PetscFunctionBegin; 1261 1262 // do parallel read unless only one processor 1263 if (numproc > 1) { 1264 str << "PARALLEL=" << mode << ";"; 1265 if (fsetid>=0) str << "PARALLEL_COMM=" << fsetid << ";"; 1266 } 1267 1268 if (dbglevel) 1269 str << "CPUTIME;DEBUG_IO=" << dbglevel << ";"; 1270 1271 if (extra_opts) 1272 str << extra_opts ; 1273 1274 *write_opts = str.str().c_str(); 1275 PetscFunctionReturn(0); 1276 } 1277 1278 1279 #undef __FUNCT__ 1280 #define __FUNCT__ "DMMoabOutput" 1281 PetscErrorCode DMMoabOutput(DM dm,const char* filename,const char* usrwriteopts) 1282 { 1283 DM_Moab *dmmoab; 1284 PetscInt dbglevel=0; 1285 const char *writeopts; 1286 1287 PetscErrorCode ierr; 1288 moab::ErrorCode merr; 1289 1290 PetscFunctionBegin; 1291 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1292 dmmoab = (DM_Moab*)(dm)->data; 1293 1294 PetscBarrier((PetscObject)dm); 1295 1296 /* TODO: Use command-line options to control by_rank, verbosity, MoabReadMode and extra options */ 1297 ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for reading/writing MOAB based meshes from file", "DMMoab"); 1298 ierr = PetscOptionsInt("-dmmb_rw_dbg", "The verbosity level for reading and writing MOAB meshes", "dmmbutil.cxx", dbglevel, &dbglevel, NULL);CHKERRQ(ierr); 1299 ierr = PetscOptionsEnd(); 1300 1301 /* add mesh loading options specific to the DM */ 1302 ierr = DMMoab_GetWriteOptions_Private(dmmoab->pcomm->get_id(), dmmoab->pcomm->size(), dmmoab->dim, MOAB_PARWOPTS_WRITE_PART, dbglevel, usrwriteopts, &writeopts);CHKERRQ(ierr); 1303 PetscInfo2(dm, "Writing file %s with options: %s\n",filename,writeopts); 1304 1305 /* output file, using parallel write */ 1306 merr = dmmoab->mbiface->write_file(filename, NULL, writeopts, &dmmoab->fileset, 1);MBERRVM(dmmoab->mbiface,"Writing output of DMMoab failed.",merr); 1307 PetscFunctionReturn(0); 1308 } 1309 1310