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