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