Actual source code: dmmbvec.cxx

petsc-3.7.7 2017-09-25
Report Typos and Errors
  1: #include <petsc/private/dmmbimpl.h> /*I  "petscdmmoab.h"   I*/
  2: #include <petsc/private/vecimpl.h>

  4: #include <petscdmmoab.h>
  5: #include <MBTagConventions.hpp>

  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);

 15: /*@C
 16:   DMMoabCreateVector - Create a Vec from either an existing tag, or a specified tag size, and a range of entities

 18:   Collective on MPI_Comm

 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

 27:   Output Parameter:
 28: . vec             - The created vector

 30:   Level: beginner

 32: .keywords: Vec, create

 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;

 41:   if(!tag && (!range || range->empty())) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Both tag and range cannot be null.");

 43:   DMCreateVector_Moab_Private(dm,tag,range,is_global_vec,destroy_tag,vec);
 44:   return(0);
 45: }


 50: /*@C
 51:   DMMoabGetVecTag - Get the MOAB tag associated with this Vec

 53:   Input Parameter:
 54: . vec - Vec being queried

 56:   Output Parameter:
 57: . tag - Tag associated with this Vec. NULL if vec is a native PETSc Vec.

 59:   Level: beginner

 61: .keywords: DMMoab, MOAB tag

 63: .seealso: DMMoabCreateVector(), DMMoabGetVecRange()
 64: @*/
 65: PetscErrorCode DMMoabGetVecTag(Vec vec,moab::Tag *tag)
 66: {
 67:   PetscContainer  moabdata;
 68:   Vec_MOAB        *vmoab;
 69:   PetscErrorCode  ierr;


 74:   /* Get the MOAB private data */
 75:   PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);
 76:   PetscContainerGetPointer(moabdata, (void**) &vmoab);

 78:   *tag = vmoab->tag;
 79:   return(0);
 80: }


 85: /*@C
 86:   DMMoabGetVecRange - Get the MOAB entities associated with this Vec

 88:   Input Parameter:
 89: . vec   - Vec being queried

 91:   Output Parameter:
 92: . range - Entities associated with this Vec. NULL if vec is a native PETSc Vec.

 94:   Level: beginner

 96: .keywords: DMMoab, MOAB range

 98: .seealso: DMMoabCreateVector(), DMMoabGetVecTag()
 99: @*/
100: PetscErrorCode DMMoabGetVecRange(Vec vec,moab::Range *range)
101: {
102:   PetscContainer  moabdata;
103:   Vec_MOAB        *vmoab;
104:   PetscErrorCode  ierr;


109:   /* Get the MOAB private data handle */
110:   PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);
111:   PetscContainerGetPointer(moabdata, (void**) &vmoab);

113:   *range = *vmoab->tag_range;
114:   return(0);
115: }


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

123:   Collective on MPI_Comm

125:   Input Parameter:
126: . dm              - The DMMoab object being set
127: . vec             - The Vector whose underlying data is requested

129:   Output Parameter:
130: . array           - The local data array

132:   Level: intermediate

134: .keywords: MOAB, distributed array

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;

154:   dmmoab=(DM_Moab*)dm->data;

156:   /* Get the Vec_MOAB struct for the original vector */
157:   PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);
158:   PetscContainerGetPointer(moabdata, (void**)&vmoab);

160:   /* Get the real scalar array handle */
161:   varray = reinterpret_cast<PetscScalar**>(array);

163:   if (vmoab->is_native_vec) {

165:     /* get the local representation of the arrays from Vectors */
166:     DMGetLocalVector(dm,&vmoab->local);
167:     DMGlobalToLocalBegin(dm,vec,INSERT_VALUES,vmoab->local);
168:     DMGlobalToLocalEnd(dm,vec,INSERT_VALUES,vmoab->local);

170:     /* Get the Vec_MOAB struct for the original vector */
171:     PetscObjectQuery((PetscObject)vmoab->local,"MOABData", (PetscObject*) &moabdata);
172:     PetscContainerGetPointer(moabdata, (void**)&xmoab);

174:     /* get the local representation of the arrays from Vectors */
175:     VecGhostGetLocalForm(vmoab->local,&xmoab->local);
176:     VecGhostUpdateBegin(vmoab->local,INSERT_VALUES,SCATTER_FORWARD);
177:     VecGhostUpdateEnd(vmoab->local,INSERT_VALUES,SCATTER_FORWARD);

179:     VecGetArray(xmoab->local, varray);
180:   }
181:   else {

183:     /* Get the MOAB private data */
184:     DMMoabGetVecTag(vec,&vtag);

186:     /* exchange the data into ghost cells first */
187:     merr = dmmoab->pcomm->exchange_tags(vtag,*dmmoab->vlocal);MBERRNM(merr);

189:     PetscMalloc1((dmmoab->nloc+dmmoab->nghost)*dmmoab->numFields,varray);

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);

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:   return(0);
202: }


