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