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