207: /*@C
208:   DMMoabVecRestoreArray - Restores the writable direct access array obtained via DMMoabVecGetArray

210:   Collective on MPI_Comm

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

217:   Level: intermediate

219: .keywords: MOAB, distributed array

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;

239:   dmmoab=(DM_Moab*)dm->data;

241:   /* Get the Vec_MOAB struct for the original vector */
242:   PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);
243:   PetscContainerGetPointer(moabdata, (void**)&vmoab);

245:   /* Get the real scalar array handle */
246:   varray = reinterpret_cast<PetscScalar**>(array);

248:   if (vmoab->is_native_vec) {

250:     /* Get the Vec_MOAB struct for the original vector */
251:     PetscObjectQuery((PetscObject)vmoab->local,"MOABData", (PetscObject*) &moabdata);
252:     PetscContainerGetPointer(moabdata, (void**)&xmoab);

254:     /* get the local representation of the arrays from Vectors */
255:     VecRestoreArray(xmoab->local, varray);
256:     VecGhostUpdateBegin(vmoab->local,ADD_VALUES,SCATTER_REVERSE);
257:     VecGhostUpdateEnd(vmoab->local,ADD_VALUES,SCATTER_REVERSE);
258:     VecGhostRestoreLocalForm(vmoab->local,&xmoab->local);

260:     /* restore local pieces */
261:     DMLocalToGlobalBegin(dm,vmoab->local,INSERT_VALUES,vec);
262:     DMLocalToGlobalEnd(dm,vmoab->local,INSERT_VALUES,vec);
263:     DMRestoreLocalVector(dm, &vmoab->local);
264:   }
265:   else {

267:     /* Get the MOAB private data */
268:     DMMoabGetVecTag(vec,&vtag);

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);

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:     }

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);

284:     PetscFree(*varray);
285:   }
286:   return(0);
287: }

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

294:   Collective on MPI_Comm

296:   Input Parameter:
297: + dm              - The DMMoab object being set
298: . vec             - The Vector whose underlying data is requested

300:   Output Parameter:
301: . array           - The local data array

303:   Level: intermediate

305: .keywords: MOAB, distributed array

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;

325:   dmmoab=(DM_Moab*)dm->data;

327:   /* Get the Vec_MOAB struct for the original vector */
328:   PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);
329:   PetscContainerGetPointer(moabdata, (void**)&vmoab);

331:   /* Get the real scalar array handle */
332:   varray = reinterpret_cast<PetscScalar**>(array);

334:   if (vmoab->is_native_vec) {
335:     /* get the local representation of the arrays from Vectors */
336:     DMGetLocalVector(dm,&vmoab->local);
337:     DMGlobalToLocalBegin(dm,vec,INSERT_VALUES,vmoab->local);
338:     DMGlobalToLocalEnd(dm,vec,INSERT_VALUES,vmoab->local);

340:     /* Get the Vec_MOAB struct for the original vector */
341:     PetscObjectQuery((PetscObject)vmoab->local,"MOABData", (PetscObject*) &moabdata);
342:     PetscContainerGetPointer(moabdata, (void**)&xmoab);

344:     /* get the local representation of the arrays from Vectors */
345:     VecGhostGetLocalForm(vmoab->local,&xmoab->local);
346:     VecGhostUpdateBegin(vmoab->local,INSERT_VALUES,SCATTER_FORWARD);
347:     VecGhostUpdateEnd(vmoab->local,INSERT_VALUES,SCATTER_FORWARD);
348:     VecGetArray(xmoab->local, varray);
349:   }
350:   else {
351:     /* Get the MOAB private data */
352:     DMMoabGetVecTag(vec,&vtag);

354:     /* exchange the data into ghost cells first */
355:     merr = dmmoab->pcomm->exchange_tags(vtag,*dmmoab->vlocal);MBERRNM(merr);

357:     PetscMalloc1((dmmoab->nloc+dmmoab->nghost)*dmmoab->numFields,varray);

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);

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:   return(0);
370: }


