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