1 #include <petsc/private/dmmbimpl.h> /*I "petscdmmoab.h" I*/ 2 #include <petsc/private/vecimpl.h> 3 4 #include <petscdmmoab.h> 5 #include <MBTagConventions.hpp> 6 7 /* declare some private DMMoab specific overrides */ 8 static PetscErrorCode DMCreateVector_Moab_Private(DM dm,moab::Tag tag,const moab::Range* userrange,PetscBool is_global_vec,PetscBool destroy_tag,Vec *vec); 9 static PetscErrorCode DMVecUserDestroy_Moab(void *user); 10 static PetscErrorCode DMVecDuplicate_Moab(Vec x,Vec *y); 11 #ifdef MOAB_HAVE_MPI 12 static PetscErrorCode DMVecCreateTagName_Moab_Private(moab::Interface *mbiface,moab::ParallelComm *pcomm,char** tag_name); 13 #else 14 static PetscErrorCode DMVecCreateTagName_Moab_Private(moab::Interface *mbiface,char** tag_name); 15 #endif 16 17 /*@C 18 DMMoabCreateVector - Create a Vec from either an existing tag, or a specified tag size, and a range of entities 19 20 Collective on MPI_Comm 21 22 Input Parameter: 23 + dm - The DMMoab object being set 24 . tag - If non-zero, block size will be taken from the tag size 25 . range - If non-empty, Vec corresponds to these entities, otherwise to the entities set on the DMMoab 26 . is_global_vec - If true, this is a local representation of the Vec (including ghosts in parallel), otherwise a truly parallel one 27 . destroy_tag - If true, MOAB tag is destroyed with Vec, otherwise it is left on MOAB 28 29 Output Parameter: 30 . vec - The created vector 31 32 Level: beginner 33 34 .keywords: Vec, create 35 36 .seealso: VecCreate() 37 @*/ 38 PetscErrorCode DMMoabCreateVector(DM dm, moab::Tag tag, const moab::Range* range, PetscBool is_global_vec, PetscBool destroy_tag, Vec *vec) 39 { 40 PetscErrorCode ierr; 41 42 PetscFunctionBegin; 43 if (!tag && (!range || range->empty())) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Both tag and range cannot be null."); 44 45 ierr = DMCreateVector_Moab_Private(dm, tag, range, is_global_vec, destroy_tag, vec); CHKERRQ(ierr); 46 PetscFunctionReturn(0); 47 } 48 49 50 /*@C 51 DMMoabGetVecTag - Get the MOAB tag associated with this Vec 52 53 Input Parameter: 54 . vec - Vec being queried 55 56 Output Parameter: 57 . tag - Tag associated with this Vec. NULL if vec is a native PETSc Vec. 58 59 Level: beginner 60 61 .keywords: DMMoab, MOAB tag 62 63 .seealso: DMMoabCreateVector(), DMMoabGetVecRange() 64 @*/ 65 PetscErrorCode DMMoabGetVecTag(Vec vec, moab::Tag *tag) 66 { 67 PetscContainer moabdata; 68 Vec_MOAB *vmoab; 69 PetscErrorCode ierr; 70 71 PetscFunctionBegin; 72 PetscValidPointer(tag, 2); 73 74 /* Get the MOAB private data */ 75 ierr = PetscObjectQuery((PetscObject)vec, "MOABData", (PetscObject*) &moabdata); CHKERRQ(ierr); 76 ierr = PetscContainerGetPointer(moabdata, (void**) &vmoab); CHKERRQ(ierr); 77 78 *tag = vmoab->tag; 79 PetscFunctionReturn(0); 80 } 81 82 83 /*@C 84 DMMoabGetVecRange - Get the MOAB entities associated with this Vec 85 86 Input Parameter: 87 . vec - Vec being queried 88 89 Output Parameter: 90 . range - Entities associated with this Vec. NULL if vec is a native PETSc Vec. 91 92 Level: beginner 93 94 .keywords: DMMoab, MOAB range 95 96 .seealso: DMMoabCreateVector(), DMMoabGetVecTag() 97 @*/ 98 PetscErrorCode DMMoabGetVecRange(Vec vec, moab::Range *range) 99 { 100 PetscContainer moabdata; 101 Vec_MOAB *vmoab; 102 PetscErrorCode ierr; 103 104 PetscFunctionBegin; 105 PetscValidPointer(range, 2); 106 107 /* Get the MOAB private data handle */ 108 ierr = PetscObjectQuery((PetscObject)vec, "MOABData", (PetscObject*) &moabdata); CHKERRQ(ierr); 109 ierr = PetscContainerGetPointer(moabdata, (void**) &vmoab); CHKERRQ(ierr); 110 111 *range = *vmoab->tag_range; 112 PetscFunctionReturn(0); 113 } 114 115 116 /*@C 117 DMMoabVecGetArray - Returns the writable direct access array to the local representation of MOAB tag data for the underlying vector using locally owned+ghosted range of entities 118 119 Collective on MPI_Comm 120 121 Input Parameter: 122 . dm - The DMMoab object being set 123 . vec - The Vector whose underlying data is requested 124 125 Output Parameter: 126 . array - The local data array 127 128 Level: intermediate 129 130 .keywords: MOAB, distributed array 131 132 .seealso: DMMoabVecRestoreArray(), DMMoabVecGetArrayRead(), DMMoabVecRestoreArrayRead() 133 @*/ 134 PetscErrorCode DMMoabVecGetArray(DM dm, Vec vec, void* array) 135 { 136 DM_Moab *dmmoab; 137 moab::ErrorCode merr; 138 PetscErrorCode ierr; 139 PetscInt count, i, f; 140 moab::Tag vtag; 141 PetscScalar **varray; 142 PetscScalar *marray; 143 PetscContainer moabdata; 144 Vec_MOAB *vmoab, *xmoab; 145 146 PetscFunctionBegin; 147 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 148 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 149 PetscValidPointer(array, 3); 150 dmmoab = (DM_Moab*)dm->data; 151 152 /* Get the Vec_MOAB struct for the original vector */ 153 ierr = PetscObjectQuery((PetscObject)vec, "MOABData", (PetscObject*) &moabdata); CHKERRQ(ierr); 154 ierr = PetscContainerGetPointer(moabdata, (void**)&vmoab); CHKERRQ(ierr); 155 156 /* Get the real scalar array handle */ 157 varray = reinterpret_cast<PetscScalar**>(array); 158 159 if (vmoab->is_native_vec) { 160 161 /* get the local representation of the arrays from Vectors */ 162 ierr = DMGetLocalVector(dm, &vmoab->local); CHKERRQ(ierr); 163 ierr = DMGlobalToLocalBegin(dm, vec, INSERT_VALUES, vmoab->local); CHKERRQ(ierr); 164 ierr = DMGlobalToLocalEnd(dm, vec, INSERT_VALUES, vmoab->local); CHKERRQ(ierr); 165 166 /* Get the Vec_MOAB struct for the original vector */ 167 ierr = PetscObjectQuery((PetscObject)vmoab->local, "MOABData", (PetscObject*) &moabdata); CHKERRQ(ierr); 168 ierr = PetscContainerGetPointer(moabdata, (void**)&xmoab); CHKERRQ(ierr); 169 170 /* get the local representation of the arrays from Vectors */ 171 ierr = VecGhostGetLocalForm(vmoab->local, &xmoab->local); CHKERRQ(ierr); 172 ierr = VecGhostUpdateBegin(vmoab->local, INSERT_VALUES, SCATTER_FORWARD); CHKERRQ(ierr); 173 ierr = VecGhostUpdateEnd(vmoab->local, INSERT_VALUES, SCATTER_FORWARD); CHKERRQ(ierr); 174 175 ierr = VecGetArray(xmoab->local, varray); CHKERRQ(ierr); 176 } 177 else { 178 179 /* Get the MOAB private data */ 180 ierr = DMMoabGetVecTag(vec, &vtag); CHKERRQ(ierr); 181 182 #ifdef MOAB_HAVE_MPI 183 /* exchange the data into ghost cells first */ 184 merr = dmmoab->pcomm->exchange_tags(vtag, *dmmoab->vlocal); MBERRNM(merr); 185 #endif 186 187 ierr = PetscMalloc1((dmmoab->nloc + dmmoab->nghost) * dmmoab->numFields, varray); CHKERRQ(ierr); 188 189 /* Get the array data for local entities */ 190 merr = dmmoab->mbiface->tag_iterate(vtag, dmmoab->vlocal->begin(), dmmoab->vlocal->end(), count, reinterpret_cast<void*&>(marray), false); MBERRNM(merr); 191 if (count != dmmoab->nloc + dmmoab->nghost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mismatch between local vertices and tag partition for Vec. %D != %D.", count, dmmoab->nloc + dmmoab->nghost); 192 193 i = 0; 194 for (moab::Range::iterator iter = dmmoab->vlocal->begin(); iter != dmmoab->vlocal->end(); iter++) { 195 for (f = 0; f < dmmoab->numFields; f++, i++) 196 (*varray)[dmmoab->lidmap[(PetscInt)*iter - dmmoab->seqstart]*dmmoab->numFields + f] = marray[i]; 197 //(*varray)[dmmoab->llmap[dmmoab->lidmap[((PetscInt)*iter-dmmoab->seqstart)]*dmmoab->numFields+f]]=marray[i]; 198 } 199 } 200 PetscFunctionReturn(0); 201 } 202 203 204 #undef __FUNCT__ 205 #define __FUNCT__ "DMMoabVecRestoreArray" 206 /*@C 207 DMMoabVecRestoreArray - Restores the writable direct access array obtained via DMMoabVecGetArray 208 209 Collective on MPI_Comm 210 211 Input Parameter: 212 + dm - The DMMoab object being set 213 . vec - The Vector whose underlying data is requested 214 . array - The local data array 215 216 Level: intermediate 217 218 .keywords: MOAB, distributed array 219 220 .seealso: DMMoabVecGetArray(), DMMoabVecGetArrayRead(), DMMoabVecRestoreArrayRead() 221 @*/ 222 PetscErrorCode DMMoabVecRestoreArray(DM dm, Vec vec, void* array) 223 { 224 DM_Moab *dmmoab; 225 moab::ErrorCode merr; 226 PetscErrorCode ierr; 227 moab::Tag vtag; 228 PetscInt count, i, f; 229 PetscScalar **varray; 230 PetscScalar *marray; 231 PetscContainer moabdata; 232 Vec_MOAB *vmoab, *xmoab; 233 234 PetscFunctionBegin; 235 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 236 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 237 PetscValidPointer(array, 3); 238 dmmoab = (DM_Moab*)dm->data; 239 240 /* Get the Vec_MOAB struct for the original vector */ 241 ierr = PetscObjectQuery((PetscObject)vec, "MOABData", (PetscObject*) &moabdata); CHKERRQ(ierr); 242 ierr = PetscContainerGetPointer(moabdata, (void**)&vmoab); CHKERRQ(ierr); 243 244 /* Get the real scalar array handle */ 245 varray = reinterpret_cast<PetscScalar**>(array); 246 247 if (vmoab->is_native_vec) { 248 249 /* Get the Vec_MOAB struct for the original vector */ 250 ierr = PetscObjectQuery((PetscObject)vmoab->local, "MOABData", (PetscObject*) &moabdata); CHKERRQ(ierr); 251 ierr = PetscContainerGetPointer(moabdata, (void**)&xmoab); CHKERRQ(ierr); 252 253 /* get the local representation of the arrays from Vectors */ 254 ierr = VecRestoreArray(xmoab->local, varray); CHKERRQ(ierr); 255 ierr = VecGhostUpdateBegin(vmoab->local, ADD_VALUES, SCATTER_REVERSE); CHKERRQ(ierr); 256 ierr = VecGhostUpdateEnd(vmoab->local, ADD_VALUES, SCATTER_REVERSE); CHKERRQ(ierr); 257 ierr = VecGhostRestoreLocalForm(vmoab->local, &xmoab->local); CHKERRQ(ierr); 258 259 /* restore local pieces */ 260 ierr = DMLocalToGlobalBegin(dm, vmoab->local, INSERT_VALUES, vec); CHKERRQ(ierr); 261 ierr = DMLocalToGlobalEnd(dm, vmoab->local, INSERT_VALUES, vec); CHKERRQ(ierr); 262 ierr = DMRestoreLocalVector(dm, &vmoab->local); CHKERRQ(ierr); 263 } 264 else { 265 266 /* Get the MOAB private data */ 267 ierr = DMMoabGetVecTag(vec, &vtag); CHKERRQ(ierr); 268 269 /* Get the array data for local entities */ 270 merr = dmmoab->mbiface->tag_iterate(vtag, dmmoab->vlocal->begin(), dmmoab->vlocal->end(), count, reinterpret_cast<void*&>(marray), false); MBERRNM(merr); 271 if (count != dmmoab->nloc + dmmoab->nghost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mismatch between local vertices and tag partition for Vec. %D != %D.", count, dmmoab->nloc + dmmoab->nghost); 272 273 i = 0; 274 for (moab::Range::iterator iter = dmmoab->vlocal->begin(); iter != dmmoab->vlocal->end(); iter++) { 275 for (f = 0; f < dmmoab->numFields; f++, i++) 276 marray[i] = (*varray)[dmmoab->lidmap[(PetscInt) * iter - dmmoab->seqstart] * dmmoab->numFields + f]; 277 //marray[i] = (*varray)[dmmoab->llmap[dmmoab->lidmap[((PetscInt)*iter-dmmoab->seqstart)]*dmmoab->numFields+f]]; 278 } 279 280 #ifdef MOAB_HAVE_MPI 281 /* reduce the tags correctly -> should probably let the user choose how to reduce in the future 282 For all FEM residual based assembly calculations, MPI_SUM should serve well */ 283 merr = dmmoab->pcomm->reduce_tags(vtag, MPI_SUM, *dmmoab->vlocal); MBERRV(dmmoab->mbiface, merr); 284 #endif 285 ierr = PetscFree(*varray); CHKERRQ(ierr); 286 } 287 PetscFunctionReturn(0); 288 } 289 290 /*@C 291 DMMoabVecGetArrayRead - Returns the read-only direct access array to the local representation of MOAB tag data for the underlying vector using locally owned+ghosted range of entities 292 293 Collective on MPI_Comm 294 295 Input Parameter: 296 + dm - The DMMoab object being set 297 . vec - The Vector whose underlying data is requested 298 299 Output Parameter: 300 . array - The local data array 301 302 Level: intermediate 303 304 .keywords: MOAB, distributed array 305 306 .seealso: DMMoabVecRestoreArrayRead(), DMMoabVecGetArray(), DMMoabVecRestoreArray() 307 @*/ 308 PetscErrorCode DMMoabVecGetArrayRead(DM dm, Vec vec, void* array) 309 { 310 DM_Moab *dmmoab; 311 moab::ErrorCode merr; 312 PetscErrorCode ierr; 313 PetscInt count, i, f; 314 moab::Tag vtag; 315 PetscScalar **varray; 316 PetscScalar *marray; 317 PetscContainer moabdata; 318 Vec_MOAB *vmoab, *xmoab; 319 320 PetscFunctionBegin; 321 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 322 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 323 PetscValidPointer(array, 3); 324 dmmoab = (DM_Moab*)dm->data; 325 326 /* Get the Vec_MOAB struct for the original vector */ 327 ierr = PetscObjectQuery((PetscObject)vec, "MOABData", (PetscObject*) &moabdata); CHKERRQ(ierr); 328 ierr = PetscContainerGetPointer(moabdata, (void**)&vmoab); CHKERRQ(ierr); 329 330 /* Get the real scalar array handle */ 331 varray = reinterpret_cast<PetscScalar**>(array); 332 333 if (vmoab->is_native_vec) { 334 /* get the local representation of the arrays from Vectors */ 335 ierr = DMGetLocalVector(dm, &vmoab->local); CHKERRQ(ierr); 336 ierr = DMGlobalToLocalBegin(dm, vec, INSERT_VALUES, vmoab->local); CHKERRQ(ierr); 337 ierr = DMGlobalToLocalEnd(dm, vec, INSERT_VALUES, vmoab->local); CHKERRQ(ierr); 338 339 /* Get the Vec_MOAB struct for the original vector */ 340 ierr = PetscObjectQuery((PetscObject)vmoab->local, "MOABData", (PetscObject*) &moabdata); CHKERRQ(ierr); 341 ierr = PetscContainerGetPointer(moabdata, (void**)&xmoab); CHKERRQ(ierr); 342 343 /* get the local representation of the arrays from Vectors */ 344 ierr = VecGhostGetLocalForm(vmoab->local, &xmoab->local); CHKERRQ(ierr); 345 ierr = VecGhostUpdateBegin(vmoab->local, INSERT_VALUES, SCATTER_FORWARD); CHKERRQ(ierr); 346 ierr = VecGhostUpdateEnd(vmoab->local, INSERT_VALUES, SCATTER_FORWARD); CHKERRQ(ierr); 347 ierr = VecGetArray(xmoab->local, varray); CHKERRQ(ierr); 348 } 349 else { 350 /* Get the MOAB private data */ 351 ierr = DMMoabGetVecTag(vec, &vtag); CHKERRQ(ierr); 352 353 #ifdef MOAB_HAVE_MPI 354 /* exchange the data into ghost cells first */ 355 merr = dmmoab->pcomm->exchange_tags(vtag, *dmmoab->vlocal); MBERRNM(merr); 356 #endif 357 ierr = PetscMalloc1((dmmoab->nloc + dmmoab->nghost) * dmmoab->numFields, varray); CHKERRQ(ierr); 358 359 /* Get the array data for local entities */ 360 merr = dmmoab->mbiface->tag_iterate(vtag, dmmoab->vlocal->begin(), dmmoab->vlocal->end(), count, reinterpret_cast<void*&>(marray), false); MBERRNM(merr); 361 if (count != dmmoab->nloc + dmmoab->nghost) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mismatch between local vertices and tag partition for Vec. %D != %D.", count, dmmoab->nloc + dmmoab->nghost); 362 363 i = 0; 364 for (moab::Range::iterator iter = dmmoab->vlocal->begin(); iter != dmmoab->vlocal->end(); iter++) { 365 for (f = 0; f < dmmoab->numFields; f++, i++) 366 (*varray)[dmmoab->lidmap[(PetscInt)*iter - dmmoab->seqstart]*dmmoab->numFields + f] = marray[i]; 367 //(*varray)[dmmoab->llmap[dmmoab->lidmap[((PetscInt)*iter-dmmoab->seqstart)]*dmmoab->numFields+f]]=marray[i]; 368 } 369 } 370 PetscFunctionReturn(0); 371 } 372 373 374 /*@C 375 DMMoabVecRestoreArray - Restores the read-only direct access array obtained via DMMoabVecGetArray 376 377 Collective on MPI_Comm 378 379 Input Parameter: 380 + dm - The DMMoab object being set 381 . vec - The Vector whose underlying data is requested 382 . array - The local data array 383 384 Level: intermediate 385 386 .keywords: MOAB, distributed array 387 388 .seealso: DMMoabVecGetArrayRead(), DMMoabVecGetArray(), DMMoabVecRestoreArray() 389 @*/ 390 PetscErrorCode DMMoabVecRestoreArrayRead(DM dm, Vec vec, void* array) 391 { 392 PetscErrorCode ierr; 393 PetscScalar **varray; 394 PetscContainer moabdata; 395 Vec_MOAB *vmoab, *xmoab; 396 397 PetscFunctionBegin; 398 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 399 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 400 PetscValidPointer(array, 3); 401 402 /* Get the Vec_MOAB struct for the original vector */ 403 ierr = PetscObjectQuery((PetscObject)vec, "MOABData", (PetscObject*) &moabdata); CHKERRQ(ierr); 404 ierr = PetscContainerGetPointer(moabdata, (void**)&vmoab); CHKERRQ(ierr); 405 406 /* Get the real scalar array handle */ 407 varray = reinterpret_cast<PetscScalar**>(array); 408 409 if (vmoab->is_native_vec) { 410 /* Get the Vec_MOAB struct for the original vector */ 411 ierr = PetscObjectQuery((PetscObject)vmoab->local, "MOABData", (PetscObject*) &moabdata); CHKERRQ(ierr); 412 ierr = PetscContainerGetPointer(moabdata, (void**)&xmoab); CHKERRQ(ierr); 413 414 /* restore the local representation of the arrays from Vectors */ 415 ierr = VecRestoreArray(xmoab->local, varray); CHKERRQ(ierr); 416 ierr = VecGhostRestoreLocalForm(vmoab->local, &xmoab->local); CHKERRQ(ierr); 417 418 /* restore local pieces */ 419 ierr = DMRestoreLocalVector(dm, &vmoab->local); CHKERRQ(ierr); 420 } 421 else { 422 /* Nothing to do but just free the memory allocated before */ 423 ierr = PetscFree(*varray); CHKERRQ(ierr); 424 425 } 426 PetscFunctionReturn(0); 427 } 428 429 430 #undef __FUNCT__ 431 #define __FUNCT__ "DMCreateVector_Moab_Private" 432 PetscErrorCode DMCreateVector_Moab_Private(DM dm, moab::Tag tag, const moab::Range* userrange, PetscBool is_global_vec, PetscBool destroy_tag, Vec *vec) 433 { 434 PetscErrorCode ierr; 435 moab::ErrorCode merr; 436 PetscBool is_newtag; 437 const moab::Range *range; 438 PetscInt count, lnative_vec, gnative_vec; 439 std::string ttname; 440 PetscScalar *data_ptr, *defaultvals; 441 442 Vec_MOAB *vmoab; 443 DM_Moab *dmmoab = (DM_Moab*)dm->data; 444 #ifdef MOAB_HAVE_MPI 445 moab::ParallelComm *pcomm = dmmoab->pcomm; 446 #endif 447 moab::Interface *mbiface = dmmoab->mbiface; 448 449 PetscFunctionBegin; 450 if (sizeof(PetscReal) != sizeof(PetscScalar)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "MOAB tags only support Real types (Complex-type unsupported)"); 451 if (!userrange) range = dmmoab->vowned; 452 else range = userrange; 453 if (!range) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Input range cannot be empty or call DMSetUp first."); 454 455 #ifndef USE_NATIVE_PETSCVEC 456 /* If the tag data is in a single sequence, use PETSc native vector since tag_iterate isn't useful anymore */ 457 lnative_vec = (range->psize() - 1); 458 #else 459 lnative_vec = 1; /* NOTE: Testing PETSc vector: will force to create native vector all the time */ 460 // lnative_vec=0; /* NOTE: Testing MOAB vector: will force to create MOAB tag_iterate based vector all the time */ 461 #endif 462 463 #ifdef MOAB_HAVE_MPI 464 ierr = MPIU_Allreduce(&lnative_vec, &gnative_vec, 1, MPI_INT, MPI_MAX, (((PetscObject)dm)->comm)); CHKERRQ(ierr); 465 #else 466 gnative_vec = lnative_vec; 467 #endif 468 469 /* Create the MOAB internal data object */ 470 ierr = PetscNew(&vmoab); CHKERRQ(ierr); 471 vmoab->is_native_vec = (gnative_vec > 0 ? PETSC_TRUE : PETSC_FALSE); 472 473 if (!vmoab->is_native_vec) { 474 merr = moab::MB_SUCCESS; 475 if (tag != 0) merr = mbiface->tag_get_name(tag, ttname); 476 if (!ttname.length() || merr != moab::MB_SUCCESS) { 477 /* get the new name for the anonymous MOABVec -> the tag_name will be destroyed along with Tag */ 478 char *tag_name = NULL; 479 #ifdef MOAB_HAVE_MPI 480 ierr = DMVecCreateTagName_Moab_Private(mbiface,pcomm,&tag_name); CHKERRQ(ierr); 481 #else 482 ierr = DMVecCreateTagName_Moab_Private(mbiface,&tag_name); CHKERRQ(ierr); 483 #endif 484 is_newtag = PETSC_TRUE; 485 486 /* Create the default value for the tag (all zeros) */ 487 ierr = PetscCalloc1(dmmoab->numFields, &defaultvals); CHKERRQ(ierr); 488 489 /* Create the tag */ 490 merr = mbiface->tag_get_handle(tag_name, dmmoab->numFields, moab::MB_TYPE_DOUBLE, tag, 491 moab::MB_TAG_DENSE | moab::MB_TAG_CREAT, defaultvals); MBERRNM(merr); 492 ierr = PetscFree(tag_name); CHKERRQ(ierr); 493 ierr = PetscFree(defaultvals); CHKERRQ(ierr); 494 } 495 else { 496 /* Make sure the tag data is of type "double" */ 497 moab::DataType tag_type; 498 merr = mbiface->tag_get_data_type(tag, tag_type); MBERRNM(merr); 499 if (tag_type != moab::MB_TYPE_DOUBLE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Tag data type must be MB_TYPE_DOUBLE"); 500 is_newtag = destroy_tag; 501 } 502 503 vmoab->tag = tag; 504 vmoab->new_tag = is_newtag; 505 } 506 vmoab->mbiface = mbiface; 507 #ifdef MOAB_HAVE_MPI 508 vmoab->pcomm = pcomm; 509 #endif 510 vmoab->is_global_vec = is_global_vec; 511 vmoab->tag_size = dmmoab->bs; 512 513 if (vmoab->is_native_vec) { 514 515 /* Create the PETSc Vector directly and attach our functions accordingly */ 516 if (!is_global_vec) { 517 /* This is an MPI Vector with ghosted padding */ 518 ierr = VecCreateGhostBlock((((PetscObject)dm)->comm), dmmoab->bs, dmmoab->numFields * dmmoab->nloc, 519 dmmoab->numFields * dmmoab->n, dmmoab->nghost, &dmmoab->gsindices[dmmoab->nloc], vec); CHKERRQ(ierr); 520 } 521 else { 522 /* This is an MPI/SEQ Vector */ 523 ierr = VecCreate((((PetscObject)dm)->comm), vec); CHKERRQ(ierr); 524 ierr = VecSetSizes(*vec, dmmoab->numFields * dmmoab->nloc, PETSC_DECIDE); CHKERRQ(ierr); 525 ierr = VecSetBlockSize(*vec, dmmoab->bs); CHKERRQ(ierr); 526 ierr = VecSetType(*vec, VECMPI); CHKERRQ(ierr); 527 } 528 } 529 else { 530 /* Call tag_iterate. This will cause MOAB to allocate memory for the 531 tag data if it hasn't already happened */ 532 merr = mbiface->tag_iterate(tag, range->begin(), range->end(), count, (void*&)data_ptr); MBERRNM(merr); 533 534 /* set the reference for vector range */ 535 vmoab->tag_range = new moab::Range(*range); 536 merr = mbiface->tag_get_length(tag, dmmoab->numFields); MBERRNM(merr); 537 538 /* Create the PETSc Vector 539 Query MOAB mesh to check if there are any ghosted entities 540 -> if we do, create a ghosted vector to map correctly to the same layout 541 -> else, create a non-ghosted parallel vector */ 542 if (!is_global_vec) { 543 /* This is an MPI Vector with ghosted padding */ 544 ierr = VecCreateGhostBlockWithArray((((PetscObject)dm)->comm), dmmoab->bs, dmmoab->numFields * dmmoab->nloc, 545 dmmoab->numFields * dmmoab->n, dmmoab->nghost, &dmmoab->gsindices[dmmoab->nloc], data_ptr, vec); CHKERRQ(ierr); 546 } 547 else { 548 /* This is an MPI Vector without ghosted padding */ 549 ierr = VecCreateMPIWithArray((((PetscObject)dm)->comm), dmmoab->bs, dmmoab->numFields * range->size(), 550 PETSC_DECIDE, data_ptr, vec); CHKERRQ(ierr); 551 } 552 } 553 ierr = VecSetFromOptions(*vec); CHKERRQ(ierr); 554 555 /* create a container and store the internal MOAB data for faster access based on Entities etc */ 556 PetscContainer moabdata; 557 ierr = PetscContainerCreate(PETSC_COMM_WORLD, &moabdata); CHKERRQ(ierr); 558 ierr = PetscContainerSetPointer(moabdata, vmoab); CHKERRQ(ierr); 559 ierr = PetscContainerSetUserDestroy(moabdata, DMVecUserDestroy_Moab); CHKERRQ(ierr); 560 ierr = PetscObjectCompose((PetscObject) * vec, "MOABData", (PetscObject)moabdata); CHKERRQ(ierr); 561 (*vec)->ops->duplicate = DMVecDuplicate_Moab; 562 ierr = PetscContainerDestroy(&moabdata); CHKERRQ(ierr); 563 564 /* Vector created, manually set local to global mapping */ 565 if (dmmoab->ltog_map) { 566 ierr = VecSetLocalToGlobalMapping(*vec, dmmoab->ltog_map); CHKERRQ(ierr); 567 } 568 569 /* set the DM reference to the vector */ 570 ierr = VecSetDM(*vec, dm); CHKERRQ(ierr); 571 PetscFunctionReturn(0); 572 } 573 574 575 /* DMVecCreateTagName_Moab_Private 576 * 577 * Creates a unique tag name that will be shared across processes. If 578 * pcomm is NULL, then this is a serial vector. A unique tag name 579 * will be returned in tag_name in either case. 580 * 581 * The tag names have the format _PETSC_VEC_N where N is some integer. 582 * 583 * NOTE: The tag_name is allocated in this routine; The user needs to free 584 * the character array. 585 */ 586 #ifdef MOAB_HAVE_MPI 587 PetscErrorCode DMVecCreateTagName_Moab_Private(moab::Interface *mbiface, moab::ParallelComm *pcomm, char** tag_name) 588 #else 589 PetscErrorCode DMVecCreateTagName_Moab_Private(moab::Interface *mbiface, char** tag_name) 590 #endif 591 { 592 moab::ErrorCode mberr; 593 PetscErrorCode ierr; 594 PetscInt n, global_n; 595 moab::Tag indexTag; 596 597 PetscFunctionBegin; 598 const char* PVEC_PREFIX = "__PETSC_VEC_"; 599 ierr = PetscMalloc1(PETSC_MAX_PATH_LEN, tag_name); CHKERRQ(ierr); 600 601 moab::EntityHandle rootset = mbiface->get_root_set(); 602 603 /* Check to see if there are any PETSc vectors defined */ 604 /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 605 mberr = mbiface->tag_get_handle("__PETSC_VECS__", 1, moab::MB_TYPE_INTEGER, indexTag, 606 moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT, 0); MBERRNM(mberr); 607 mberr = mbiface->tag_get_data(indexTag, &rootset, 1, &n); 608 if (mberr == moab::MB_TAG_NOT_FOUND) n = 0; /* this is the first temporary vector */ 609 else MBERRNM(mberr); 610 611 /* increment the new value of n */ 612 ++n; 613 614 #ifdef MOAB_HAVE_MPI 615 /* Make sure that n is consistent across all processes */ 616 ierr = MPIU_Allreduce(&n, &global_n, 1, MPI_INT, MPI_MAX, pcomm->comm()); CHKERRQ(ierr); 617 #else 618 global_n = n; 619 #endif 620 621 /* Set the new name accordingly and return */ 622 ierr = PetscSNPrintf(*tag_name, PETSC_MAX_PATH_LEN - 1, "%s_%D", PVEC_PREFIX, global_n); CHKERRQ(ierr); 623 mberr = mbiface->tag_set_data(indexTag, &rootset, 1, (const void*)&global_n); MBERRNM(mberr); 624 PetscFunctionReturn(0); 625 } 626 627 628 #undef __FUNCT__ 629 #define __FUNCT__ "DMCreateGlobalVector_Moab" 630 PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Moab(DM dm, Vec *gvec) 631 { 632 PetscErrorCode ierr; 633 DM_Moab *dmmoab = (DM_Moab*)dm->data; 634 635 PetscFunctionBegin; 636 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 637 PetscValidPointer(gvec,2); 638 ierr = DMCreateVector_Moab_Private(dm,NULL,dmmoab->vowned,PETSC_TRUE,PETSC_TRUE,gvec);CHKERRQ(ierr); 639 PetscFunctionReturn(0); 640 } 641 642 643 #undef __FUNCT__ 644 #define __FUNCT__ "DMCreateLocalVector_Moab" 645 PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Moab(DM dm, Vec *lvec) 646 { 647 PetscErrorCode ierr; 648 DM_Moab *dmmoab = (DM_Moab*)dm->data; 649 650 PetscFunctionBegin; 651 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 652 PetscValidPointer(lvec,2); 653 ierr = DMCreateVector_Moab_Private(dm,NULL,dmmoab->vlocal,PETSC_FALSE,PETSC_TRUE,lvec);CHKERRQ(ierr); 654 PetscFunctionReturn(0); 655 } 656 657 658 #undef __FUNCT__ 659 #define __FUNCT__ "DMVecDuplicate_Moab" 660 PetscErrorCode DMVecDuplicate_Moab(Vec x, Vec *y) 661 { 662 PetscErrorCode ierr; 663 DM dm; 664 PetscContainer moabdata; 665 Vec_MOAB *vmoab; 666 667 PetscFunctionBegin; 668 PetscValidHeaderSpecific(x, VEC_CLASSID, 1); 669 PetscValidPointer(y, 2); 670 671 /* Get the Vec_MOAB struct for the original vector */ 672 ierr = PetscObjectQuery((PetscObject)x, "MOABData", (PetscObject*) &moabdata); CHKERRQ(ierr); 673 ierr = PetscContainerGetPointer(moabdata, (void**)&vmoab); CHKERRQ(ierr); 674 675 ierr = VecGetDM(x, &dm); CHKERRQ(ierr); 676 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 677 678 ierr = DMCreateVector_Moab_Private(dm,NULL,vmoab->tag_range,vmoab->is_global_vec,PETSC_TRUE,y);CHKERRQ(ierr); 679 ierr = VecSetDM(*y, dm);CHKERRQ(ierr); 680 PetscFunctionReturn(0); 681 } 682 683 684 PetscErrorCode DMVecUserDestroy_Moab(void *user) 685 { 686 Vec_MOAB *vmoab = (Vec_MOAB*)user; 687 PetscErrorCode ierr; 688 moab::ErrorCode merr; 689 690 PetscFunctionBegin; 691 if (vmoab->new_tag && vmoab->tag) { 692 /* Tag was created via a call to VecDuplicate, delete the underlying tag in MOAB */ 693 merr = vmoab->mbiface->tag_delete(vmoab->tag); MBERRNM(merr); 694 } 695 delete vmoab->tag_range; 696 vmoab->tag = NULL; 697 vmoab->mbiface = NULL; 698 #ifdef MOAB_HAVE_MPI 699 vmoab->pcomm = NULL; 700 #endif 701 ierr = PetscFree(vmoab); CHKERRQ(ierr); 702 PetscFunctionReturn(0); 703 } 704 705 706 #undef __FUNCT__ 707 #define __FUNCT__ "DMGlobalToLocalBegin_Moab" 708 PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin_Moab(DM dm, Vec g, InsertMode mode, Vec l) 709 { 710 PetscErrorCode ierr; 711 DM_Moab *dmmoab = (DM_Moab*)dm->data; 712 713 PetscFunctionBegin; 714 ierr = VecScatterBegin(dmmoab->ltog_sendrecv, g, l, mode, SCATTER_REVERSE); CHKERRQ(ierr); 715 PetscFunctionReturn(0); 716 } 717 718 719 #undef __FUNCT__ 720 #define __FUNCT__ "DMGlobalToLocalEnd_Moab" 721 PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd_Moab(DM dm, Vec g, InsertMode mode, Vec l) 722 { 723 PetscErrorCode ierr; 724 DM_Moab *dmmoab = (DM_Moab*)dm->data; 725 726 PetscFunctionBegin; 727 ierr = VecScatterEnd(dmmoab->ltog_sendrecv, g, l, mode, SCATTER_REVERSE); CHKERRQ(ierr); 728 PetscFunctionReturn(0); 729 } 730 731 732 #undef __FUNCT__ 733 #define __FUNCT__ "DMLocalToGlobalBegin_Moab" 734 PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin_Moab(DM dm, Vec l, InsertMode mode, Vec g) 735 { 736 PetscErrorCode ierr; 737 DM_Moab *dmmoab = (DM_Moab*)dm->data; 738 739 PetscFunctionBegin; 740 ierr = VecScatterBegin(dmmoab->ltog_sendrecv, l, g, mode, SCATTER_FORWARD); CHKERRQ(ierr); 741 PetscFunctionReturn(0); 742 } 743 744 745 #undef __FUNCT__ 746 #define __FUNCT__ "DMLocalToGlobalEnd_Moab" 747 PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd_Moab(DM dm, Vec l, InsertMode mode, Vec g) 748 { 749 PetscErrorCode ierr; 750 DM_Moab *dmmoab = (DM_Moab*)dm->data; 751 752 PetscFunctionBegin; 753 ierr = VecScatterEnd(dmmoab->ltog_sendrecv, l, g, mode, SCATTER_FORWARD); CHKERRQ(ierr); 754 PetscFunctionReturn(0); 755 } 756 757