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