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)->n = 0; 227 ((DM_Moab*)dm->data)->nloc = 0; 228 ((DM_Moab*)dm->data)->nele = 0; 229 ((DM_Moab*)dm->data)->neleloc = 0; 230 ((DM_Moab*)dm->data)->nghost = 0; 231 ((DM_Moab*)dm->data)->ltog_map = PETSC_NULL; 232 ((DM_Moab*)dm->data)->ltog_sendrecv = PETSC_NULL; 233 234 ((DM_Moab*)dm->data)->vlocal = new moab::Range(); 235 ((DM_Moab*)dm->data)->vowned = new moab::Range(); 236 ((DM_Moab*)dm->data)->vghost = new moab::Range(); 237 ((DM_Moab*)dm->data)->elocal = new moab::Range(); 238 ((DM_Moab*)dm->data)->eghost = new moab::Range(); 239 240 dm->ops->createglobalvector = DMCreateGlobalVector_Moab; 241 dm->ops->createlocalvector = DMCreateLocalVector_Moab; 242 dm->ops->creatematrix = DMCreateMatrix_Moab; 243 dm->ops->setup = DMSetUp_Moab; 244 dm->ops->destroy = DMDestroy_Moab; 245 dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Moab; 246 dm->ops->globaltolocalend = DMGlobalToLocalEnd_Moab; 247 dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Moab; 248 dm->ops->localtoglobalend = DMLocalToGlobalEnd_Moab; 249 PetscFunctionReturn(0); 250 } 251 252 #undef __FUNCT__ 253 #define __FUNCT__ "DMMoabCreate" 254 /*@ 255 DMMoabCreate - Creates a DMMoab object, which encapsulates a moab instance 256 257 Collective on MPI_Comm 258 259 Input Parameter: 260 . comm - The communicator for the DMMoab object 261 262 Output Parameter: 263 . dmb - The DMMoab object 264 265 Level: beginner 266 267 .keywords: DMMoab, create 268 @*/ 269 PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *dmb) 270 { 271 PetscErrorCode ierr; 272 273 PetscFunctionBegin; 274 PetscValidPointer(dmb,2); 275 ierr = DMCreate(comm, dmb);CHKERRQ(ierr); 276 ierr = DMSetType(*dmb, DMMOAB);CHKERRQ(ierr); 277 PetscFunctionReturn(0); 278 } 279 280 #undef __FUNCT__ 281 #define __FUNCT__ "DMMoabCreateMoab" 282 /*@ 283 DMMoabCreate - Creates a DMMoab object, optionally from an instance and other data 284 285 Collective on MPI_Comm 286 287 Input Parameter: 288 . comm - The communicator for the DMMoab object 289 . mbiface - (ptr to) the MOAB Instance; if passed in NULL, MOAB instance is created inside PETSc, and destroyed 290 along with the DMMoab 291 . pcomm - (ptr to) a ParallelComm; if NULL, creates one internally for the whole communicator 292 . ltog_tag - A tag to use to retrieve global id for an entity; if 0, will use GLOBAL_ID_TAG_NAME/tag 293 . range - If non-NULL, contains range of entities to which DOFs will be assigned 294 295 Output Parameter: 296 . dmb - The DMMoab object 297 298 Level: intermediate 299 300 .keywords: DMMoab, create 301 @*/ 302 PetscErrorCode DMMoabCreateMoab(MPI_Comm comm, moab::Interface *mbiface, moab::ParallelComm *pcomm, moab::Tag *ltog_tag, moab::Range *range, DM *dmb) 303 { 304 PetscErrorCode ierr; 305 moab::ErrorCode merr; 306 moab::EntityHandle partnset; 307 PetscInt rank, nprocs; 308 DM_Moab *dmmoab; 309 310 PetscFunctionBegin; 311 PetscValidPointer(dmb,6); 312 ierr = DMMoabCreate(comm, dmb);CHKERRQ(ierr); 313 dmmoab = (DM_Moab*)(*dmb)->data; 314 315 if (!mbiface) { 316 dmmoab->mbiface = new moab::Core(); 317 dmmoab->icreatedinstance = PETSC_TRUE; 318 } 319 else { 320 dmmoab->mbiface = mbiface; 321 dmmoab->icreatedinstance = PETSC_FALSE; 322 } 323 324 /* create a fileset to store the hierarchy of entities belonging to current DM */ 325 merr = dmmoab->mbiface->create_meshset(moab::MESHSET_ORDERED, dmmoab->fileset);MBERR("Creating file set failed", merr); 326 327 if (!pcomm) { 328 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 329 ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr); 330 331 /* Create root sets for each mesh. Then pass these 332 to the load_file functions to be populated. */ 333 merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, partnset); 334 MBERR("Creating partition set failed", merr); 335 336 /* Create the parallel communicator object with the partition handle associated with MOAB */ 337 dmmoab->pcomm = moab::ParallelComm::get_pcomm(dmmoab->mbiface, partnset, &comm); 338 } 339 else { 340 ierr = DMMoabSetParallelComm(*dmb, pcomm);CHKERRQ(ierr); 341 } 342 343 /* do the remaining initializations for DMMoab */ 344 dmmoab->bs = 1; 345 346 /* set global ID tag handle */ 347 if (!ltog_tag) { 348 merr = dmmoab->mbiface->tag_get_handle(GLOBAL_ID_TAG_NAME, dmmoab->ltog_tag);MBERRNM(merr); 349 } 350 else { 351 ierr = DMMoabSetLocalToGlobalTag(*dmb, *ltog_tag);CHKERRQ(ierr); 352 } 353 354 /* set the local range of entities (vertices) of interest */ 355 if (range) { 356 ierr = DMMoabSetLocalVertices(*dmb, range);CHKERRQ(ierr); 357 } 358 PetscFunctionReturn(0); 359 } 360 361 #undef __FUNCT__ 362 #define __FUNCT__ "DMMoabSetParallelComm" 363 /*@ 364 DMMoabSetParallelComm - Set the ParallelComm used with this DMMoab 365 366 Collective on MPI_Comm 367 368 Input Parameter: 369 . dm - The DMMoab object being set 370 . pcomm - The ParallelComm being set on the DMMoab 371 372 Level: beginner 373 374 .keywords: DMMoab, create 375 @*/ 376 PetscErrorCode DMMoabSetParallelComm(DM dm,moab::ParallelComm *pcomm) 377 { 378 DM_Moab *dmmoab = (DM_Moab*)(dm)->data; 379 380 PetscFunctionBegin; 381 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 382 PetscValidPointer(pcomm,2); 383 dmmoab->pcomm = pcomm; 384 dmmoab->mbiface = pcomm->get_moab(); 385 dmmoab->icreatedinstance = PETSC_FALSE; 386 PetscFunctionReturn(0); 387 } 388 389 390 #undef __FUNCT__ 391 #define __FUNCT__ "DMMoabGetParallelComm" 392 /*@ 393 DMMoabGetParallelComm - Get the ParallelComm used with this DMMoab 394 395 Collective on MPI_Comm 396 397 Input Parameter: 398 . dm - The DMMoab object being set 399 400 Output Parameter: 401 . pcomm - The ParallelComm for the DMMoab 402 403 Level: beginner 404 405 .keywords: DMMoab, create 406 @*/ 407 PetscErrorCode DMMoabGetParallelComm(DM dm,moab::ParallelComm **pcomm) 408 { 409 PetscFunctionBegin; 410 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 411 *pcomm = ((DM_Moab*)(dm)->data)->pcomm; 412 PetscFunctionReturn(0); 413 } 414 415 416 #undef __FUNCT__ 417 #define __FUNCT__ "DMMoabSetInterface" 418 /*@ 419 DMMoabSetInterface - Set the MOAB instance used with this DMMoab 420 421 Collective on MPI_Comm 422 423 Input Parameter: 424 . dm - The DMMoab object being set 425 . mbiface - The MOAB instance being set on this DMMoab 426 427 Level: beginner 428 429 .keywords: DMMoab, create 430 @*/ 431 PetscErrorCode DMMoabSetInterface(DM dm,moab::Interface *mbiface) 432 { 433 DM_Moab *dmmoab = (DM_Moab*)(dm)->data; 434 435 PetscFunctionBegin; 436 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 437 PetscValidPointer(mbiface,2); 438 dmmoab->pcomm = NULL; 439 dmmoab->mbiface = mbiface; 440 dmmoab->icreatedinstance = PETSC_FALSE; 441 PetscFunctionReturn(0); 442 } 443 444 445 #undef __FUNCT__ 446 #define __FUNCT__ "DMMoabGetInterface" 447 /*@ 448 DMMoabGetInterface - Get the MOAB instance used with this DMMoab 449 450 Collective on MPI_Comm 451 452 Input Parameter: 453 . dm - The DMMoab object being set 454 455 Output Parameter: 456 . mbiface - The MOAB instance set on this DMMoab 457 458 Level: beginner 459 460 .keywords: DMMoab, create 461 @*/ 462 PetscErrorCode DMMoabGetInterface(DM dm,moab::Interface **mbiface) 463 { 464 PetscErrorCode ierr; 465 static PetscBool cite = PETSC_FALSE; 466 467 PetscFunctionBegin; 468 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 469 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); 470 *mbiface = ((DM_Moab*)dm->data)->mbiface; 471 PetscFunctionReturn(0); 472 } 473 474 475 #undef __FUNCT__ 476 #define __FUNCT__ "DMMoabSetLocalVertices" 477 /*@ 478 DMMoabSetLocalVertices - Set the entities having DOFs on this DMMoab 479 480 Collective on MPI_Comm 481 482 Input Parameter: 483 . dm - The DMMoab object being set 484 . range - The entities treated by this DMMoab 485 486 Level: beginner 487 488 .keywords: DMMoab, create 489 @*/ 490 PetscErrorCode DMMoabSetLocalVertices(DM dm,moab::Range *range) 491 { 492 moab::ErrorCode merr; 493 PetscErrorCode ierr; 494 DM_Moab *dmmoab = (DM_Moab*)(dm)->data; 495 496 PetscFunctionBegin; 497 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 498 dmmoab->vlocal->clear(); 499 dmmoab->vowned->clear(); 500 dmmoab->vlocal->insert(range->begin(), range->end()); 501 *dmmoab->vowned = *dmmoab->vlocal; 502 merr = dmmoab->pcomm->filter_pstatus(*dmmoab->vowned,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr); 503 *dmmoab->vghost = moab::subtract(*range, *dmmoab->vowned); 504 dmmoab->nloc=dmmoab->vowned->size(); 505 dmmoab->nghost=dmmoab->vghost->size(); 506 ierr = MPI_Allreduce(&dmmoab->nloc, &dmmoab->n, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr); 507 PetscFunctionReturn(0); 508 } 509 510 511 #undef __FUNCT__ 512 #define __FUNCT__ "DMMoabGetLocalVertices" 513 /*@ 514 DMMoabGetLocalVertices - Get the entities having DOFs on this DMMoab 515 516 Collective on MPI_Comm 517 518 Input Parameter: 519 . dm - The DMMoab object being set 520 521 Output Parameter: 522 . owned - The owned vertex entities in this DMMoab 523 . ghost - The ghosted entities (non-owned) stored locally in this partition 524 525 Level: beginner 526 527 .keywords: DMMoab, create 528 @*/ 529 PetscErrorCode DMMoabGetLocalVertices(DM dm,moab::Range *owned,moab::Range *ghost) 530 { 531 PetscFunctionBegin; 532 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 533 if (owned) *owned = *((DM_Moab*)dm->data)->vowned; 534 if (ghost) *ghost = *((DM_Moab*)dm->data)->vghost; 535 PetscFunctionReturn(0); 536 } 537 538 #undef __FUNCT__ 539 #define __FUNCT__ "DMMoabGetLocalElements" 540 /*@ 541 DMMoabGetLocalElements - Get the higher-dimensional entities that are locally owned 542 543 Collective on MPI_Comm 544 545 Input Parameter: 546 . dm - The DMMoab object being set 547 548 Output Parameter: 549 . range - The entities owned locally 550 551 Level: beginner 552 553 .keywords: DMMoab, create 554 @*/ 555 PetscErrorCode DMMoabGetLocalElements(DM dm,moab::Range *range) 556 { 557 PetscFunctionBegin; 558 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 559 if (range) *range = *((DM_Moab*)dm->data)->elocal; 560 PetscFunctionReturn(0); 561 } 562 563 564 #undef __FUNCT__ 565 #define __FUNCT__ "DMMoabSetLocalElements" 566 /*@ 567 DMMoabSetLocalElements - Set the entities having DOFs on this DMMoab 568 569 Collective on MPI_Comm 570 571 Input Parameter: 572 . dm - The DMMoab object being set 573 . range - The entities treated by this DMMoab 574 575 Level: beginner 576 577 .keywords: DMMoab, create 578 @*/ 579 PetscErrorCode DMMoabSetLocalElements(DM dm,moab::Range *range) 580 { 581 moab::ErrorCode merr; 582 PetscErrorCode ierr; 583 DM_Moab *dmmoab = (DM_Moab*)(dm)->data; 584 585 PetscFunctionBegin; 586 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 587 dmmoab->elocal->clear(); 588 dmmoab->eghost->clear(); 589 dmmoab->elocal->insert(range->begin(), range->end()); 590 PetscInfo2(dm, "Range size = %D; elocal size = %D.\n", range->size(), dmmoab->elocal->size()); 591 merr = dmmoab->pcomm->filter_pstatus(*dmmoab->elocal,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr); 592 *dmmoab->eghost = moab::subtract(*range, *dmmoab->elocal); 593 dmmoab->neleloc=dmmoab->elocal->size(); 594 ierr = MPI_Allreduce(&dmmoab->nele, &dmmoab->neleloc, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr); 595 PetscInfo2(dm, "Created %D local and %D glocal elements.\n", dmmoab->neleloc, dmmoab->nele); 596 PetscFunctionReturn(0); 597 } 598 599 600 #undef __FUNCT__ 601 #define __FUNCT__ "DMMoabSetLocalToGlobalTag" 602 /*@ 603 DMMoabSetLocalToGlobalTag - Set the tag used for local to global numbering 604 605 Collective on MPI_Comm 606 607 Input Parameter: 608 . dm - The DMMoab object being set 609 . ltogtag - The MOAB tag used for local to global ids 610 611 Level: beginner 612 613 .keywords: DMMoab, create 614 @*/ 615 PetscErrorCode DMMoabSetLocalToGlobalTag(DM dm,moab::Tag ltogtag) 616 { 617 PetscFunctionBegin; 618 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 619 ((DM_Moab*)dm->data)->ltog_tag = ltogtag; 620 PetscFunctionReturn(0); 621 } 622 623 624 #undef __FUNCT__ 625 #define __FUNCT__ "DMMoabGetLocalToGlobalTag" 626 /*@ 627 DMMoabGetLocalToGlobalTag - Get the tag used for local to global numbering 628 629 Collective on MPI_Comm 630 631 Input Parameter: 632 . dm - The DMMoab object being set 633 634 Output Parameter: 635 . ltogtag - The MOAB tag used for local to global ids 636 637 Level: beginner 638 639 .keywords: DMMoab, create 640 @*/ 641 PetscErrorCode DMMoabGetLocalToGlobalTag(DM dm,moab::Tag *ltog_tag) 642 { 643 PetscFunctionBegin; 644 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 645 *ltog_tag = ((DM_Moab*)dm->data)->ltog_tag; 646 PetscFunctionReturn(0); 647 } 648 649 650 #undef __FUNCT__ 651 #define __FUNCT__ "DMMoabSetBlockSize" 652 /*@ 653 DMMoabSetBlockSize - Set the block size used with this DMMoab 654 655 Collective on MPI_Comm 656 657 Input Parameter: 658 . dm - The DMMoab object being set 659 . bs - The block size used with this DMMoab 660 661 Level: beginner 662 663 .keywords: DMMoab, create 664 @*/ 665 PetscErrorCode DMMoabSetBlockSize(DM dm,PetscInt bs) 666 { 667 PetscFunctionBegin; 668 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 669 ((DM_Moab*)dm->data)->bs = bs; 670 PetscFunctionReturn(0); 671 } 672 673 674 #undef __FUNCT__ 675 #define __FUNCT__ "DMMoabGetBlockSize" 676 /*@ 677 DMMoabGetBlockSize - Get the block size used with this DMMoab 678 679 Collective on MPI_Comm 680 681 Input Parameter: 682 . dm - The DMMoab object being set 683 684 Output Parameter: 685 . bs - The block size used with this DMMoab 686 687 Level: beginner 688 689 .keywords: DMMoab, create 690 @*/ 691 PetscErrorCode DMMoabGetBlockSize(DM dm,PetscInt *bs) 692 { 693 PetscFunctionBegin; 694 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 695 *bs = ((DM_Moab*)dm->data)->bs; 696 PetscFunctionReturn(0); 697 } 698 699 700 #undef __FUNCT__ 701 #define __FUNCT__ "DMMoabGetSize" 702 /*@ 703 DMMoabGetSize - Get the global vertex size used with this DMMoab 704 705 Collective on MPI_Comm 706 707 Input Parameter: 708 . dm - The DMMoab object being set 709 710 Output Parameter: 711 . ng - The global size of the DMMoab instance 712 713 Level: beginner 714 715 .keywords: DMMoab, create 716 @*/ 717 PetscErrorCode DMMoabGetSize(DM dm,PetscInt *ng) 718 { 719 PetscFunctionBegin; 720 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 721 if(ng) *ng = ((DM_Moab*)dm->data)->n; 722 PetscFunctionReturn(0); 723 } 724 725 726 #undef __FUNCT__ 727 #define __FUNCT__ "DMMoabGetLocalSize" 728 /*@ 729 DMMoabGetLocalSize - Get the local and ghosted vertex size used with this DMMoab 730 731 Collective on MPI_Comm 732 733 Input Parameter: 734 . dm - The DMMoab object being set 735 736 Output Parameter: 737 . nl - The local size of the DMMoab instance 738 . ng - The ghosted size of the DMMoab instance 739 740 Level: beginner 741 742 .keywords: DMMoab, create 743 @*/ 744 PetscErrorCode DMMoabGetLocalSize(DM dm,PetscInt *nl,PetscInt *ng) 745 { 746 PetscFunctionBegin; 747 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 748 if(nl) *nl = ((DM_Moab*)dm->data)->nloc; 749 if(ng) *ng = ((DM_Moab*)dm->data)->nghost; 750 PetscFunctionReturn(0); 751 } 752 753 754 #undef __FUNCT__ 755 #define __FUNCT__ "DMMoabSetFieldVector" 756 PetscErrorCode DMMoabSetFieldVector(DM dm, PetscInt ifield, Vec fvec) 757 { 758 DM_Moab *dmmoab; 759 moab::Tag vtag,ntag; 760 const PetscScalar *varray; 761 PetscScalar *farray; 762 moab::ErrorCode merr; 763 PetscErrorCode ierr; 764 PetscInt doff; 765 std::string tag_name; 766 moab::Range::iterator iter; 767 768 PetscFunctionBegin; 769 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 770 dmmoab = (DM_Moab*)(dm)->data; 771 772 /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 773 merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 774 moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr); 775 776 ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr); 777 778 merr = dmmoab->mbiface->tag_get_name(vtag, tag_name); 779 if (!tag_name.length() && merr !=moab::MB_SUCCESS) { 780 ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr); 781 782 for(iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++) { 783 moab::EntityHandle vtx = (*iter); 784 785 /* get field dof index */ 786 ierr = DMMoabGetFieldDof(dm, vtx, ifield, &doff); 787 788 /* use the entity handle and the Dof index to set the right value */ 789 merr = dmmoab->mbiface->tag_set_data(ntag, &vtx, 1, (const void*)&varray[doff]);MBERRNM(merr); 790 } 791 ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr); 792 } 793 else { 794 ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&farray);CHKERRQ(ierr); 795 /* we are using a MOAB Vec - directly copy the tag data to new one */ 796 merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)farray);MBERRNM(merr); 797 merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr); 798 /* make sure the parallel exchange for ghosts are done appropriately */ 799 merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr); 800 ierr = PetscFree(farray);CHKERRQ(ierr); 801 } 802 PetscFunctionReturn(0); 803 } 804 805 806 #undef __FUNCT__ 807 #define __FUNCT__ "DMMoabSetGlobalFieldVector" 808 PetscErrorCode DMMoabSetGlobalFieldVector(DM dm, Vec fvec) 809 { 810 DM_Moab *dmmoab; 811 moab::Tag vtag,ntag; 812 const PetscScalar *varray; 813 PetscScalar *farray; 814 moab::ErrorCode merr; 815 PetscErrorCode ierr; 816 PetscSection section; 817 PetscInt i,doff,ifield; 818 std::string tag_name; 819 moab::Range::iterator iter; 820 821 PetscFunctionBegin; 822 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 823 dmmoab = (DM_Moab*)(dm)->data; 824 825 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 826 827 /* get the Tag corresponding to the global vector - possible that there is no tag associated.. */ 828 ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr); 829 merr = dmmoab->mbiface->tag_get_name(vtag, tag_name); 830 if (!tag_name.length() && merr !=moab::MB_SUCCESS) { 831 /* not a MOAB vector - use VecGetSubVector to get the parts as needed */ 832 833 ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr); 834 for (ifield=0; ifield<dmmoab->nfields; ++ifield) { 835 836 /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 837 merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 838 moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr); 839 840 for(iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++) { 841 moab::EntityHandle vtx = (*iter); 842 843 /* get field dof index */ 844 ierr = DMMoabGetFieldDof(dm, vtx, ifield, &doff); 845 846 /* use the entity handle and the Dof index to set the right value */ 847 merr = dmmoab->mbiface->tag_set_data(ntag, &vtx, 1, (const void*)&varray[doff]);MBERRNM(merr); 848 } 849 } 850 ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr); 851 } 852 else { 853 ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&farray);CHKERRQ(ierr); 854 ierr = PetscMalloc(dmmoab->nloc*dmmoab->bs*sizeof(PetscScalar),&varray);CHKERRQ(ierr); 855 856 /* we are using a MOAB Vec - directly copy the tag data to new one */ 857 merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)varray);MBERRNM(merr); 858 for (ifield=0; ifield<dmmoab->nfields; ++ifield) { 859 860 /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 861 merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 862 moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr); 863 864 /* we are using a MOAB Vec - directly copy the tag data to new one */ 865 for(i=0; i < dmmoab->nloc; i++) { 866 farray[i] = varray[i*dmmoab->bs+ifield]; 867 } 868 869 merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr); 870 /* make sure the parallel exchange for ghosts are done appropriately */ 871 merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr); 872 } 873 ierr = PetscFree(farray);CHKERRQ(ierr); 874 ierr = PetscFree(varray);CHKERRQ(ierr); 875 } 876 PetscFunctionReturn(0); 877 } 878 879 880 881 #undef __FUNCT__ 882 #define __FUNCT__ "DMMoabGetVertexCoordinates" 883 PetscErrorCode DMMoabGetVertexCoordinates(DM dm,PetscInt nconn,const moab::EntityHandle *conn,PetscScalar *vpos) 884 { 885 DM_Moab *dmmoab; 886 PetscErrorCode ierr; 887 moab::ErrorCode merr; 888 889 PetscFunctionBegin; 890 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 891 PetscValidPointer(conn,3); 892 dmmoab = (DM_Moab*)(dm)->data; 893 894 if (!vpos) { 895 ierr = PetscMalloc(sizeof(PetscScalar)*nconn*3, &vpos);CHKERRQ(ierr); 896 } 897 898 /* Get connectivity information in MOAB canonical ordering */ 899 merr = dmmoab->mbiface->get_coords(conn, nconn, vpos);MBERRNM(merr); 900 PetscFunctionReturn(0); 901 } 902 903 904 #undef __FUNCT__ 905 #define __FUNCT__ "DMMoabGetElementConnectivity" 906 PetscErrorCode DMMoabGetElementConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn,const moab::EntityHandle **conn) 907 { 908 DM_Moab *dmmoab; 909 const moab::EntityHandle *connect; 910 moab::ErrorCode merr; 911 PetscInt nnodes; 912 913 PetscFunctionBegin; 914 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 915 PetscValidPointer(conn,4); 916 dmmoab = (DM_Moab*)(dm)->data; 917 918 /* Get connectivity information in MOAB canonical ordering */ 919 merr = dmmoab->mbiface->get_connectivity(ehandle, connect, nnodes);MBERRNM(merr); 920 if (conn) *conn=connect; 921 if (nconn) *nconn=nnodes; 922 PetscFunctionReturn(0); 923 } 924 925 926 #undef __FUNCT__ 927 #define __FUNCT__ "DMMoabIsEntityOnBoundary" 928 PetscErrorCode DMMoabIsEntityOnBoundary(DM dm,const moab::EntityHandle ent,PetscBool* ent_on_boundary) 929 { 930 moab::EntityType etype; 931 DM_Moab *dmmoab; 932 PetscInt edim; 933 934 PetscFunctionBegin; 935 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 936 PetscValidPointer(ent_on_boundary,3); 937 dmmoab = (DM_Moab*)(dm)->data; 938 939 /* get the entity type and handle accordingly */ 940 etype=dmmoab->mbiface->type_from_handle(ent); 941 if(etype >= moab::MBPOLYHEDRON) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Entity type on the boundary skin is invalid. EntityType = %D\n",etype); 942 943 /* get the entity dimension */ 944 edim=dmmoab->mbiface->dimension_from_handle(ent); 945 946 *ent_on_boundary=PETSC_FALSE; 947 if(etype == moab::MBVERTEX && edim == 0) { 948 moab::Range::const_iterator giter = dmmoab->bndyvtx->find(ent); 949 if (giter != dmmoab->bndyvtx->end()) *ent_on_boundary=PETSC_TRUE; 950 } 951 else { 952 if (edim == dmmoab->dim) { /* check the higher-dimensional elements first */ 953 moab::Range::const_iterator geiter = dmmoab->bndyelems->find(ent); 954 if (geiter != dmmoab->bndyelems->end()) *ent_on_boundary=PETSC_TRUE; 955 } 956 else { /* next check the lower-dimensional faces */ 957 moab::Range::const_iterator gfiter = dmmoab->bndyfaces->find(ent); 958 if (gfiter != dmmoab->bndyfaces->end()) *ent_on_boundary=PETSC_TRUE; 959 } 960 } 961 PetscFunctionReturn(0); 962 } 963 964 965 #undef __FUNCT__ 966 #define __FUNCT__ "DMMoabCheckBoundaryVertices" 967 PetscErrorCode DMMoabCheckBoundaryVertices(DM dm,PetscInt nconn,const moab::EntityHandle *cnt,PetscBool* isbdvtx) 968 { 969 DM_Moab *dmmoab; 970 PetscInt i; 971 972 PetscFunctionBegin; 973 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 974 PetscValidPointer(cnt,3); 975 PetscValidPointer(isbdvtx,4); 976 dmmoab = (DM_Moab*)(dm)->data; 977 978 for (i=0; i < nconn; ++i) { 979 moab::Range::const_iterator giter = dmmoab->bndyvtx->find(cnt[i]); 980 if (giter != dmmoab->bndyvtx->end()) isbdvtx[i] = PETSC_TRUE; 981 else isbdvtx[i] = PETSC_FALSE; 982 } 983 PetscFunctionReturn(0); 984 } 985 986 987 #undef __FUNCT__ 988 #define __FUNCT__ "DMMoabGetBoundaryEntities" 989 PetscErrorCode DMMoabGetBoundaryEntities(DM dm,moab::Range *bdvtx,moab::Range* bdfaces,moab::Range* bdelems) 990 { 991 DM_Moab *dmmoab; 992 993 PetscFunctionBegin; 994 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 995 dmmoab = (DM_Moab*)(dm)->data; 996 997 if (bdvtx) *bdvtx = *dmmoab->bndyvtx; 998 if (bdfaces) *bdfaces = *dmmoab->bndyfaces; 999 if (bdelems) *bdfaces = *dmmoab->bndyelems; 1000 PetscFunctionReturn(0); 1001 } 1002 1003 1004 #undef __FUNCT__ 1005 #define __FUNCT__ "DMMoabSetFields" 1006 PetscErrorCode DMMoabSetFields(DM dm,PetscInt nfields,const char** fields) 1007 { 1008 DM_Moab *dmmoab; 1009 1010 PetscFunctionBegin; 1011 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1012 dmmoab = (DM_Moab*)(dm)->data; 1013 1014 dmmoab->fields = fields; 1015 dmmoab->nfields = nfields; 1016 PetscFunctionReturn(0); 1017 } 1018 1019 1020 #undef __FUNCT__ 1021 #define __FUNCT__ "DMMoabGetFieldDof" 1022 PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point,PetscInt field,PetscInt* dof) 1023 { 1024 PetscSection section; 1025 PetscInt gid; 1026 PetscErrorCode ierr; 1027 moab::ErrorCode merr; 1028 DM_Moab *dmmoab; 1029 1030 PetscFunctionBegin; 1031 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1032 dmmoab = (DM_Moab*)(dm)->data; 1033 1034 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1035 1036 /* first get the global ID for the point */ 1037 merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&point,1,&gid);MBERRNM(merr); 1038 1039 /* get the dof value for the field */ 1040 ierr = PetscSectionGetFieldDof(section, gid, field, dof);CHKERRQ(ierr); 1041 PetscFunctionReturn(0); 1042 } 1043 1044 1045 #undef __FUNCT__ 1046 #define __FUNCT__ "DMMoabGetFieldDofs" 1047 PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof) 1048 { 1049 PetscInt i,gid; 1050 PetscSection section; 1051 PetscErrorCode ierr; 1052 moab::ErrorCode merr; 1053 DM_Moab *dmmoab; 1054 1055 PetscFunctionBegin; 1056 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1057 dmmoab = (DM_Moab*)(dm)->data; 1058 1059 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1060 if (!dof) { 1061 ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr); 1062 } 1063 1064 /* first get the local indices */ 1065 merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr); 1066 1067 for (i=0; i<npoints; ++i) { 1068 gid=dof[i]; 1069 ierr = PetscSectionGetFieldDof(section, gid, field, &dof[i]);CHKERRQ(ierr); 1070 } 1071 PetscFunctionReturn(0); 1072 } 1073 1074 1075 #undef __FUNCT__ 1076 #define __FUNCT__ "DMMoabGetFieldDofsLocal" 1077 PetscErrorCode DMMoabGetFieldDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof) 1078 { 1079 PetscInt i,offset; 1080 PetscErrorCode ierr; 1081 DM_Moab *dmmoab; 1082 1083 PetscFunctionBegin; 1084 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1085 dmmoab = (DM_Moab*)(dm)->data; 1086 1087 if (!dof) { 1088 ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr); 1089 } 1090 1091 if (dmmoab->bs > 1) { 1092 for (i=0; i<npoints; ++i) 1093 dof[i] = (points[i]-1)*dmmoab->bs+field; 1094 } 1095 else { 1096 offset = field*dmmoab->n; /* assume all fields have equal distribution */ 1097 for (i=0; i<npoints; ++i) 1098 dof[i] = offset+points[i]-1; 1099 } 1100 PetscFunctionReturn(0); 1101 } 1102 1103 1104 #undef __FUNCT__ 1105 #define __FUNCT__ "DMMoabGetDofs" 1106 PetscErrorCode DMMoabGetDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 1107 { 1108 PetscInt i,f,gid; 1109 PetscSection section; 1110 PetscErrorCode ierr; 1111 moab::ErrorCode merr; 1112 DM_Moab *dmmoab; 1113 1114 PetscFunctionBegin; 1115 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1116 dmmoab = (DM_Moab*)(dm)->data; 1117 1118 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1119 if (!dof) { 1120 ierr = PetscMalloc(sizeof(PetscScalar)*dmmoab->nfields*npoints, &dof);CHKERRQ(ierr); 1121 } 1122 1123 /* first get the local indices */ 1124 merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr); 1125 1126 for (i=0; i<npoints; ++i) { 1127 gid=dof[i]; 1128 for (f=0; f<dmmoab->nfields; ++f) { 1129 ierr = PetscSectionGetFieldDof(section, gid, f, &dof[i*dmmoab->nfields+f]);CHKERRQ(ierr); 1130 } 1131 } 1132 PetscFunctionReturn(0); 1133 } 1134 1135 1136 #undef __FUNCT__ 1137 #define __FUNCT__ "DMMoabGetDofsLocal" 1138 PetscErrorCode DMMoabGetDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 1139 { 1140 PetscInt i,f,offset; 1141 PetscErrorCode ierr; 1142 DM_Moab *dmmoab; 1143 1144 PetscFunctionBegin; 1145 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1146 dmmoab = (DM_Moab*)(dm)->data; 1147 1148 if (!dof) { 1149 ierr = PetscMalloc(sizeof(PetscScalar)*dmmoab->nfields*npoints, &dof);CHKERRQ(ierr); 1150 } 1151 1152 if (dmmoab->bs > 1) { 1153 for (f=0; f<dmmoab->nfields; ++f) 1154 for (i=0; i<npoints; ++i) 1155 dof[i*dmmoab->nfields+f] = (points[i]-1)*dmmoab->bs+f; 1156 } 1157 else { 1158 for (f=0; f<dmmoab->nfields; ++f) { 1159 offset = f*dmmoab->n; /* assume all fields have equal distribution - say all vertex based */ 1160 for (i=0; i<npoints; ++i) 1161 dof[i*dmmoab->nfields+f] = offset+points[i]-1; 1162 } 1163 } 1164 PetscFunctionReturn(0); 1165 } 1166 1167 1168 #undef __FUNCT__ 1169 #define __FUNCT__ "DMMoabGetDofsBlocked" 1170 PetscErrorCode DMMoabGetDofsBlocked(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 1171 { 1172 PetscInt i,gid,dofindx; 1173 PetscSection section; 1174 PetscErrorCode ierr; 1175 moab::ErrorCode merr; 1176 DM_Moab *dmmoab; 1177 1178 PetscFunctionBegin; 1179 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1180 dmmoab = (DM_Moab*)(dm)->data; 1181 1182 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1183 if (!dof) { 1184 ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr); 1185 } 1186 1187 /* first get the local indices */ 1188 merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr); 1189 1190 for (i=0; i<npoints; ++i) { 1191 gid=dof[i]; 1192 ierr = PetscSectionGetFieldDof(section, gid, 0, &dofindx);CHKERRQ(ierr); 1193 if (dmmoab->bs > 1) dof[i]=dofindx/dmmoab->bs; 1194 else dof[i]=dofindx; 1195 // PetscPrintf(PETSC_COMM_SELF, "I=%D, GID=%D, DOFINDX=%D, DOF=%D\n", i, gid, dofindx, dof[i]); 1196 } 1197 PetscFunctionReturn(0); 1198 } 1199 1200 1201 #undef __FUNCT__ 1202 #define __FUNCT__ "DMMoabGetDofsBlockedLocal" 1203 PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 1204 { 1205 PetscInt i; 1206 PetscErrorCode ierr; 1207 1208 PetscFunctionBegin; 1209 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1210 1211 if (!dof) { 1212 ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr); 1213 } 1214 1215 for (i=0; i<npoints; ++i) 1216 dof[i] = points[i]-1; 1217 PetscFunctionReturn(0); 1218 } 1219 1220 #undef __FUNCT__ 1221 #define __FUNCT__ "DMMoab_GetWriteOptions_Private" 1222 PetscErrorCode DMMoab_GetWriteOptions_Private(PetscInt fsetid, PetscInt numproc, PetscInt dim, MoabWriteMode mode, PetscInt dbglevel, const char* extra_opts, const char** write_opts) 1223 { 1224 std::ostringstream str; 1225 1226 PetscFunctionBegin; 1227 1228 // do parallel read unless only one processor 1229 if (numproc > 1) { 1230 str << "PARALLEL=" << mode << ";"; 1231 if (fsetid>=0) str << "PARALLEL_COMM=" << fsetid << ";"; 1232 } 1233 1234 if (dbglevel) 1235 str << "CPUTIME;DEBUG_IO=" << dbglevel << ";"; 1236 1237 if (extra_opts) 1238 str << extra_opts ; 1239 1240 *write_opts = str.str().c_str(); 1241 PetscFunctionReturn(0); 1242 } 1243 1244 1245 #undef __FUNCT__ 1246 #define __FUNCT__ "DMMoabOutput" 1247 PetscErrorCode DMMoabOutput(DM dm,const char* filename,const char* usrwriteopts) 1248 { 1249 DM_Moab *dmmoab; 1250 PetscInt dbglevel=0; 1251 const char *writeopts; 1252 1253 PetscErrorCode ierr; 1254 moab::ErrorCode merr; 1255 1256 PetscFunctionBegin; 1257 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1258 dmmoab = (DM_Moab*)(dm)->data; 1259 1260 PetscBarrier((PetscObject)dm); 1261 1262 /* TODO: Use command-line options to control by_rank, verbosity, MoabReadMode and extra options */ 1263 ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for reading/writing MOAB based meshes from file", "DMMoab"); 1264 ierr = PetscOptionsInt("-dmmb_rw_dbg", "The verbosity level for reading and writing MOAB meshes", "dmmbutil.cxx", dbglevel, &dbglevel, NULL);CHKERRQ(ierr); 1265 ierr = PetscOptionsEnd(); 1266 1267 /* add mesh loading options specific to the DM */ 1268 ierr = DMMoab_GetWriteOptions_Private(dmmoab->pcomm->get_id(), dmmoab->pcomm->size(), dmmoab->dim, MOAB_PARWOPTS_WRITE_PART, dbglevel, usrwriteopts, &writeopts);CHKERRQ(ierr); 1269 PetscInfo2(dm, "Writing file %s with options: %s\n",filename,writeopts); 1270 1271 /* output file, using parallel write */ 1272 merr = dmmoab->mbiface->write_file(filename, NULL, writeopts, &dmmoab->fileset, 1);MBERRVM(dmmoab->mbiface,"Writing output of DMMoab failed.",merr); 1273 PetscFunctionReturn(0); 1274 } 1275 1276