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