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