1 #include <petsc/private/dmmbimpl.h> /*I "petscdmmoab.h" I*/ 2 3 #include <petscdmmoab.h> 4 5 /*@C 6 DMMoabSetFieldVector - Sets the vector reference that represents the solution associated 7 with a particular field component. 8 9 Not Collective 10 11 Input Parameters: 12 + dm - the discretization manager object 13 . ifield - the index of the field as set before via DMMoabSetFieldName. 14 . fvec - the Vector solution corresponding to the field (component) 15 16 Level: intermediate 17 18 .keywords: discretization manager, set, component solution 19 20 .seealso: DMMoabGetFieldName(), DMMoabSetGlobalFieldVector() 21 @*/ 22 PetscErrorCode DMMoabSetFieldVector(DM dm, PetscInt ifield, Vec fvec) 23 { 24 DM_Moab *dmmoab; 25 moab::Tag vtag,ntag; 26 const PetscScalar *varray; 27 PetscScalar *farray; 28 moab::ErrorCode merr; 29 PetscErrorCode ierr; 30 std::string tag_name; 31 32 PetscFunctionBegin; 33 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 34 dmmoab = (DM_Moab*)(dm)->data; 35 36 if ((ifield < 0) || (ifield >= dmmoab->numFields)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "The field %d should be positive and less than %d.", ifield, dmmoab->numFields); 37 38 /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 39 merr = dmmoab->mbiface->tag_get_handle(dmmoab->fieldNames[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 40 moab::MB_TAG_DENSE|moab::MB_TAG_CREAT);MBERRNM(merr); 41 42 ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr); 43 44 merr = dmmoab->mbiface->tag_get_name(vtag, tag_name); 45 if (!tag_name.length() && merr !=moab::MB_SUCCESS) { 46 ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr); 47 /* use the entity handle and the Dof index to set the right value */ 48 merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)varray);MBERRNM(merr); 49 ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr); 50 } 51 else { 52 ierr = PetscMalloc1(dmmoab->nloc,&farray);CHKERRQ(ierr); 53 /* we are using a MOAB Vec - directly copy the tag data to new one */ 54 merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)farray);MBERRNM(merr); 55 merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr); 56 /* make sure the parallel exchange for ghosts are done appropriately */ 57 ierr = PetscFree(farray);CHKERRQ(ierr); 58 } 59 merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vowned);MBERRNM(merr); 60 PetscFunctionReturn(0); 61 } 62 63 64 /*@C 65 DMMoabSetGlobalFieldVector - Sets the vector reference that represents the global solution associated 66 with all fields (components) managed by DM. 67 A useful utility when updating the DM solution after a solve, to be serialized with the mesh for 68 checkpointing purposes. 69 70 Not Collective 71 72 Input Parameters: 73 + dm - the discretization manager object 74 . fvec - the global Vector solution corresponding to all the fields managed by DM 75 76 Level: intermediate 77 78 .keywords: discretization manager, set, component solution 79 80 .seealso: DMMoabGetFieldName(), DMMoabSetFieldVector() 81 @*/ 82 PetscErrorCode DMMoabSetGlobalFieldVector(DM dm, Vec fvec) 83 { 84 DM_Moab *dmmoab; 85 moab::Tag vtag,ntag; 86 const PetscScalar *rarray; 87 PetscScalar *varray,*farray; 88 moab::ErrorCode merr; 89 PetscErrorCode ierr; 90 PetscInt i,ifield; 91 std::string tag_name; 92 moab::Range::iterator iter; 93 94 PetscFunctionBegin; 95 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 96 dmmoab = (DM_Moab*)(dm)->data; 97 98 /* get the Tag corresponding to the global vector - possible that there is no tag associated.. */ 99 ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr); 100 merr = dmmoab->mbiface->tag_get_name(vtag, tag_name); 101 ierr = PetscMalloc1(dmmoab->nloc,&farray);CHKERRQ(ierr); 102 if (!tag_name.length() && merr !=moab::MB_SUCCESS) { 103 /* not a MOAB vector - use VecGetSubVector to get the parts as needed */ 104 105 ierr = VecGetArrayRead(fvec,&rarray);CHKERRQ(ierr); 106 for (ifield=0; ifield<dmmoab->numFields; ++ifield) { 107 108 /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 109 merr = dmmoab->mbiface->tag_get_handle(dmmoab->fieldNames[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 110 moab::MB_TAG_DENSE|moab::MB_TAG_CREAT);MBERRNM(merr); 111 112 for(i=0;i<dmmoab->nloc;i++) { 113 if (dmmoab->bs == 1) 114 farray[i]=rarray[ifield*dmmoab->nloc+i]; 115 else 116 farray[i]=rarray[i*dmmoab->numFields+ifield]; 117 } 118 119 /* use the entity handle and the Dof index to set the right value */ 120 merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr); 121 } 122 ierr = VecRestoreArrayRead(fvec,&rarray);CHKERRQ(ierr); 123 } 124 else { 125 ierr = PetscMalloc1(dmmoab->nloc*dmmoab->numFields,&varray);CHKERRQ(ierr); 126 127 /* we are using a MOAB Vec - directly copy the tag data to new one */ 128 merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)varray);MBERRNM(merr); 129 for (ifield=0; ifield<dmmoab->numFields; ++ifield) { 130 131 /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 132 merr = dmmoab->mbiface->tag_get_handle(dmmoab->fieldNames[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 133 moab::MB_TAG_DENSE|moab::MB_TAG_CREAT);MBERRNM(merr); 134 135 /* we are using a MOAB Vec - directly copy the tag data to new one */ 136 for(i=0; i < dmmoab->nloc; i++) { 137 if (dmmoab->bs == 1) 138 farray[i]=varray[ifield*dmmoab->nloc+i]; 139 else 140 farray[i]=varray[i*dmmoab->numFields+ifield]; 141 } 142 143 merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr); 144 /* make sure the parallel exchange for ghosts are done appropriately */ 145 merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr); 146 } 147 ierr = PetscFree(varray);CHKERRQ(ierr); 148 } 149 ierr = PetscFree(farray);CHKERRQ(ierr); 150 PetscFunctionReturn(0); 151 } 152 153 154 /*@C 155 DMMoabSetFieldNames - Sets the number of fields and their names to be managed by the DM 156 157 Not Collective 158 159 Input Parameters: 160 + dm - the discretization manager object 161 . numFields - the total number of fields 162 . fields - the array containing the names of each field (component); Can be NULL. 163 164 Level: intermediate 165 166 .keywords: discretization manager, set, component name 167 168 .seealso: DMMoabGetFieldName(), DMMoabSetFieldName() 169 @*/ 170 PetscErrorCode DMMoabSetFieldNames(DM dm,PetscInt numFields,const char* fields[]) 171 { 172 PetscErrorCode ierr; 173 PetscInt i; 174 DM_Moab *dmmoab; 175 176 PetscFunctionBegin; 177 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 178 dmmoab = (DM_Moab*)(dm)->data; 179 180 /* first deallocate the existing field structure */ 181 if (dmmoab->fieldNames) { 182 for(i=0; i<dmmoab->numFields; i++) { 183 ierr = PetscFree(dmmoab->fieldNames[i]);CHKERRQ(ierr); 184 } 185 ierr = PetscFree(dmmoab->fieldNames);CHKERRQ(ierr); 186 } 187 188 /* now re-allocate and assign field names */ 189 dmmoab->numFields = numFields; 190 ierr = PetscMalloc1(numFields,&dmmoab->fieldNames);CHKERRQ(ierr); 191 if (fields) { 192 for(i=0; i<dmmoab->numFields; i++) { 193 ierr = PetscStrallocpy(fields[i], (char**) &dmmoab->fieldNames[i]);CHKERRQ(ierr); 194 } 195 } 196 PetscFunctionReturn(0); 197 } 198 199 200 /*@C 201 DMMoabGetFieldName - Gets the names of individual field components in multicomponent 202 vectors associated with a DMDA. 203 204 Not Collective 205 206 Input Parameter: 207 + dm - the discretization manager object 208 . field - field number for the DMMoab (0, 1, ... dof-1), where dof indicates the 209 number of degrees of freedom per node within the DMMoab 210 211 Output Parameter: 212 . fieldName - the name of the field (component) 213 214 Level: intermediate 215 216 .keywords: discretization manager, get, component name 217 218 .seealso: DMMoabSetFieldName(), DMMoabSetFields() 219 @*/ 220 PetscErrorCode DMMoabGetFieldName(DM dm, PetscInt field, const char **fieldName) 221 { 222 DM_Moab *dmmoab; 223 224 PetscFunctionBegin; 225 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 226 dmmoab = (DM_Moab*)(dm)->data; 227 if ((field < 0) || (field >= dmmoab->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "DM field %d should be in [%d, %d)", field, 0, dmmoab->numFields); 228 229 *fieldName = dmmoab->fieldNames[field]; 230 PetscFunctionReturn(0); 231 } 232 233 234 /*@C 235 DMMoabSetFieldName - Sets the name of a field (component) managed by the DM 236 237 Not Collective 238 239 Input Parameters: 240 + dm - the discretization manager object 241 . field - the field number 242 . fieldName - the field (component) name 243 244 Level: intermediate 245 Notes: Can only be called after DMMoabSetFields supplied with correct numFields 246 247 .keywords: discretization manager, set, component name 248 249 .seealso: DMMoabGetFieldName(), DMMoabSetFields() 250 @*/ 251 PetscErrorCode DMMoabSetFieldName(DM dm, PetscInt field, const char *fieldName) 252 { 253 PetscErrorCode ierr; 254 DM_Moab *dmmoab; 255 256 PetscFunctionBegin; 257 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 258 PetscValidCharPointer(fieldName,3); 259 260 dmmoab = (DM_Moab*)(dm)->data; 261 if ((field < 0) || (field >= dmmoab->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "DM field %d should be in [%d, %d)", field, 0, dmmoab->numFields); 262 263 if (dmmoab->fieldNames[field]) { 264 ierr = PetscFree(dmmoab->fieldNames[field]);CHKERRQ(ierr); 265 } 266 ierr = PetscStrallocpy(fieldName, (char**) &dmmoab->fieldNames[field]);CHKERRQ(ierr); 267 PetscFunctionReturn(0); 268 } 269 270 271 /*@C 272 DMMoabGetFieldDof - Gets the global degree-of-freedom of a field (component) defined on a 273 particular MOAB EntityHandle. 274 275 Not Collective 276 277 Input Parameters: 278 + dm - the discretization manager object 279 . point - the MOAB EntityHandle container which holds the field degree-of-freedom values 280 . field - the field (component) index 281 282 Output Parameter: 283 + dof - the global degree-of-freedom index corresponding to the field in the discrete representation (Vec, Mat) 284 285 Level: beginner 286 287 .keywords: discretization manager, get, global degree of freedom 288 289 .seealso: DMMoabGetFieldDofs(), DMMoabGetFieldDofsLocal() 290 @*/ 291 PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point,PetscInt field,PetscInt* dof) 292 { 293 DM_Moab *dmmoab; 294 295 PetscFunctionBegin; 296 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 297 dmmoab = (DM_Moab*)(dm)->data; 298 299 *dof=dmmoab->gidmap[(PetscInt)point]; 300 PetscFunctionReturn(0); 301 } 302 303 304 /*@C 305 DMMoabGetFieldDofs - Gets the global degree-of-freedom of a field (component) defined on an 306 array of MOAB EntityHandles. 307 308 Not Collective 309 310 Input Parameters: 311 + dm - the discretization manager object 312 . npoints - the total number of Entities in the points array 313 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values 314 . field - the field (component) index 315 316 Output Parameter: 317 + dof - the global degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat) 318 319 Level: intermediate 320 321 .keywords: discretization manager, get, global degrees of freedom 322 323 .seealso: DMMoabGetFieldDof(), DMMoabGetFieldDofsLocal() 324 @*/ 325 PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof) 326 { 327 PetscInt i; 328 PetscErrorCode ierr; 329 DM_Moab *dmmoab; 330 331 PetscFunctionBegin; 332 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 333 PetscValidPointer(points,3); 334 dmmoab = (DM_Moab*)(dm)->data; 335 336 if (!dof) { 337 ierr = PetscMalloc1(npoints,&dof);CHKERRQ(ierr); 338 } 339 340 if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */ 341 for (i=0; i<npoints; ++i) 342 dof[i] = dmmoab->gidmap[(PetscInt)points[i]]*dmmoab->numFields+field; 343 } 344 else { 345 /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */ 346 /* TODO: eliminate the limitation using PetscSection to manage DOFs */ 347 for (i=0; i<npoints; ++i) 348 dof[i] = dmmoab->gidmap[(PetscInt)points[i]]+field*dmmoab->n; 349 } 350 PetscFunctionReturn(0); 351 } 352 353 354 /*@C 355 DMMoabGetFieldDofsLocal - Gets the local degrees-of-freedom of a field (component) defined on an 356 array of MOAB EntityHandles. 357 358 Not Collective 359 360 Input Parameters: 361 + dm - the discretization manager object 362 . npoints - the total number of Entities in the points array 363 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values 364 . field - the field (component) index 365 366 Output Parameter: 367 + dof - the local degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat) 368 369 Level: intermediate 370 371 .keywords: discretization manager, get, local degrees of freedom 372 373 .seealso: DMMoabGetFieldDof(), DMMoabGetFieldDofs() 374 @*/ 375 PetscErrorCode DMMoabGetFieldDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof) 376 { 377 PetscInt i,offset; 378 PetscErrorCode ierr; 379 DM_Moab *dmmoab; 380 381 PetscFunctionBegin; 382 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 383 PetscValidPointer(points,3); 384 dmmoab = (DM_Moab*)(dm)->data; 385 386 if (!dof) { 387 ierr = PetscMalloc1(npoints,&dof);CHKERRQ(ierr); 388 } 389 390 if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */ 391 for (i=0; i<npoints; ++i) 392 dof[i] = dmmoab->lidmap[(PetscInt)points[i]]*dmmoab->numFields+field; 393 } 394 else { 395 /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */ 396 /* TODO: eliminate the limitation using PetscSection to manage DOFs */ 397 offset = field*dmmoab->n; 398 for (i=0; i<npoints; ++i) 399 dof[i] = dmmoab->lidmap[(PetscInt)points[i]]+offset; 400 } 401 PetscFunctionReturn(0); 402 } 403 404 405 /*@C 406 DMMoabGetDofs - Gets the global degree-of-freedom for all fields (components) defined on an 407 array of MOAB EntityHandles. 408 409 Not Collective 410 411 Input Parameters: 412 + dm - the discretization manager object 413 . npoints - the total number of Entities in the points array 414 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values 415 416 Output Parameter: 417 + dof - the global degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat) 418 419 Level: intermediate 420 421 .keywords: discretization manager, get, global degrees of freedom 422 423 .seealso: DMMoabGetFieldDofs(), DMMoabGetDofsLocal(), DMMoabGetDofsBlocked() 424 @*/ 425 PetscErrorCode DMMoabGetDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 426 { 427 PetscInt i,field,offset; 428 PetscErrorCode ierr; 429 DM_Moab *dmmoab; 430 431 PetscFunctionBegin; 432 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 433 PetscValidPointer(points,3); 434 dmmoab = (DM_Moab*)(dm)->data; 435 436 if (!dof) { 437 ierr = PetscMalloc1(dmmoab->numFields*npoints,&dof);CHKERRQ(ierr); 438 } 439 440 if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */ 441 for (field=0; field<dmmoab->numFields; ++field) { 442 for (i=0; i<npoints; ++i) 443 dof[i*dmmoab->numFields+field] = dmmoab->gidmap[(PetscInt)points[i]]*dmmoab->numFields+field; 444 } 445 } 446 else { 447 /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */ 448 /* TODO: eliminate the limitation using PetscSection to manage DOFs */ 449 for (field=0; field<dmmoab->numFields; ++field) { 450 offset = field*dmmoab->n; 451 for (i=0; i<npoints; ++i) 452 dof[i*dmmoab->numFields+field] = dmmoab->gidmap[(PetscInt)points[i]]+offset; 453 } 454 } 455 PetscFunctionReturn(0); 456 } 457 458 459 /*@C 460 DMMoabGetDofsLocal - Gets the local degree-of-freedom for all fields (components) defined on an 461 array of MOAB EntityHandles. 462 463 Not Collective 464 465 Input Parameters: 466 + dm - the discretization manager object 467 . npoints - the total number of Entities in the points array 468 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values 469 470 Output Parameter: 471 + dof - the local degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat) 472 473 Level: intermediate 474 475 .keywords: discretization manager, get, global degrees of freedom 476 477 .seealso: DMMoabGetFieldDofs(), DMMoabGetDofs(), DMMoabGetDofsBlocked() 478 @*/ 479 PetscErrorCode DMMoabGetDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 480 { 481 PetscInt i,field,offset; 482 PetscErrorCode ierr; 483 DM_Moab *dmmoab; 484 485 PetscFunctionBegin; 486 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 487 PetscValidPointer(points,3); 488 dmmoab = (DM_Moab*)(dm)->data; 489 490 if (!dof) { 491 ierr = PetscMalloc1(dmmoab->numFields*npoints,&dof);CHKERRQ(ierr); 492 } 493 494 if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */ 495 for (field=0; field<dmmoab->numFields; ++field) { 496 for (i=0; i<npoints; ++i) 497 dof[i*dmmoab->numFields+field] = dmmoab->lidmap[(PetscInt)points[i]]*dmmoab->numFields+field; 498 } 499 } 500 else { 501 /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */ 502 /* TODO: eliminate the limitation using PetscSection to manage DOFs */ 503 for (field=0; field<dmmoab->numFields; ++field) { 504 offset = field*dmmoab->n; 505 for (i=0; i<npoints; ++i) 506 dof[i*dmmoab->numFields+field] = dmmoab->lidmap[(PetscInt)points[i]]+offset; 507 } 508 } 509 PetscFunctionReturn(0); 510 } 511 512 513 /*@C 514 DMMoabGetDofsBlocked - Gets the global degree-of-freedom for the first field (component) defined on an 515 array of MOAB EntityHandles. It is useful when performing Blocked(Get/Set) methods in computation 516 of element residuals and assembly of the discrete systems when all fields are co-located. 517 518 Not Collective 519 520 Input Parameters: 521 + dm - the discretization manager object 522 . npoints - the total number of Entities in the points array 523 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values 524 525 Output Parameter: 526 + dof - the blocked global degree-of-freedom index array in the discrete representation (Vec, Mat) 527 528 Level: intermediate 529 530 .keywords: discretization manager, get, global degrees of freedom 531 532 .seealso: DMMoabGetDofsLocal(), DMMoabGetDofs(), DMMoabGetDofsBlockedLocal() 533 @*/ 534 PetscErrorCode DMMoabGetDofsBlocked(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 535 { 536 PetscInt i; 537 DM_Moab *dmmoab; 538 PetscErrorCode ierr; 539 540 PetscFunctionBegin; 541 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 542 PetscValidPointer(points,2); 543 dmmoab = (DM_Moab*)(dm)->data; 544 545 if (!dof) { 546 ierr = PetscMalloc1(npoints,&dof);CHKERRQ(ierr); 547 } 548 549 for (i=0; i<npoints; ++i) { 550 dof[i]=dmmoab->gidmap[(PetscInt)points[i]]; 551 } 552 PetscFunctionReturn(0); 553 } 554 555 556 /*@C 557 DMMoabGetDofsBlockedLocal - Gets the local degree-of-freedom for the first field (component) defined on an 558 array of MOAB EntityHandles. It is useful when performing local Blocked(Get/Set) methods in computation 559 of element residuals and assembly of the discrete systems when all fields are co-located. 560 561 Not Collective 562 563 Input Parameters: 564 + dm - the discretization manager object 565 . npoints - the total number of Entities in the points array 566 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values 567 568 Output Parameter: 569 + dof - the blocked local degree-of-freedom index array in the discrete representation (Vec, Mat) 570 571 Level: intermediate 572 573 .keywords: discretization manager, get, global degrees of freedom 574 575 .seealso: DMMoabGetDofsLocal(), DMMoabGetDofs(), DMMoabGetDofsBlockedLocal() 576 @*/ 577 PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 578 { 579 PetscInt i; 580 DM_Moab *dmmoab; 581 PetscErrorCode ierr; 582 583 PetscFunctionBegin; 584 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 585 PetscValidPointer(points,2); 586 dmmoab = (DM_Moab*)(dm)->data; 587 588 if (!dof) { 589 ierr = PetscMalloc1(npoints,&dof);CHKERRQ(ierr); 590 } 591 592 for (i=0; i<npoints; ++i) 593 dof[i] = dmmoab->lidmap[(PetscInt)points[i]]; 594 PetscFunctionReturn(0); 595 } 596 597 598 /*@C 599 DMMoabGetVertexDofsBlocked - Gets the global degree-of-freedom for the first field (component) defined on an 600 array of locally owned MOAB mesh vertices. It's utility is when performing Finite-Difference type calculations 601 where vertex traversal is faster than element-wise assembly that is typically done in FEM calculations. 602 603 Not Collective 604 605 Input Parameters: 606 + dm - the discretization manager object 607 608 Output Parameter: 609 + dof - the blocked global degree-of-freedom index array in the discrete representation (Vec, Mat) that is vertex-based based on local numbering 610 611 Level: intermediate 612 613 .keywords: discretization manager, get, blocked degrees of freedom 614 615 .seealso: DMMoabGetVertexDofsBlockedLocal(), DMMoabGetDofsBlocked(), DMMoabGetDofsBlockedLocal() 616 @*/ 617 PetscErrorCode DMMoabGetVertexDofsBlocked(DM dm,PetscInt** dof) 618 { 619 DM_Moab *dmmoab; 620 621 PetscFunctionBegin; 622 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 623 dmmoab = (DM_Moab*)(dm)->data; 624 625 *dof = dmmoab->gidmap; 626 PetscFunctionReturn(0); 627 } 628 629 630 /*@C 631 DMMoabGetVertexDofsBlockedLocal - Gets the local degree-of-freedom for the first field (component) defined on an 632 array of locally owned MOAB mesh vertices. It's utility is when performing Finite-Difference type calculations 633 where vertex traversal is faster than element-wise assembly that is typically done in FEM calculations. 634 635 Not Collective 636 637 Input Parameters: 638 + dm - the discretization manager object 639 640 Output Parameter: 641 + dof - the blocked local degree-of-freedom index array in the discrete representation (Vec, Mat) that is vertex-based based on local numbering 642 643 Level: intermediate 644 645 .keywords: discretization manager, get, blocked degrees of freedom 646 647 .seealso: DMMoabGetVertexDofsBlocked(), DMMoabGetDofsBlocked(), DMMoabGetDofsBlockedLocal() 648 @*/ 649 PetscErrorCode DMMoabGetVertexDofsBlockedLocal(DM dm,PetscInt** dof) 650 { 651 DM_Moab *dmmoab; 652 653 PetscFunctionBegin; 654 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 655 PetscValidPointer(dof,2); 656 dmmoab = (DM_Moab*)(dm)->data; 657 658 *dof = dmmoab->lidmap; 659 PetscFunctionReturn(0); 660 } 661 662