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