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