375: /*@C
376:   DMMoabVecRestoreArray - Restores the read-only direct access array obtained via DMMoabVecGetArray

378:   Collective on MPI_Comm

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

385:   Level: intermediate

387: .keywords: MOAB, distributed array

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;


403:   /* Get the Vec_MOAB struct for the original vector */
404:   PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);
405:   PetscContainerGetPointer(moabdata, (void**)&vmoab);

407:   /* Get the real scalar array handle */
408:   varray = reinterpret_cast<PetscScalar**>(array);

410:   if (vmoab->is_native_vec) {
411:     /* Get the Vec_MOAB struct for the original vector */
412:     PetscObjectQuery((PetscObject)vmoab->local,"MOABData", (PetscObject*) &moabdata);
413:     PetscContainerGetPointer(moabdata, (void**)&xmoab);

415:     /* restore the local representation of the arrays from Vectors */
416:     VecRestoreArray(xmoab->local, varray);
417:     VecGhostRestoreLocalForm(vmoab->local,&xmoab->local);

419:     /* restore local pieces */
420:     DMRestoreLocalVector(dm, &vmoab->local);
421:   }
422:   else {
423:     /* Nothing to do but just free the memory allocated before */
424:     PetscFree(*varray);

426:   }
427:   return(0);
428: }


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;

443:   Vec_MOAB *vmoab;
444:   DM_Moab *dmmoab = (DM_Moab*)dm->data;
445:   moab::ParallelComm *pcomm = dmmoab->pcomm;
446:   moab::Interface *mbiface = dmmoab->mbiface;

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.");

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);

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:   MPIU_Allreduce(&lnative_vec, &gnative_vec, 1, MPI_INT, MPI_MAX, pcomm->comm());

463:   /* Create the MOAB internal data object */
464:   PetscNew(&vmoab);
465:   vmoab->is_native_vec=(gnative_vec>0?PETSC_TRUE:PETSC_FALSE);

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:       DMVecCreateTagName_Moab_Private(pcomm,&tag_name);
473:       is_newtag = PETSC_TRUE;

475:       /* Create the default value for the tag (all zeros) */
476:       PetscCalloc1(dmmoab->numFields,&defaultvals);

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:       PetscFree(tag_name);
482:       PetscFree(defaultvals);
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:     }

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) {

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:       VecCreateGhostBlock(pcomm->comm(),dmmoab->bs,dmmoab->numFields*dmmoab->nloc,
506:                                  dmmoab->numFields*dmmoab->n,dmmoab->nghost,&dmmoab->gsindices[dmmoab->nloc],vec);
507:     }
508:     else {
509:       /* This is an MPI/SEQ Vector */
510:       VecCreate(pcomm->comm(),vec);
511:       VecSetSizes(*vec,dmmoab->numFields*dmmoab->nloc,PETSC_DECIDE);
512:       VecSetBlockSize(*vec,dmmoab->bs);
513:       VecSetType(*vec, VECMPI);
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);

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);

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:       VecCreateGhostBlockWithArray(pcomm->comm(),dmmoab->bs,dmmoab->numFields*dmmoab->nloc,
532:                                 dmmoab->numFields*dmmoab->n,dmmoab->nghost,&dmmoab->gsindices[dmmoab->nloc],data_ptr,vec);
533:     }
534:     else {
535:       /* This is an MPI Vector without ghosted padding */
536:       VecCreateMPIWithArray(pcomm->comm(),dmmoab->bs,dmmoab->numFields*range->size(),
537:                                 PETSC_DECIDE,data_ptr,vec);
538:     }
539:   }
540:   VecSetFromOptions(*vec);

542:   /* create a container and store the internal MOAB data for faster access based on Entities etc */
543:   PetscContainer moabdata;
544:   PetscContainerCreate(PETSC_COMM_WORLD,&moabdata);
545:   PetscContainerSetPointer(moabdata,vmoab);
546:   PetscContainerSetUserDestroy(moabdata,DMVecUserDestroy_Moab);
547:   PetscObjectCompose((PetscObject)*vec,"MOABData",(PetscObject)moabdata);
548:   (*vec)->ops->duplicate = DMVecDuplicate_Moab;
549:   PetscContainerDestroy(&moabdata);

551:   /* Vector created, manually set local to global mapping */
552:   if (dmmoab->ltog_map) {
553:     VecSetLocalToGlobalMapping(*vec,dmmoab->ltog_map);
554:   }

