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,PetscInt tag_size,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 *gindices,*gsindices; 21 PetscInt i,count,icount,dof; 22 PetscInt size,rank; 23 PetscScalar *data_ptr; 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 moab::Tag ltog_tag = dmmoab->ltog_tag; 30 31 PetscFunctionBegin; 32 if(!userrange) range = dmmoab->vowned; 33 else range = userrange; 34 if(!range) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Input range cannot be empty or call DMSetUp first."); 35 36 if (!tag) { 37 /* get the new name for the anonymous MOABVec -> the tag_name will be destroyed along with Tag */ 38 char *tag_name = PETSC_NULL; 39 ierr = DMMoab_VecCreateTagName_Private(pcomm,&tag_name);CHKERRQ(ierr); 40 is_newtag = PETSC_TRUE; 41 42 /* Create the default value for the tag (all zeros) */ 43 std::vector<PetscScalar> default_value(tag_size, 0.0); 44 45 /* Create the tag */ 46 merr = mbiface->tag_get_handle(tag_name,tag_size,moab::MB_TYPE_DOUBLE,tag, 47 moab::MB_TAG_DENSE|moab::MB_TAG_CREAT,default_value.data());MBERRNM(merr); 48 ierr = PetscFree(tag_name);CHKERRQ(ierr); 49 } 50 else { 51 /* Make sure the tag data is of type "double" */ 52 moab::DataType tag_type; 53 merr = mbiface->tag_get_data_type(tag, tag_type);MBERRNM(merr); 54 if(tag_type != moab::MB_TYPE_DOUBLE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Tag data type must be MB_TYPE_DOUBLE"); 55 is_newtag = destroy_tag; 56 } 57 58 /* Create the MOAB internal data object */ 59 ierr = PetscNew(Vec_MOAB,&vmoab);CHKERRQ(ierr); 60 vmoab->tag = tag; 61 vmoab->mbiface = mbiface; 62 vmoab->pcomm = pcomm; 63 vmoab->new_tag = is_newtag; 64 vmoab->is_global_vec = is_global_vec; 65 merr = mbiface->tag_get_length(tag,vmoab->tag_size);MBERR("tag_get_size", merr); 66 67 /* set the reference for vector range */ 68 vmoab->tag_range = new moab::Range(*range); 69 70 /* Call tag_iterate. This will cause MOAB to allocate memory for the 71 tag data if it hasn't already happened */ 72 merr = mbiface->tag_iterate(tag,range->begin(),range->end(),count,(void*&)data_ptr);MBERRNM(merr); 73 74 /* Check to make sure the tag data is in a single sequence */ 75 if ((unsigned)count != range->size()) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only create MOAB Vector for single sequence"); 76 77 ierr = MPI_Comm_size(((PetscObject)dm)->comm, &size);CHKERRQ(ierr); 78 ierr = MPI_Comm_rank(((PetscObject)dm)->comm, &rank);CHKERRQ(ierr); 79 80 /* Create the PETSc Vector 81 Query MOAB mesh to check if there are any ghosted entities 82 -> if we do, create a ghosted vector to map correctly to the same layout 83 -> else, create a non-ghosted parallel vector */ 84 if (!is_global_vec && (size>1) && dmmoab->nghost) { 85 moab::Range::iterator iter; 86 ierr = PetscMalloc(dmmoab->nghost*sizeof(PetscInt), &gsindices);CHKERRQ(ierr); 87 88 for(iter = dmmoab->vghost->begin(),icount=0; iter != dmmoab->vghost->end(); iter++) { 89 merr = mbiface->tag_get_data(ltog_tag,&(*iter),1,&dof);MBERRNM(merr); 90 gsindices[icount++] = dof; 91 } 92 93 /* This is an MPI Vector with ghosted padding */ 94 ierr = VecCreateGhostBlockWithArray(vmoab->pcomm->comm(),vmoab->tag_size,vmoab->tag_size*dmmoab->nloc, 95 vmoab->tag_size*dmmoab->n,dmmoab->nghost,gsindices,data_ptr,vec);CHKERRQ(ierr); 96 97 ierr = PetscFree(gsindices);CHKERRQ(ierr); 98 } 99 else if (size>1) { 100 /* This is an MPI Vector without ghosted padding */ 101 ierr = VecCreateMPIWithArray(vmoab->pcomm->comm(),vmoab->tag_size,vmoab->tag_size*range->size(), 102 PETSC_DECIDE,data_ptr,vec);CHKERRQ(ierr); 103 } 104 else { 105 /* This is a sequential vector - valid only for the single processor case since MOAB tags are always partitioned 106 and we cannot define a Vec using the Tag array with size>1 will be of full length */ 107 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,vmoab->tag_size,vmoab->tag_size*dmmoab->n,data_ptr,vec);CHKERRQ(ierr); 108 } 109 110 PetscContainer moabdata; 111 ierr = PetscContainerCreate(PETSC_COMM_SELF,&moabdata);CHKERRQ(ierr); 112 ierr = PetscContainerSetPointer(moabdata,vmoab);CHKERRQ(ierr); 113 ierr = PetscContainerSetUserDestroy(moabdata,DMMoab_VecUserDestroy);CHKERRQ(ierr); 114 ierr = PetscObjectCompose((PetscObject)*vec,"MOABData",(PetscObject)moabdata);CHKERRQ(ierr); 115 (*vec)->ops->duplicate = DMMoab_VecDuplicate; 116 ierr = PetscContainerDestroy(&moabdata);CHKERRQ(ierr); 117 118 if (!dmmoab->ltog_map) { 119 /* Vector created, manually set local to global mapping */ 120 ierr = PetscMalloc(range->size()*sizeof(PetscInt)*vmoab->tag_size, &gindices);CHKERRQ(ierr); 121 moab::Range::iterator iter; 122 for(iter = range->begin(),count=0; iter != range->end(); iter++,count+=vmoab->tag_size) { 123 merr = mbiface->tag_get_data(ltog_tag,&(*iter),1,&dof);MBERRNM(merr); 124 for(i=0; i<vmoab->tag_size; ++i) 125 gindices[count+i] = (dof)*vmoab->tag_size+i; 126 } 127 128 ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,range->size(),gindices, 129 PETSC_COPY_VALUES,&dmmoab->ltog_map);CHKERRQ(ierr); 130 131 ierr = VecSetLocalToGlobalMappingBlock(*vec,dmmoab->ltog_map);CHKERRQ(ierr); 132 133 /* Clean up */ 134 ierr = PetscFree(gindices);CHKERRQ(ierr); 135 } 136 else { 137 ierr = VecSetLocalToGlobalMappingBlock(*vec,dmmoab->ltog_map);CHKERRQ(ierr); 138 } 139 140 /* set the DM reference to the vector */ 141 ierr = VecSetDM(*vec, dm);CHKERRQ(ierr); 142 PetscFunctionReturn(0); 143 } 144 145 146 #undef __FUNCT__ 147 #define __FUNCT__ "DMMoabCreateVector" 148 /*@ 149 DMMoabCreateVector - Create a Vec from either an existing tag, or a specified tag size, and a range of entities 150 151 Collective on MPI_Comm 152 153 Input Parameter: 154 . dm - The DMMoab object being set 155 . tag - If non-zero, block size will be taken from the tag size 156 . tag_size - If tag was zero, this parameter specifies the block size; unique tag name will be generated automatically 157 . range - If non-empty, Vec corresponds to these entities, otherwise to the entities set on the DMMoab 158 . is_global_vec - If true, this is a local representation of the Vec (including ghosts in parallel), otherwise a truly parallel one 159 . destroy_tag - If true, MOAB tag is destroyed with Vec, otherwise it is left on MOAB 160 161 Output Parameter: 162 . vec - The created vector 163 164 Level: beginner 165 166 .keywords: DMMoab, create 167 @*/ 168 PetscErrorCode DMMoabCreateVector(DM dm,moab::Tag tag,PetscInt tag_size,moab::Range* range,PetscBool is_global_vec,PetscBool destroy_tag,Vec *vec) 169 { 170 PetscErrorCode ierr; 171 172 PetscFunctionBegin; 173 if(!tag && !tag_size) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Both tag_size and tag cannot be null."); 174 175 ierr = DMMoab_CreateVector_Private(dm,tag,tag_size,range,is_global_vec,destroy_tag,vec);CHKERRQ(ierr); 176 PetscFunctionReturn(0); 177 } 178 179 180 #undef __FUNCT__ 181 #define __FUNCT__ "DMCreateGlobalVector_Moab" 182 PetscErrorCode DMCreateGlobalVector_Moab(DM dm,Vec *gvec) 183 { 184 PetscErrorCode ierr; 185 DM_Moab *dmmoab = (DM_Moab*)dm->data; 186 187 PetscFunctionBegin; 188 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 189 PetscValidPointer(gvec,2); 190 ierr = DMMoab_CreateVector_Private(dm,PETSC_NULL,dmmoab->bs,dmmoab->vowned,PETSC_TRUE,PETSC_TRUE,gvec);CHKERRQ(ierr); 191 PetscFunctionReturn(0); 192 } 193 194 195 #undef __FUNCT__ 196 #define __FUNCT__ "DMCreateLocalVector_Moab" 197 PetscErrorCode DMCreateLocalVector_Moab(DM dm,Vec *lvec) 198 { 199 PetscErrorCode ierr; 200 moab::Range vlocal; 201 DM_Moab *dmmoab = (DM_Moab*)dm->data; 202 203 PetscFunctionBegin; 204 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 205 PetscValidPointer(lvec,2); 206 vlocal = *dmmoab->vowned; 207 vlocal.merge(*dmmoab->vghost); 208 ierr = DMMoab_CreateVector_Private(dm,PETSC_NULL,dmmoab->bs,dmmoab->vowned,PETSC_FALSE,PETSC_TRUE,lvec);CHKERRQ(ierr); 209 PetscFunctionReturn(0); 210 } 211 212 213 #undef __FUNCT__ 214 #define __FUNCT__ "DMMoabGetVecTag" 215 /*@ 216 DMMoabGetVecTag - Get the MOAB tag associated with this Vec 217 218 Input Parameter: 219 . vec - Vec being queried 220 221 Output Parameter: 222 . tag - Tag associated with this Vec 223 224 Level: beginner 225 226 .keywords: DMMoab, create 227 @*/ 228 PetscErrorCode DMMoabGetVecTag(Vec vec,moab::Tag *tag) 229 { 230 PetscContainer moabdata; 231 Vec_MOAB *vmoab; 232 PetscErrorCode ierr; 233 234 PetscFunctionBegin; 235 PetscValidPointer(tag,2); 236 237 /* Get the MOAB private data */ 238 ierr = PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr); 239 ierr = PetscContainerGetPointer(moabdata, (void**) &vmoab);CHKERRQ(ierr); 240 241 *tag = vmoab->tag; 242 PetscFunctionReturn(0); 243 } 244 245 246 #undef __FUNCT__ 247 #define __FUNCT__ "DMMoabGetVecRange" 248 /*@ 249 DMMoabGetVecRange - Get the MOAB entities associated with this Vec 250 251 Input Parameter: 252 . vec - Vec being queried 253 254 Output Parameter: 255 . range - Entities associated with this Vec 256 257 Level: beginner 258 259 .keywords: DMMoab, create 260 @*/ 261 PetscErrorCode DMMoabGetVecRange(Vec vec,moab::Range *range) 262 { 263 PetscContainer moabdata; 264 Vec_MOAB *vmoab; 265 PetscErrorCode ierr; 266 267 PetscFunctionBegin; 268 PetscValidPointer(range,2); 269 270 /* Get the MOAB private data handle */ 271 ierr = PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr); 272 ierr = PetscContainerGetPointer(moabdata, (void**) &vmoab);CHKERRQ(ierr); 273 274 *range = *vmoab->tag_range; 275 PetscFunctionReturn(0); 276 } 277 278 279 #undef __FUNCT__ 280 #define __FUNCT__ "DMMoab_VecDuplicate" 281 PetscErrorCode DMMoab_VecDuplicate(Vec x,Vec *y) 282 { 283 PetscErrorCode ierr; 284 DM dm; 285 PetscContainer moabdata; 286 Vec_MOAB *vmoab; 287 288 PetscFunctionBegin; 289 PetscValidHeaderSpecific(x,VEC_CLASSID,1); 290 PetscValidPointer(y,2); 291 292 /* Get the Vec_MOAB struct for the original vector */ 293 ierr = PetscObjectQuery((PetscObject)x,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr); 294 ierr = PetscContainerGetPointer(moabdata, (void**)&vmoab);CHKERRQ(ierr); 295 296 ierr = VecGetDM(x, &dm);CHKERRQ(ierr); 297 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 298 299 ierr = DMMoab_CreateVector_Private(dm,PETSC_NULL,vmoab->tag_size,vmoab->tag_range,vmoab->is_global_vec,PETSC_TRUE,y);CHKERRQ(ierr); 300 ierr = VecSetDM(*y, dm);CHKERRQ(ierr); 301 PetscFunctionReturn(0); 302 } 303 304 305 #undef __FUNCT__ 306 #define __FUNCT__ "DMMoab_VecCreateTagName_Private" 307 /* DMMoab_VecCreateTagName_Private 308 * 309 * Creates a unique tag name that will be shared across processes. If 310 * pcomm is NULL, then this is a serial vector. A unique tag name 311 * will be returned in tag_name in either case. 312 * 313 * The tag names have the format _PETSC_VEC_N where N is some integer. 314 * 315 * NOTE: The tag_name is allocated in this routine; The user needs to free 316 * the character array. 317 */ 318 PetscErrorCode DMMoab_VecCreateTagName_Private(moab::ParallelComm *pcomm,char** tag_name) 319 { 320 moab::ErrorCode mberr; 321 PetscErrorCode ierr; 322 PetscInt n,global_n; 323 moab::Tag indexTag; 324 325 PetscFunctionBegin; 326 const char* PVEC_PREFIX = "__PETSC_VEC_"; 327 ierr = PetscMalloc(PETSC_MAX_PATH_LEN, tag_name);CHKERRQ(ierr); 328 329 /* Check to see if there are any PETSc vectors defined */ 330 moab::Interface *mbiface = pcomm->get_moab(); 331 moab::EntityHandle rootset = mbiface->get_root_set(); 332 333 /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 334 mberr = mbiface->tag_get_handle("__PETSC_VECS__",1,moab::MB_TYPE_INTEGER,indexTag, 335 moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT,0);MBERRNM(mberr); 336 mberr = mbiface->tag_get_data(indexTag, &rootset, 1, &n); 337 if (mberr == moab::MB_TAG_NOT_FOUND) n=0; /* this is the first temporary vector */ 338 else MBERRNM(mberr); 339 340 /* increment the new value of n */ 341 ++n; 342 343 /* Make sure that n is consistent across all processes */ 344 ierr = MPI_Allreduce(&n,&global_n,1,MPI_INT,MPI_MAX,pcomm->comm());CHKERRQ(ierr); 345 346 /* Set the new name accordingly and return */ 347 ierr = PetscSNPrintf(*tag_name, PETSC_MAX_PATH_LEN-1, "%s_%D", PVEC_PREFIX, global_n);CHKERRQ(ierr); 348 mberr = mbiface->tag_set_data(indexTag, &rootset, 1, (const void*)&global_n);MBERRNM(mberr); 349 PetscFunctionReturn(0); 350 } 351 352 353 #undef __FUNCT__ 354 #define __FUNCT__ "DMMoab_VecUserDestroy" 355 PetscErrorCode DMMoab_VecUserDestroy(void *user) 356 { 357 Vec_MOAB *vmoab=(Vec_MOAB*)user; 358 PetscErrorCode ierr; 359 moab::ErrorCode merr; 360 361 PetscFunctionBegin; 362 if(vmoab->new_tag && vmoab->tag) { 363 /* Tag was created via a call to VecDuplicate, delete the underlying tag in MOAB */ 364 merr = vmoab->mbiface->tag_delete(vmoab->tag);MBERRNM(merr); 365 } 366 delete vmoab->tag_range; 367 vmoab->tag = PETSC_NULL; 368 vmoab->mbiface = PETSC_NULL; 369 vmoab->pcomm = PETSC_NULL; 370 ierr = PetscFree(vmoab);CHKERRQ(ierr); 371 PetscFunctionReturn(0); 372 } 373 374 #undef __FUNCT__ 375 #define __FUNCT__ "DMMoabVecGetArray" 376 /*@ 377 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 378 379 Collective on MPI_Comm 380 381 Input Parameter: 382 . dm - The DMMoab object being set 383 . vec - The Vector whose underlying data is requested 384 385 Output Parameter: 386 . array - The local data array 387 388 Level: intermediate 389 390 .keywords: MOAB, distributed array 391 392 .seealso: DMMoabVecRestoreArray(), DMMoabVecGetArrayRead(), DMMoabVecRestoreArrayRead() 393 @*/ 394 PetscErrorCode DMMoabVecGetArray(DM dm,Vec vec,void* array) 395 { 396 DM_Moab *moab; 397 moab::ErrorCode merr; 398 PetscErrorCode ierr; 399 PetscInt count; 400 moab::Tag vtag; 401 PetscScalar **varray; 402 403 PetscFunctionBegin; 404 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 405 PetscValidHeaderSpecific(vec,VEC_CLASSID,2); 406 moab=(DM_Moab*)dm->data; 407 408 /* Get the MOAB private data */ 409 ierr = DMMoabGetVecTag(vec,&vtag);CHKERRQ(ierr); 410 411 /* Get the real scalar array handle */ 412 varray = reinterpret_cast<PetscScalar**>(array); 413 414 /* Get the array data for local entities */ 415 merr = moab->mbiface->tag_iterate(vtag,moab->vlocal->begin(),moab->vlocal->end(),count,reinterpret_cast<void*&>(*varray),true);MBERRNM(merr); 416 if (count!=(PetscInt)moab->vlocal->size()) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mismatch between local vertices and tag partition for Vec. %D != %D.",count,moab->vlocal->size()); 417 418 merr = moab->pcomm->exchange_tags(vtag,*moab->vlocal);MBERRNM(merr); 419 PetscFunctionReturn(0); 420 } 421 422 423 #undef __FUNCT__ 424 #define __FUNCT__ "DMMoabVecRestoreArray" 425 /*@ 426 DMMoabVecRestoreArray - Restores the writable direct access array obtained via DMMoabVecGetArray 427 428 Collective on MPI_Comm 429 430 Input Parameter: 431 + dm - The DMMoab object being set 432 . vec - The Vector whose underlying data is requested 433 - array - The local data array 434 435 Level: intermediate 436 437 .keywords: MOAB, distributed array 438 439 .seealso: DMMoabVecGetArray(), DMMoabVecGetArrayRead(), DMMoabVecRestoreArrayRead() 440 @*/ 441 PetscErrorCode DMMoabVecRestoreArray(DM dm,Vec v,void* array) 442 { 443 DM_Moab *moab; 444 moab::ErrorCode merr; 445 PetscErrorCode ierr; 446 moab::Tag vtag; 447 448 PetscFunctionBegin; 449 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 450 PetscValidHeaderSpecific(v,VEC_CLASSID,2); 451 moab=(DM_Moab*)dm->data; 452 453 /* Get the MOAB private data */ 454 ierr = DMMoabGetVecTag(v,&vtag);CHKERRQ(ierr); 455 456 /* reduce the tags correctly -> should probably let the user choose how to reduce in the future 457 For all FEM residual based assembly calculations, MPI_SUM should serve well 458 */ 459 merr = moab->pcomm->reduce_tags(vtag,MPI_SUM,*moab->vghost);MBERRNM(merr); 460 PetscFunctionReturn(0); 461 } 462 463 #undef __FUNCT__ 464 #define __FUNCT__ "DMMoabVecGetArrayRead" 465 /*@ 466 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 467 468 Collective on MPI_Comm 469 470 Input Parameter: 471 + dm - The DMMoab object being set 472 . vec - The Vector whose underlying data is requested 473 474 Output Parameter: 475 . array - The local data array 476 477 Level: intermediate 478 479 .keywords: MOAB, distributed array 480 481 .seealso: DMMoabVecRestoreArrayRead(), DMMoabVecGetArray(), DMMoabVecRestoreArray() 482 @*/ 483 PetscErrorCode DMMoabVecGetArrayRead(DM dm,Vec vec,void* array) 484 { 485 DM_Moab *moab; 486 moab::ErrorCode merr; 487 PetscErrorCode ierr; 488 PetscInt count; 489 moab::Tag vtag; 490 PetscScalar **varray; 491 492 PetscFunctionBegin; 493 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 494 PetscValidHeaderSpecific(vec,VEC_CLASSID,2); 495 moab=(DM_Moab*)dm->data; 496 497 /* Get the MOAB private data */ 498 ierr = DMMoabGetVecTag(vec,&vtag);CHKERRQ(ierr); 499 500 /* Get the real scalar array handle */ 501 varray = reinterpret_cast<PetscScalar**>(array); 502 503 /* Get the array data for local entities */ 504 merr = moab->mbiface->tag_iterate(vtag,moab->vlocal->begin(),moab->vlocal->end(),count,reinterpret_cast<void*&>(*varray),true);MBERRNM(merr); 505 if (count!=(PetscInt)moab->vlocal->size()) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mismatch between local vertices and tag partition for Vec. %D != %D.",count,moab->vlocal->size()); 506 507 merr = moab->pcomm->exchange_tags(vtag,*moab->vlocal);MBERRNM(merr); 508 PetscFunctionReturn(0); 509 } 510 511 512 #undef __FUNCT__ 513 #define __FUNCT__ "DMMoabVecRestoreArrayRead" 514 /*@ 515 DMMoabVecRestoreArray - Restores the read-only direct access array obtained via DMMoabVecGetArray 516 517 Collective on MPI_Comm 518 519 Input Parameter: 520 + dm - The DMMoab object being set 521 . vec - The Vector whose underlying data is requested 522 - array - The local data array 523 524 Level: intermediate 525 526 .keywords: MOAB, distributed array 527 528 .seealso: DMMoabVecGetArrayRead(), DMMoabVecGetArray(), DMMoabVecRestoreArray() 529 @*/ 530 PetscErrorCode DMMoabVecRestoreArrayRead(DM dm,Vec v,void* array) 531 { 532 PetscFunctionBegin; 533 /* Nothing to do -> do not free the array memory obtained from tag_iterate */ 534 PetscFunctionReturn(0); 535 } 536 537 538 #undef __FUNCT__ 539 #define __FUNCT__ "DMGlobalToLocalBegin_Moab" 540 PetscErrorCode DMGlobalToLocalBegin_Moab(DM dm,Vec g,InsertMode mode,Vec l) 541 { 542 PetscErrorCode ierr; 543 DM_Moab *dmmoab = (DM_Moab*)dm->data; 544 545 PetscFunctionBegin; 546 ierr = VecScatterBegin(dmmoab->ltog_sendrecv,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 547 PetscFunctionReturn(0); 548 } 549 550 551 #undef __FUNCT__ 552 #define __FUNCT__ "DMGlobalToLocalEnd_Moab" 553 PetscErrorCode DMGlobalToLocalEnd_Moab(DM dm,Vec g,InsertMode mode,Vec l) 554 { 555 PetscErrorCode ierr; 556 DM_Moab *dmmoab = (DM_Moab*)dm->data; 557 558 PetscFunctionBegin; 559 ierr = VecScatterEnd(dmmoab->ltog_sendrecv,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 560 PetscFunctionReturn(0); 561 } 562 563 564 #undef __FUNCT__ 565 #define __FUNCT__ "DMLocalToGlobalBegin_Moab" 566 PetscErrorCode DMLocalToGlobalBegin_Moab(DM dm,Vec l,InsertMode mode,Vec g) 567 { 568 PetscErrorCode ierr; 569 DM_Moab *dmmoab = (DM_Moab*)dm->data; 570 571 PetscFunctionBegin; 572 ierr = VecScatterBegin(dmmoab->ltog_sendrecv,l,g,mode,SCATTER_REVERSE);CHKERRQ(ierr); 573 PetscFunctionReturn(0); 574 } 575 576 577 #undef __FUNCT__ 578 #define __FUNCT__ "DMLocalToGlobalEnd_Moab" 579 PetscErrorCode DMLocalToGlobalEnd_Moab(DM dm,Vec l,InsertMode mode,Vec g) 580 { 581 PetscErrorCode ierr; 582 DM_Moab *dmmoab = (DM_Moab*)dm->data; 583 584 PetscFunctionBegin; 585 ierr = VecScatterEnd(dmmoab->ltog_sendrecv,l,g,mode,SCATTER_REVERSE);CHKERRQ(ierr); 586 PetscFunctionReturn(0); 587 } 588 589 590