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