556:   /* set the DM reference to the vector */
557:   VecSetDM(*vec, dm);
558:   return(0);
559: }


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;

583:   const char*       PVEC_PREFIX      = "__PETSC_VEC_";
584:   PetscMalloc(PETSC_MAX_PATH_LEN, tag_name);

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);

597:   /* increment the new value of n */
598:   ++n;

600:   /* Make sure that n is consistent across all processes */
601:   MPIU_Allreduce(&n,&global_n,1,MPI_INT,MPI_MAX,pcomm->comm());

603:   /* Set the new name accordingly and return */
604:   PetscSNPrintf(*tag_name, PETSC_MAX_PATH_LEN-1, "%s_%D", PVEC_PREFIX, global_n);
605:   mberr = mbiface->tag_set_data(indexTag, &rootset, 1, (const void*)&global_n);MBERRNM(mberr);
606:   return(0);
607: }


612: PetscErrorCode DMCreateGlobalVector_Moab(DM dm,Vec *gvec)
613: {
614:   PetscErrorCode  ierr;
615:   DM_Moab         *dmmoab = (DM_Moab*)dm->data;

620:   DMCreateVector_Moab_Private(dm,PETSC_NULL,dmmoab->vowned,PETSC_TRUE,PETSC_TRUE,gvec);
621:   return(0);
622: }


627: PetscErrorCode DMCreateLocalVector_Moab(DM dm,Vec *lvec)
628: {
629:   PetscErrorCode  ierr;
630:   DM_Moab         *dmmoab = (DM_Moab*)dm->data;

635:   DMCreateVector_Moab_Private(dm,PETSC_NULL,dmmoab->vlocal,PETSC_FALSE,PETSC_TRUE,lvec);
636:   return(0);
637: }


642: PetscErrorCode DMVecDuplicate_Moab(Vec x,Vec *y)
643: {
645:   DM             dm;
646:   PetscContainer  moabdata;
647:   Vec_MOAB        *vmoab;


653:   /* Get the Vec_MOAB struct for the original vector */
654:   PetscObjectQuery((PetscObject)x,"MOABData", (PetscObject*) &moabdata);
655:   PetscContainerGetPointer(moabdata, (void**)&vmoab);

657:   VecGetDM(x, &dm);

660:   DMCreateVector_Moab_Private(dm,PETSC_NULL,vmoab->tag_range,vmoab->is_global_vec,PETSC_TRUE,y);
661:   VecSetDM(*y, dm);
662:   return(0);
663: }


668: PetscErrorCode DMVecUserDestroy_Moab(void *user)
669: {
670:   Vec_MOAB        *vmoab=(Vec_MOAB*)user;
671:   PetscErrorCode  ierr;
672:   moab::ErrorCode merr;

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:   PetscFree(vmoab);
684:   return(0);
685: }


690: PetscErrorCode  DMGlobalToLocalBegin_Moab(DM dm,Vec g,InsertMode mode,Vec l)
691: {
692:   PetscErrorCode    ierr;
693:   DM_Moab         *dmmoab = (DM_Moab*)dm->data;

696:   VecScatterBegin(dmmoab->ltog_sendrecv,g,l,mode,SCATTER_REVERSE);
697:   return(0);
698: }


703: PetscErrorCode  DMGlobalToLocalEnd_Moab(DM dm,Vec g,InsertMode mode,Vec l)
704: {
705:   PetscErrorCode    ierr;
706:   DM_Moab         *dmmoab = (DM_Moab*)dm->data;

709:   VecScatterEnd(dmmoab->ltog_sendrecv,g,l,mode,SCATTER_REVERSE);
710:   return(0);
711: }


716: PetscErrorCode  DMLocalToGlobalBegin_Moab(DM dm,Vec l,InsertMode mode,Vec g)
717: {
718:   PetscErrorCode    ierr;
719:   DM_Moab         *dmmoab = (DM_Moab*)dm->data;

722:   VecScatterBegin(dmmoab->ltog_sendrecv,l,g,mode,SCATTER_FORWARD);
723:   return(0);
724: }


729: PetscErrorCode  DMLocalToGlobalEnd_Moab(DM dm,Vec l,InsertMode mode,Vec g)
730: {
731:   PetscErrorCode    ierr;
732:   DM_Moab         *dmmoab = (DM_Moab*)dm->data;

735:   VecScatterEnd(dmmoab->ltog_sendrecv,l,g,mode,SCATTER_FORWARD);
736:   return(0);
737: }