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