xref: /petsc/src/dm/impls/moab/dmmbvec.cxx (revision f6829af0cc05a47053f5ab6972ac976581ca15ed)
1032b8ab6SVijay Mahadevan #include <petsc-private/dmmbimpl.h> /*I  "petscdm.h"   I*/
2032b8ab6SVijay Mahadevan #include <petsc-private/vecimpl.h> /*I  "petscdm.h"   I*/
3032b8ab6SVijay Mahadevan 
4032b8ab6SVijay Mahadevan #include <petscdmmoab.h>
5032b8ab6SVijay Mahadevan #include <MBTagConventions.hpp>
6032b8ab6SVijay Mahadevan 
7*f6829af0SVijay Mahadevan /* declare some private DMMoab specific overrides */
8*f6829af0SVijay Mahadevan static PetscErrorCode DMCreateVector_Moab_Private(DM dm,moab::Tag tag,const moab::Range* userrange,PetscBool is_global_vec,PetscBool destroy_tag,Vec *vec);
9*f6829af0SVijay Mahadevan static PetscErrorCode DMVecUserDestroy_Moab(void *user);
10*f6829af0SVijay Mahadevan static PetscErrorCode DMVecDuplicate_Moab(Vec x,Vec *y);
11*f6829af0SVijay Mahadevan static PetscErrorCode DMVecCreateTagName_Moab_Private(moab::ParallelComm *pcomm,char** tag_name);
12032b8ab6SVijay Mahadevan 
13032b8ab6SVijay Mahadevan #undef __FUNCT__
14*f6829af0SVijay Mahadevan #define __FUNCT__ "DMMoabCreateVector"
15*f6829af0SVijay Mahadevan /*@
16*f6829af0SVijay Mahadevan   DMMoabCreateVector - Create a Vec from either an existing tag, or a specified tag size, and a range of entities
17*f6829af0SVijay Mahadevan 
18*f6829af0SVijay Mahadevan   Collective on MPI_Comm
19*f6829af0SVijay Mahadevan 
20*f6829af0SVijay Mahadevan   Input Parameter:
21*f6829af0SVijay Mahadevan . dm              - The DMMoab object being set
22*f6829af0SVijay Mahadevan . tag             - If non-zero, block size will be taken from the tag size
23*f6829af0SVijay Mahadevan . range           - If non-empty, Vec corresponds to these entities, otherwise to the entities set on the DMMoab
24*f6829af0SVijay Mahadevan . is_global_vec   - If true, this is a local representation of the Vec (including ghosts in parallel), otherwise a truly parallel one
25*f6829af0SVijay Mahadevan . destroy_tag     - If true, MOAB tag is destroyed with Vec, otherwise it is left on MOAB
26*f6829af0SVijay Mahadevan 
27*f6829af0SVijay Mahadevan   Output Parameter:
28*f6829af0SVijay Mahadevan . vec             - The created vector
29*f6829af0SVijay Mahadevan 
30*f6829af0SVijay Mahadevan   Level: beginner
31*f6829af0SVijay Mahadevan 
32*f6829af0SVijay Mahadevan .keywords: DMMoab, create
33*f6829af0SVijay Mahadevan @*/
34*f6829af0SVijay Mahadevan PetscErrorCode DMMoabCreateVector(DM dm,moab::Tag tag,const moab::Range* range,PetscBool is_global_vec,PetscBool destroy_tag,Vec *vec)
35*f6829af0SVijay Mahadevan {
36*f6829af0SVijay Mahadevan   PetscErrorCode     ierr;
37*f6829af0SVijay Mahadevan 
38*f6829af0SVijay Mahadevan   PetscFunctionBegin;
39*f6829af0SVijay Mahadevan   if(!tag && (!range || range->empty())) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Both tag and range cannot be null.");
40*f6829af0SVijay Mahadevan 
41*f6829af0SVijay Mahadevan   ierr = DMCreateVector_Moab_Private(dm,tag,range,is_global_vec,destroy_tag,vec);CHKERRQ(ierr);
42*f6829af0SVijay Mahadevan   PetscFunctionReturn(0);
43*f6829af0SVijay Mahadevan }
44*f6829af0SVijay Mahadevan 
45*f6829af0SVijay Mahadevan 
46*f6829af0SVijay Mahadevan #undef __FUNCT__
47*f6829af0SVijay Mahadevan #define __FUNCT__ "DMMoabGetVecTag"
48*f6829af0SVijay Mahadevan /*@
49*f6829af0SVijay Mahadevan   DMMoabGetVecTag - Get the MOAB tag associated with this Vec
50*f6829af0SVijay Mahadevan 
51*f6829af0SVijay Mahadevan   Input Parameter:
52*f6829af0SVijay Mahadevan . vec - Vec being queried
53*f6829af0SVijay Mahadevan 
54*f6829af0SVijay Mahadevan   Output Parameter:
55*f6829af0SVijay Mahadevan . tag - Tag associated with this Vec. NULL if vec is a native PETSc Vec.
56*f6829af0SVijay Mahadevan 
57*f6829af0SVijay Mahadevan   Level: beginner
58*f6829af0SVijay Mahadevan 
59*f6829af0SVijay Mahadevan .keywords: DMMoab, create
60*f6829af0SVijay Mahadevan @*/
61*f6829af0SVijay Mahadevan PetscErrorCode DMMoabGetVecTag(Vec vec,moab::Tag *tag)
62*f6829af0SVijay Mahadevan {
63*f6829af0SVijay Mahadevan   PetscContainer  moabdata;
64*f6829af0SVijay Mahadevan   Vec_MOAB        *vmoab;
65*f6829af0SVijay Mahadevan   PetscErrorCode  ierr;
66*f6829af0SVijay Mahadevan 
67*f6829af0SVijay Mahadevan   PetscFunctionBegin;
68*f6829af0SVijay Mahadevan   PetscValidPointer(tag,2);
69*f6829af0SVijay Mahadevan 
70*f6829af0SVijay Mahadevan   /* Get the MOAB private data */
71*f6829af0SVijay Mahadevan   ierr = PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr);
72*f6829af0SVijay Mahadevan   ierr = PetscContainerGetPointer(moabdata, (void**) &vmoab);CHKERRQ(ierr);
73*f6829af0SVijay Mahadevan 
74*f6829af0SVijay Mahadevan   *tag = vmoab->tag;
75*f6829af0SVijay Mahadevan   PetscFunctionReturn(0);
76*f6829af0SVijay Mahadevan }
77*f6829af0SVijay Mahadevan 
78*f6829af0SVijay Mahadevan 
79*f6829af0SVijay Mahadevan #undef __FUNCT__
80*f6829af0SVijay Mahadevan #define __FUNCT__ "DMMoabGetVecRange"
81*f6829af0SVijay Mahadevan /*@
82*f6829af0SVijay Mahadevan   DMMoabGetVecRange - Get the MOAB entities associated with this Vec
83*f6829af0SVijay Mahadevan 
84*f6829af0SVijay Mahadevan   Input Parameter:
85*f6829af0SVijay Mahadevan . vec   - Vec being queried
86*f6829af0SVijay Mahadevan 
87*f6829af0SVijay Mahadevan   Output Parameter:
88*f6829af0SVijay Mahadevan . range - Entities associated with this Vec. NULL if vec is a native PETSc Vec.
89*f6829af0SVijay Mahadevan 
90*f6829af0SVijay Mahadevan   Level: beginner
91*f6829af0SVijay Mahadevan 
92*f6829af0SVijay Mahadevan .keywords: DMMoab, create
93*f6829af0SVijay Mahadevan @*/
94*f6829af0SVijay Mahadevan PetscErrorCode DMMoabGetVecRange(Vec vec,moab::Range *range)
95*f6829af0SVijay Mahadevan {
96*f6829af0SVijay Mahadevan   PetscContainer  moabdata;
97*f6829af0SVijay Mahadevan   Vec_MOAB        *vmoab;
98*f6829af0SVijay Mahadevan   PetscErrorCode  ierr;
99*f6829af0SVijay Mahadevan 
100*f6829af0SVijay Mahadevan   PetscFunctionBegin;
101*f6829af0SVijay Mahadevan   PetscValidPointer(range,2);
102*f6829af0SVijay Mahadevan 
103*f6829af0SVijay Mahadevan   /* Get the MOAB private data handle */
104*f6829af0SVijay Mahadevan   ierr = PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr);
105*f6829af0SVijay Mahadevan   ierr = PetscContainerGetPointer(moabdata, (void**) &vmoab);CHKERRQ(ierr);
106*f6829af0SVijay Mahadevan 
107*f6829af0SVijay Mahadevan   *range = *vmoab->tag_range;
108*f6829af0SVijay Mahadevan   PetscFunctionReturn(0);
109*f6829af0SVijay Mahadevan }
110*f6829af0SVijay Mahadevan 
111*f6829af0SVijay Mahadevan 
112*f6829af0SVijay Mahadevan #undef __FUNCT__
113*f6829af0SVijay Mahadevan #define __FUNCT__ "DMMoabVecGetArray"
114*f6829af0SVijay Mahadevan /*@
115*f6829af0SVijay Mahadevan   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
116*f6829af0SVijay Mahadevan 
117*f6829af0SVijay Mahadevan   Collective on MPI_Comm
118*f6829af0SVijay Mahadevan 
119*f6829af0SVijay Mahadevan   Input Parameter:
120*f6829af0SVijay Mahadevan . dm              - The DMMoab object being set
121*f6829af0SVijay Mahadevan . vec             - The Vector whose underlying data is requested
122*f6829af0SVijay Mahadevan 
123*f6829af0SVijay Mahadevan   Output Parameter:
124*f6829af0SVijay Mahadevan . array           - The local data array
125*f6829af0SVijay Mahadevan 
126*f6829af0SVijay Mahadevan   Level: intermediate
127*f6829af0SVijay Mahadevan 
128*f6829af0SVijay Mahadevan .keywords: MOAB, distributed array
129*f6829af0SVijay Mahadevan 
130*f6829af0SVijay Mahadevan .seealso: DMMoabVecRestoreArray(), DMMoabVecGetArrayRead(), DMMoabVecRestoreArrayRead()
131*f6829af0SVijay Mahadevan @*/
132*f6829af0SVijay Mahadevan PetscErrorCode  DMMoabVecGetArray(DM dm,Vec vec,void* array)
133*f6829af0SVijay Mahadevan {
134*f6829af0SVijay Mahadevan   DM_Moab        *dmmoab;
135*f6829af0SVijay Mahadevan   moab::ErrorCode merr;
136*f6829af0SVijay Mahadevan   PetscErrorCode  ierr;
137*f6829af0SVijay Mahadevan   PetscInt        count,i,f;
138*f6829af0SVijay Mahadevan   moab::Tag       vtag;
139*f6829af0SVijay Mahadevan   PetscScalar     **varray;
140*f6829af0SVijay Mahadevan   PetscScalar     *marray;
141*f6829af0SVijay Mahadevan   PetscContainer  moabdata;
142*f6829af0SVijay Mahadevan   Vec_MOAB        *vmoab,*xmoab;
143*f6829af0SVijay Mahadevan 
144*f6829af0SVijay Mahadevan   PetscFunctionBegin;
145*f6829af0SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
146*f6829af0SVijay Mahadevan   PetscValidHeaderSpecific(vec,VEC_CLASSID,2);
147*f6829af0SVijay Mahadevan   PetscValidPointer(array,3);
148*f6829af0SVijay Mahadevan   dmmoab=(DM_Moab*)dm->data;
149*f6829af0SVijay Mahadevan 
150*f6829af0SVijay Mahadevan   /* Get the Vec_MOAB struct for the original vector */
151*f6829af0SVijay Mahadevan   ierr = PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr);
152*f6829af0SVijay Mahadevan   ierr = PetscContainerGetPointer(moabdata, (void**)&vmoab);CHKERRQ(ierr);
153*f6829af0SVijay Mahadevan 
154*f6829af0SVijay Mahadevan   /* Get the real scalar array handle */
155*f6829af0SVijay Mahadevan   varray = reinterpret_cast<PetscScalar**>(array);
156*f6829af0SVijay Mahadevan 
157*f6829af0SVijay Mahadevan   if (vmoab->is_native_vec) {
158*f6829af0SVijay Mahadevan 
159*f6829af0SVijay Mahadevan     /* get the local representation of the arrays from Vectors */
160*f6829af0SVijay Mahadevan     ierr = DMGetLocalVector(dm,&vmoab->local);CHKERRQ(ierr);
161*f6829af0SVijay Mahadevan     ierr = DMGlobalToLocalBegin(dm,vec,INSERT_VALUES,vmoab->local);CHKERRQ(ierr);
162*f6829af0SVijay Mahadevan     ierr = DMGlobalToLocalEnd(dm,vec,INSERT_VALUES,vmoab->local);CHKERRQ(ierr);
163*f6829af0SVijay Mahadevan 
164*f6829af0SVijay Mahadevan     /* Get the Vec_MOAB struct for the original vector */
165*f6829af0SVijay Mahadevan     ierr = PetscObjectQuery((PetscObject)vmoab->local,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr);
166*f6829af0SVijay Mahadevan     ierr = PetscContainerGetPointer(moabdata, (void**)&xmoab);CHKERRQ(ierr);
167*f6829af0SVijay Mahadevan 
168*f6829af0SVijay Mahadevan     /* get the local representation of the arrays from Vectors */
169*f6829af0SVijay Mahadevan     ierr = VecGhostGetLocalForm(vmoab->local,&xmoab->local);CHKERRQ(ierr);
170*f6829af0SVijay Mahadevan     ierr = VecGhostUpdateBegin(vmoab->local,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
171*f6829af0SVijay Mahadevan     ierr = VecGhostUpdateEnd(vmoab->local,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
172*f6829af0SVijay Mahadevan 
173*f6829af0SVijay Mahadevan     ierr = VecGetArray(xmoab->local, varray);CHKERRQ(ierr);
174*f6829af0SVijay Mahadevan   }
175*f6829af0SVijay Mahadevan   else {
176*f6829af0SVijay Mahadevan 
177*f6829af0SVijay Mahadevan     /* Get the MOAB private data */
178*f6829af0SVijay Mahadevan     ierr = DMMoabGetVecTag(vec,&vtag);CHKERRQ(ierr);
179*f6829af0SVijay Mahadevan 
180*f6829af0SVijay Mahadevan     /* exchange the data into ghost cells first */
181*f6829af0SVijay Mahadevan     merr = dmmoab->pcomm->exchange_tags(vtag,*dmmoab->vlocal);MBERRNM(merr);
182*f6829af0SVijay Mahadevan 
183*f6829af0SVijay Mahadevan     ierr = PetscMalloc(sizeof(PetscScalar)*(dmmoab->nloc+dmmoab->nghost)*dmmoab->numFields,varray);CHKERRQ(ierr);
184*f6829af0SVijay Mahadevan 
185*f6829af0SVijay Mahadevan     /* Get the array data for local entities */
186*f6829af0SVijay Mahadevan     merr = dmmoab->mbiface->tag_iterate(vtag,dmmoab->vlocal->begin(),dmmoab->vlocal->end(),count,reinterpret_cast<void*&>(marray),false);MBERRNM(merr);
187*f6829af0SVijay Mahadevan     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);
188*f6829af0SVijay Mahadevan 
189*f6829af0SVijay Mahadevan     i=0;
190*f6829af0SVijay Mahadevan     for(moab::Range::iterator iter = dmmoab->vlocal->begin(); iter != dmmoab->vlocal->end(); iter++) {
191*f6829af0SVijay Mahadevan       for (f=0;f<dmmoab->numFields;f++,i++)
192*f6829af0SVijay Mahadevan         (*varray)[dmmoab->lidmap[(PetscInt)*iter]*dmmoab->numFields+f]=marray[i];
193*f6829af0SVijay Mahadevan     }
194*f6829af0SVijay Mahadevan   }
195*f6829af0SVijay Mahadevan   PetscFunctionReturn(0);
196*f6829af0SVijay Mahadevan }
197*f6829af0SVijay Mahadevan 
198*f6829af0SVijay Mahadevan 
199*f6829af0SVijay Mahadevan #undef __FUNCT__
200*f6829af0SVijay Mahadevan #define __FUNCT__ "DMMoabVecRestoreArray"
201*f6829af0SVijay Mahadevan /*@
202*f6829af0SVijay Mahadevan   DMMoabVecRestoreArray - Restores the writable direct access array obtained via DMMoabVecGetArray
203*f6829af0SVijay Mahadevan 
204*f6829af0SVijay Mahadevan   Collective on MPI_Comm
205*f6829af0SVijay Mahadevan 
206*f6829af0SVijay Mahadevan   Input Parameter:
207*f6829af0SVijay Mahadevan + dm              - The DMMoab object being set
208*f6829af0SVijay Mahadevan . vec             - The Vector whose underlying data is requested
209*f6829af0SVijay Mahadevan - array           - The local data array
210*f6829af0SVijay Mahadevan 
211*f6829af0SVijay Mahadevan   Level: intermediate
212*f6829af0SVijay Mahadevan 
213*f6829af0SVijay Mahadevan .keywords: MOAB, distributed array
214*f6829af0SVijay Mahadevan 
215*f6829af0SVijay Mahadevan .seealso: DMMoabVecGetArray(), DMMoabVecGetArrayRead(), DMMoabVecRestoreArrayRead()
216*f6829af0SVijay Mahadevan @*/
217*f6829af0SVijay Mahadevan PetscErrorCode  DMMoabVecRestoreArray(DM dm,Vec vec,void* array)
218*f6829af0SVijay Mahadevan {
219*f6829af0SVijay Mahadevan   DM_Moab        *dmmoab;
220*f6829af0SVijay Mahadevan   moab::ErrorCode merr;
221*f6829af0SVijay Mahadevan   PetscErrorCode  ierr;
222*f6829af0SVijay Mahadevan   moab::Tag       vtag;
223*f6829af0SVijay Mahadevan   PetscInt        count,i,f;
224*f6829af0SVijay Mahadevan   PetscScalar     **varray;
225*f6829af0SVijay Mahadevan   PetscScalar     *marray;
226*f6829af0SVijay Mahadevan   PetscContainer  moabdata;
227*f6829af0SVijay Mahadevan   Vec_MOAB        *vmoab,*xmoab;
228*f6829af0SVijay Mahadevan 
229*f6829af0SVijay Mahadevan   PetscFunctionBegin;
230*f6829af0SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
231*f6829af0SVijay Mahadevan   PetscValidHeaderSpecific(vec,VEC_CLASSID,2);
232*f6829af0SVijay Mahadevan   PetscValidPointer(array,3);
233*f6829af0SVijay Mahadevan   dmmoab=(DM_Moab*)dm->data;
234*f6829af0SVijay Mahadevan 
235*f6829af0SVijay Mahadevan   /* Get the Vec_MOAB struct for the original vector */
236*f6829af0SVijay Mahadevan   ierr = PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr);
237*f6829af0SVijay Mahadevan   ierr = PetscContainerGetPointer(moabdata, (void**)&vmoab);CHKERRQ(ierr);
238*f6829af0SVijay Mahadevan 
239*f6829af0SVijay Mahadevan   /* Get the real scalar array handle */
240*f6829af0SVijay Mahadevan   varray = reinterpret_cast<PetscScalar**>(array);
241*f6829af0SVijay Mahadevan 
242*f6829af0SVijay Mahadevan   if (vmoab->is_native_vec) {
243*f6829af0SVijay Mahadevan 
244*f6829af0SVijay Mahadevan     /* Get the Vec_MOAB struct for the original vector */
245*f6829af0SVijay Mahadevan     ierr = PetscObjectQuery((PetscObject)vmoab->local,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr);
246*f6829af0SVijay Mahadevan     ierr = PetscContainerGetPointer(moabdata, (void**)&xmoab);CHKERRQ(ierr);
247*f6829af0SVijay Mahadevan 
248*f6829af0SVijay Mahadevan     /* get the local representation of the arrays from Vectors */
249*f6829af0SVijay Mahadevan     ierr = VecRestoreArray(xmoab->local, varray);CHKERRQ(ierr);
250*f6829af0SVijay Mahadevan     ierr = VecGhostUpdateBegin(vmoab->local,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
251*f6829af0SVijay Mahadevan     ierr = VecGhostUpdateEnd(vmoab->local,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
252*f6829af0SVijay Mahadevan     ierr = VecGhostRestoreLocalForm(vmoab->local,&xmoab->local);CHKERRQ(ierr);
253*f6829af0SVijay Mahadevan 
254*f6829af0SVijay Mahadevan     /* restore local pieces */
255*f6829af0SVijay Mahadevan     ierr = DMLocalToGlobalBegin(dm,vmoab->local,INSERT_VALUES,vec);CHKERRQ(ierr);
256*f6829af0SVijay Mahadevan     ierr = DMLocalToGlobalEnd(dm,vmoab->local,INSERT_VALUES,vec);CHKERRQ(ierr);
257*f6829af0SVijay Mahadevan     ierr = DMRestoreLocalVector(dm, &vmoab->local);CHKERRQ(ierr);
258*f6829af0SVijay Mahadevan   }
259*f6829af0SVijay Mahadevan   else {
260*f6829af0SVijay Mahadevan 
261*f6829af0SVijay Mahadevan     /* Get the MOAB private data */
262*f6829af0SVijay Mahadevan     ierr = DMMoabGetVecTag(vec,&vtag);CHKERRQ(ierr);
263*f6829af0SVijay Mahadevan 
264*f6829af0SVijay Mahadevan     /* Get the array data for local entities */
265*f6829af0SVijay Mahadevan     merr = dmmoab->mbiface->tag_iterate(vtag,dmmoab->vlocal->begin(),dmmoab->vlocal->end(),count,reinterpret_cast<void*&>(marray),false);MBERRNM(merr);
266*f6829af0SVijay Mahadevan     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);
267*f6829af0SVijay Mahadevan 
268*f6829af0SVijay Mahadevan     i=0;
269*f6829af0SVijay Mahadevan     for(moab::Range::iterator iter = dmmoab->vlocal->begin(); iter != dmmoab->vlocal->end(); iter++) {
270*f6829af0SVijay Mahadevan       for (f=0;f<dmmoab->numFields;f++,i++)
271*f6829af0SVijay Mahadevan       marray[i] = (*varray)[dmmoab->lidmap[(PetscInt)*iter]*dmmoab->numFields+f];
272*f6829af0SVijay Mahadevan     }
273*f6829af0SVijay Mahadevan 
274*f6829af0SVijay Mahadevan     /* reduce the tags correctly -> should probably let the user choose how to reduce in the future
275*f6829af0SVijay Mahadevan       For all FEM residual based assembly calculations, MPI_SUM should serve well */
276*f6829af0SVijay Mahadevan     merr = dmmoab->pcomm->reduce_tags(vtag,MPI_SUM,*dmmoab->vlocal);MBERRV(dmmoab->mbiface,merr);
277*f6829af0SVijay Mahadevan 
278*f6829af0SVijay Mahadevan     ierr = PetscFree(*varray);CHKERRQ(ierr);
279*f6829af0SVijay Mahadevan   }
280*f6829af0SVijay Mahadevan   PetscFunctionReturn(0);
281*f6829af0SVijay Mahadevan }
282*f6829af0SVijay Mahadevan 
283*f6829af0SVijay Mahadevan #undef __FUNCT__
284*f6829af0SVijay Mahadevan #define __FUNCT__ "DMMoabVecGetArrayRead"
285*f6829af0SVijay Mahadevan /*@
286*f6829af0SVijay Mahadevan   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
287*f6829af0SVijay Mahadevan 
288*f6829af0SVijay Mahadevan   Collective on MPI_Comm
289*f6829af0SVijay Mahadevan 
290*f6829af0SVijay Mahadevan   Input Parameter:
291*f6829af0SVijay Mahadevan + dm              - The DMMoab object being set
292*f6829af0SVijay Mahadevan . vec             - The Vector whose underlying data is requested
293*f6829af0SVijay Mahadevan 
294*f6829af0SVijay Mahadevan   Output Parameter:
295*f6829af0SVijay Mahadevan . array           - The local data array
296*f6829af0SVijay Mahadevan 
297*f6829af0SVijay Mahadevan   Level: intermediate
298*f6829af0SVijay Mahadevan 
299*f6829af0SVijay Mahadevan .keywords: MOAB, distributed array
300*f6829af0SVijay Mahadevan 
301*f6829af0SVijay Mahadevan .seealso: DMMoabVecRestoreArrayRead(), DMMoabVecGetArray(), DMMoabVecRestoreArray()
302*f6829af0SVijay Mahadevan @*/
303*f6829af0SVijay Mahadevan PetscErrorCode  DMMoabVecGetArrayRead(DM dm,Vec vec,void* array)
304*f6829af0SVijay Mahadevan {
305*f6829af0SVijay Mahadevan   DM_Moab        *dmmoab;
306*f6829af0SVijay Mahadevan   moab::ErrorCode merr;
307*f6829af0SVijay Mahadevan   PetscErrorCode  ierr;
308*f6829af0SVijay Mahadevan   PetscInt        count,i,f;
309*f6829af0SVijay Mahadevan   moab::Tag       vtag;
310*f6829af0SVijay Mahadevan   PetscScalar     **varray;
311*f6829af0SVijay Mahadevan   PetscScalar     *marray;
312*f6829af0SVijay Mahadevan   PetscContainer  moabdata;
313*f6829af0SVijay Mahadevan   Vec_MOAB        *vmoab,*xmoab;
314*f6829af0SVijay Mahadevan 
315*f6829af0SVijay Mahadevan   PetscFunctionBegin;
316*f6829af0SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
317*f6829af0SVijay Mahadevan   PetscValidHeaderSpecific(vec,VEC_CLASSID,2);
318*f6829af0SVijay Mahadevan   PetscValidPointer(array,3);
319*f6829af0SVijay Mahadevan   dmmoab=(DM_Moab*)dm->data;
320*f6829af0SVijay Mahadevan 
321*f6829af0SVijay Mahadevan   /* Get the Vec_MOAB struct for the original vector */
322*f6829af0SVijay Mahadevan   ierr = PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr);
323*f6829af0SVijay Mahadevan   ierr = PetscContainerGetPointer(moabdata, (void**)&vmoab);CHKERRQ(ierr);
324*f6829af0SVijay Mahadevan 
325*f6829af0SVijay Mahadevan   /* Get the real scalar array handle */
326*f6829af0SVijay Mahadevan   varray = reinterpret_cast<PetscScalar**>(array);
327*f6829af0SVijay Mahadevan 
328*f6829af0SVijay Mahadevan   if (vmoab->is_native_vec) {
329*f6829af0SVijay Mahadevan     /* get the local representation of the arrays from Vectors */
330*f6829af0SVijay Mahadevan     ierr = DMGetLocalVector(dm,&vmoab->local);CHKERRQ(ierr);
331*f6829af0SVijay Mahadevan     ierr = DMGlobalToLocalBegin(dm,vec,INSERT_VALUES,vmoab->local);CHKERRQ(ierr);
332*f6829af0SVijay Mahadevan     ierr = DMGlobalToLocalEnd(dm,vec,INSERT_VALUES,vmoab->local);CHKERRQ(ierr);
333*f6829af0SVijay Mahadevan 
334*f6829af0SVijay Mahadevan     /* Get the Vec_MOAB struct for the original vector */
335*f6829af0SVijay Mahadevan     ierr = PetscObjectQuery((PetscObject)vmoab->local,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr);
336*f6829af0SVijay Mahadevan     ierr = PetscContainerGetPointer(moabdata, (void**)&xmoab);CHKERRQ(ierr);
337*f6829af0SVijay Mahadevan 
338*f6829af0SVijay Mahadevan     /* get the local representation of the arrays from Vectors */
339*f6829af0SVijay Mahadevan     ierr = VecGhostGetLocalForm(vmoab->local,&xmoab->local);CHKERRQ(ierr);
340*f6829af0SVijay Mahadevan     ierr = VecGhostUpdateBegin(vmoab->local,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
341*f6829af0SVijay Mahadevan     ierr = VecGhostUpdateEnd(vmoab->local,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
342*f6829af0SVijay Mahadevan     ierr = VecGetArray(xmoab->local, varray);CHKERRQ(ierr);
343*f6829af0SVijay Mahadevan   }
344*f6829af0SVijay Mahadevan   else {
345*f6829af0SVijay Mahadevan     /* Get the MOAB private data */
346*f6829af0SVijay Mahadevan     ierr = DMMoabGetVecTag(vec,&vtag);CHKERRQ(ierr);
347*f6829af0SVijay Mahadevan 
348*f6829af0SVijay Mahadevan     /* exchange the data into ghost cells first */
349*f6829af0SVijay Mahadevan     merr = dmmoab->pcomm->exchange_tags(vtag,*dmmoab->vlocal);MBERRNM(merr);
350*f6829af0SVijay Mahadevan 
351*f6829af0SVijay Mahadevan     ierr = PetscMalloc(sizeof(PetscScalar)*(dmmoab->nloc+dmmoab->nghost)*dmmoab->numFields,varray);CHKERRQ(ierr);
352*f6829af0SVijay Mahadevan 
353*f6829af0SVijay Mahadevan     /* Get the array data for local entities */
354*f6829af0SVijay Mahadevan     merr = dmmoab->mbiface->tag_iterate(vtag,dmmoab->vlocal->begin(),dmmoab->vlocal->end(),count,reinterpret_cast<void*&>(marray),false);MBERRNM(merr);
355*f6829af0SVijay Mahadevan     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);
356*f6829af0SVijay Mahadevan 
357*f6829af0SVijay Mahadevan     i=0;
358*f6829af0SVijay Mahadevan     for(moab::Range::iterator iter = dmmoab->vlocal->begin(); iter != dmmoab->vlocal->end(); iter++) {
359*f6829af0SVijay Mahadevan       for (f=0;f<dmmoab->numFields;f++,i++)
360*f6829af0SVijay Mahadevan         (*varray)[dmmoab->lidmap[(PetscInt)*iter]*dmmoab->numFields+f]=marray[i];
361*f6829af0SVijay Mahadevan     }
362*f6829af0SVijay Mahadevan   }
363*f6829af0SVijay Mahadevan   PetscFunctionReturn(0);
364*f6829af0SVijay Mahadevan }
365*f6829af0SVijay Mahadevan 
366*f6829af0SVijay Mahadevan 
367*f6829af0SVijay Mahadevan #undef __FUNCT__
368*f6829af0SVijay Mahadevan #define __FUNCT__ "DMMoabVecRestoreArrayRead"
369*f6829af0SVijay Mahadevan /*@
370*f6829af0SVijay Mahadevan   DMMoabVecRestoreArray - Restores the read-only direct access array obtained via DMMoabVecGetArray
371*f6829af0SVijay Mahadevan 
372*f6829af0SVijay Mahadevan   Collective on MPI_Comm
373*f6829af0SVijay Mahadevan 
374*f6829af0SVijay Mahadevan   Input Parameter:
375*f6829af0SVijay Mahadevan + dm              - The DMMoab object being set
376*f6829af0SVijay Mahadevan . vec             - The Vector whose underlying data is requested
377*f6829af0SVijay Mahadevan - array           - The local data array
378*f6829af0SVijay Mahadevan 
379*f6829af0SVijay Mahadevan   Level: intermediate
380*f6829af0SVijay Mahadevan 
381*f6829af0SVijay Mahadevan .keywords: MOAB, distributed array
382*f6829af0SVijay Mahadevan 
383*f6829af0SVijay Mahadevan .seealso: DMMoabVecGetArrayRead(), DMMoabVecGetArray(), DMMoabVecRestoreArray()
384*f6829af0SVijay Mahadevan @*/
385*f6829af0SVijay Mahadevan PetscErrorCode  DMMoabVecRestoreArrayRead(DM dm,Vec vec,void* array)
386*f6829af0SVijay Mahadevan {
387*f6829af0SVijay Mahadevan   PetscErrorCode  ierr;
388*f6829af0SVijay Mahadevan   PetscScalar     **varray;
389*f6829af0SVijay Mahadevan   PetscContainer  moabdata;
390*f6829af0SVijay Mahadevan   Vec_MOAB        *vmoab,*xmoab;
391*f6829af0SVijay Mahadevan 
392*f6829af0SVijay Mahadevan   PetscFunctionBegin;
393*f6829af0SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
394*f6829af0SVijay Mahadevan   PetscValidHeaderSpecific(vec,VEC_CLASSID,2);
395*f6829af0SVijay Mahadevan   PetscValidPointer(array,3);
396*f6829af0SVijay Mahadevan 
397*f6829af0SVijay Mahadevan   /* Get the Vec_MOAB struct for the original vector */
398*f6829af0SVijay Mahadevan   ierr = PetscObjectQuery((PetscObject)vec,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr);
399*f6829af0SVijay Mahadevan   ierr = PetscContainerGetPointer(moabdata, (void**)&vmoab);CHKERRQ(ierr);
400*f6829af0SVijay Mahadevan 
401*f6829af0SVijay Mahadevan   /* Get the real scalar array handle */
402*f6829af0SVijay Mahadevan   varray = reinterpret_cast<PetscScalar**>(array);
403*f6829af0SVijay Mahadevan 
404*f6829af0SVijay Mahadevan   if (vmoab->is_native_vec) {
405*f6829af0SVijay Mahadevan     /* Get the Vec_MOAB struct for the original vector */
406*f6829af0SVijay Mahadevan     ierr = PetscObjectQuery((PetscObject)vmoab->local,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr);
407*f6829af0SVijay Mahadevan     ierr = PetscContainerGetPointer(moabdata, (void**)&xmoab);CHKERRQ(ierr);
408*f6829af0SVijay Mahadevan 
409*f6829af0SVijay Mahadevan     /* restore the local representation of the arrays from Vectors */
410*f6829af0SVijay Mahadevan     ierr = VecRestoreArray(xmoab->local, varray);CHKERRQ(ierr);
411*f6829af0SVijay Mahadevan     ierr = VecGhostRestoreLocalForm(vmoab->local,&xmoab->local);CHKERRQ(ierr);
412*f6829af0SVijay Mahadevan 
413*f6829af0SVijay Mahadevan     /* restore local pieces */
414*f6829af0SVijay Mahadevan     ierr = DMRestoreLocalVector(dm, &vmoab->local);CHKERRQ(ierr);
415*f6829af0SVijay Mahadevan   }
416*f6829af0SVijay Mahadevan   else {
417*f6829af0SVijay Mahadevan     /* Nothing to do but just free the memory allocated before */
418*f6829af0SVijay Mahadevan     ierr = PetscFree(*varray);CHKERRQ(ierr);
419*f6829af0SVijay Mahadevan 
420*f6829af0SVijay Mahadevan   }
421*f6829af0SVijay Mahadevan   PetscFunctionReturn(0);
422*f6829af0SVijay Mahadevan }
423*f6829af0SVijay Mahadevan 
424*f6829af0SVijay Mahadevan 
425*f6829af0SVijay Mahadevan #undef __FUNCT__
426*f6829af0SVijay Mahadevan #define __FUNCT__ "DMCreateVector_Moab_Private"
427*f6829af0SVijay Mahadevan PetscErrorCode DMCreateVector_Moab_Private(DM dm,moab::Tag tag,const moab::Range* userrange,PetscBool is_global_vec,PetscBool destroy_tag,Vec *vec)
428032b8ab6SVijay Mahadevan {
429032b8ab6SVijay Mahadevan   PetscErrorCode         ierr;
430032b8ab6SVijay Mahadevan   moab::ErrorCode        merr;
431032b8ab6SVijay Mahadevan   PetscBool              is_newtag;
432351b8a77SVijay Mahadevan   const moab::Range      *range;
433bb8f3634SVijay Mahadevan   PetscInt               count,lnative_vec,gnative_vec;
4349088682fSVijay Mahadevan   std::string ttname;
435032b8ab6SVijay Mahadevan   PetscScalar *data_ptr;
436bb8f3634SVijay Mahadevan   ISLocalToGlobalMapping ltogb;
437032b8ab6SVijay Mahadevan 
438db66d124SVijay Mahadevan   Vec_MOAB *vmoab;
439032b8ab6SVijay Mahadevan   DM_Moab *dmmoab = (DM_Moab*)dm->data;
440032b8ab6SVijay Mahadevan   moab::ParallelComm *pcomm = dmmoab->pcomm;
441032b8ab6SVijay Mahadevan   moab::Interface *mbiface = dmmoab->mbiface;
442032b8ab6SVijay Mahadevan 
443032b8ab6SVijay Mahadevan   PetscFunctionBegin;
444032b8ab6SVijay Mahadevan   if(!userrange) range = dmmoab->vowned;
445032b8ab6SVijay Mahadevan   else range = userrange;
446db66d124SVijay Mahadevan   if(!range) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Input range cannot be empty or call DMSetUp first.");
447032b8ab6SVijay Mahadevan 
448bb8f3634SVijay Mahadevan   /* If the tag data is in a single sequence, use PETSc native vector since tag_iterate isn't useful anymore */
449bb8f3634SVijay Mahadevan   lnative_vec=(range->psize()-1);
450bb8f3634SVijay Mahadevan 
451*f6829af0SVijay Mahadevan #if 0
4526465f021SVijay Mahadevan //  lnative_vec=1; /* NOTE: Testing PETSc vector: will force to create native vector all the time */
4536465f021SVijay Mahadevan //  lnative_vec=0; /* NOTE: Testing MOAB vector: will force to create MOAB tag_iterate based vector all the time */
454*f6829af0SVijay Mahadevan #endif
455bb8f3634SVijay Mahadevan   ierr = MPI_Allreduce(&lnative_vec, &gnative_vec, 1, MPI_INT, MPI_MAX, pcomm->comm());CHKERRQ(ierr);
456bb8f3634SVijay Mahadevan 
457bb8f3634SVijay Mahadevan   /* Create the MOAB internal data object */
458bb8f3634SVijay Mahadevan   ierr = PetscNew(Vec_MOAB,&vmoab);CHKERRQ(ierr);
459bb8f3634SVijay Mahadevan   vmoab->is_native_vec=(gnative_vec>0?PETSC_TRUE:PETSC_FALSE);
460bb8f3634SVijay Mahadevan 
461bb8f3634SVijay Mahadevan   if (!vmoab->is_native_vec) {
4629088682fSVijay Mahadevan     merr = mbiface->tag_get_name(tag, ttname);
4639088682fSVijay Mahadevan     if (!ttname.length() && merr !=moab::MB_SUCCESS) {
464db66d124SVijay Mahadevan       /* get the new name for the anonymous MOABVec -> the tag_name will be destroyed along with Tag */
465032b8ab6SVijay Mahadevan       char *tag_name = PETSC_NULL;
466*f6829af0SVijay Mahadevan       ierr = DMVecCreateTagName_Moab_Private(pcomm,&tag_name);CHKERRQ(ierr);
467032b8ab6SVijay Mahadevan       is_newtag = PETSC_TRUE;
468032b8ab6SVijay Mahadevan 
469db66d124SVijay Mahadevan       /* Create the default value for the tag (all zeros) */
470addae81cSVijay Mahadevan       std::vector<PetscScalar> default_value(dmmoab->numFields, 0.0);
471032b8ab6SVijay Mahadevan 
472db66d124SVijay Mahadevan       /* Create the tag */
473addae81cSVijay Mahadevan       merr = mbiface->tag_get_handle(tag_name,dmmoab->numFields,moab::MB_TYPE_DOUBLE,tag,
474032b8ab6SVijay Mahadevan                                     moab::MB_TAG_DENSE|moab::MB_TAG_CREAT,default_value.data());MBERRNM(merr);
475032b8ab6SVijay Mahadevan       ierr = PetscFree(tag_name);CHKERRQ(ierr);
476032b8ab6SVijay Mahadevan     }
477032b8ab6SVijay Mahadevan     else {
478db66d124SVijay Mahadevan       /* Make sure the tag data is of type "double" */
479032b8ab6SVijay Mahadevan       moab::DataType tag_type;
480032b8ab6SVijay Mahadevan       merr = mbiface->tag_get_data_type(tag, tag_type);MBERRNM(merr);
481032b8ab6SVijay Mahadevan       if(tag_type != moab::MB_TYPE_DOUBLE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Tag data type must be MB_TYPE_DOUBLE");
482032b8ab6SVijay Mahadevan       is_newtag = destroy_tag;
483032b8ab6SVijay Mahadevan     }
484032b8ab6SVijay Mahadevan 
485032b8ab6SVijay Mahadevan     vmoab->tag = tag;
486bb8f3634SVijay Mahadevan     vmoab->new_tag = is_newtag;
487bb8f3634SVijay Mahadevan   }
488032b8ab6SVijay Mahadevan   vmoab->mbiface = mbiface;
489032b8ab6SVijay Mahadevan   vmoab->pcomm = pcomm;
490032b8ab6SVijay Mahadevan   vmoab->is_global_vec = is_global_vec;
491bb8f3634SVijay Mahadevan   vmoab->tag_size=dmmoab->bs;
492032b8ab6SVijay Mahadevan 
493bb8f3634SVijay Mahadevan   if (vmoab->is_native_vec) {
494fc418013SVijay Mahadevan 
495bb8f3634SVijay Mahadevan     /* Create the PETSc Vector directly and attach our functions accordingly */
496bb8f3634SVijay Mahadevan     if (!is_global_vec) {
497bb8f3634SVijay Mahadevan       /* This is an MPI Vector with ghosted padding */
498addae81cSVijay Mahadevan       ierr = VecCreateGhostBlock(pcomm->comm(),dmmoab->bs,dmmoab->numFields*dmmoab->nloc,
499addae81cSVijay Mahadevan                                  dmmoab->numFields*dmmoab->n,dmmoab->nghost,&dmmoab->gsindices[dmmoab->nloc],vec);CHKERRQ(ierr);
500bb8f3634SVijay Mahadevan     }
501bb8f3634SVijay Mahadevan     else {
502bb8f3634SVijay Mahadevan       /* This is an MPI/SEQ Vector */
503bb8f3634SVijay Mahadevan       ierr = VecCreate(pcomm->comm(),vec);CHKERRQ(ierr);
504addae81cSVijay Mahadevan       ierr = VecSetSizes(*vec,dmmoab->numFields*dmmoab->nloc,PETSC_DECIDE);CHKERRQ(ierr);
505bb8f3634SVijay Mahadevan       ierr = VecSetBlockSize(*vec,dmmoab->bs);CHKERRQ(ierr);
506bb8f3634SVijay Mahadevan       ierr = VecSetType(*vec, VECMPI);CHKERRQ(ierr);
507bb8f3634SVijay Mahadevan     }
508bb8f3634SVijay Mahadevan   }
509bb8f3634SVijay Mahadevan   else {
510db66d124SVijay Mahadevan     /* Call tag_iterate. This will cause MOAB to allocate memory for the
511db66d124SVijay Mahadevan        tag data if it hasn't already happened */
512032b8ab6SVijay Mahadevan     merr = mbiface->tag_iterate(tag,range->begin(),range->end(),count,(void*&)data_ptr);MBERRNM(merr);
513032b8ab6SVijay Mahadevan 
514bb8f3634SVijay Mahadevan     /* set the reference for vector range */
515bb8f3634SVijay Mahadevan     vmoab->tag_range = new moab::Range(*range);
516addae81cSVijay Mahadevan     merr = mbiface->tag_get_length(tag,dmmoab->numFields);MBERRNM(merr);
517032b8ab6SVijay Mahadevan 
518db66d124SVijay Mahadevan     /* Create the PETSc Vector
519db66d124SVijay Mahadevan       Query MOAB mesh to check if there are any ghosted entities
520032b8ab6SVijay Mahadevan         -> if we do, create a ghosted vector to map correctly to the same layout
521032b8ab6SVijay Mahadevan         -> else, create a non-ghosted parallel vector */
522bb8f3634SVijay Mahadevan     if (!is_global_vec) {
523db66d124SVijay Mahadevan       /* This is an MPI Vector with ghosted padding */
524addae81cSVijay Mahadevan       ierr = VecCreateGhostBlockWithArray(pcomm->comm(),dmmoab->bs,dmmoab->numFields*dmmoab->nloc,
525addae81cSVijay Mahadevan                                 dmmoab->numFields*dmmoab->n,dmmoab->nghost,&dmmoab->gsindices[dmmoab->nloc],data_ptr,vec);CHKERRQ(ierr);
526032b8ab6SVijay Mahadevan     }
527032b8ab6SVijay Mahadevan     else {
528bb8f3634SVijay Mahadevan       /* This is an MPI Vector without ghosted padding */
529addae81cSVijay Mahadevan       ierr = VecCreateMPIWithArray(pcomm->comm(),dmmoab->bs,dmmoab->numFields*range->size(),
530bb8f3634SVijay Mahadevan                                 PETSC_DECIDE,data_ptr,vec);CHKERRQ(ierr);
531032b8ab6SVijay Mahadevan     }
532bb8f3634SVijay Mahadevan   }
533bb8f3634SVijay Mahadevan   ierr = VecSetFromOptions(*vec);CHKERRQ(ierr);
534032b8ab6SVijay Mahadevan 
535*f6829af0SVijay Mahadevan   /* create a container and store the internal MOAB data for faster access based on Entities etc */
536032b8ab6SVijay Mahadevan   PetscContainer moabdata;
537bb8f3634SVijay Mahadevan   ierr = PetscContainerCreate(PETSC_COMM_WORLD,&moabdata);CHKERRQ(ierr);
538032b8ab6SVijay Mahadevan   ierr = PetscContainerSetPointer(moabdata,vmoab);CHKERRQ(ierr);
539*f6829af0SVijay Mahadevan   ierr = PetscContainerSetUserDestroy(moabdata,DMVecUserDestroy_Moab);CHKERRQ(ierr);
540032b8ab6SVijay Mahadevan   ierr = PetscObjectCompose((PetscObject)*vec,"MOABData",(PetscObject)moabdata);CHKERRQ(ierr);
541*f6829af0SVijay Mahadevan   (*vec)->ops->duplicate = DMVecDuplicate_Moab;
542032b8ab6SVijay Mahadevan   ierr = PetscContainerDestroy(&moabdata);CHKERRQ(ierr);
543032b8ab6SVijay Mahadevan 
544db66d124SVijay Mahadevan   /* Vector created, manually set local to global mapping */
5456d9eb265SVijay Mahadevan   if (dmmoab->ltog_map) {
546bb8f3634SVijay Mahadevan     ierr = VecSetLocalToGlobalMapping(*vec,dmmoab->ltog_map);CHKERRQ(ierr);
547bb8f3634SVijay Mahadevan     ierr = ISLocalToGlobalMappingBlock(dmmoab->ltog_map,dmmoab->bs,&ltogb);
548bb8f3634SVijay Mahadevan     ierr = VecSetLocalToGlobalMappingBlock(*vec,ltogb);CHKERRQ(ierr);
549bb8f3634SVijay Mahadevan     ierr = ISLocalToGlobalMappingDestroy(&ltogb);CHKERRQ(ierr);
5506d9eb265SVijay Mahadevan   }
551032b8ab6SVijay Mahadevan 
552032b8ab6SVijay Mahadevan   /* set the DM reference to the vector */
553032b8ab6SVijay Mahadevan   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
554032b8ab6SVijay Mahadevan   PetscFunctionReturn(0);
555032b8ab6SVijay Mahadevan }
556032b8ab6SVijay Mahadevan 
557032b8ab6SVijay Mahadevan 
558032b8ab6SVijay Mahadevan #undef __FUNCT__
559*f6829af0SVijay Mahadevan #define __FUNCT__ "DMVecCreateTagName_Moab_Private"
560*f6829af0SVijay Mahadevan /*  DMVecCreateTagName_Moab_Private
561032b8ab6SVijay Mahadevan  *
562032b8ab6SVijay Mahadevan  *  Creates a unique tag name that will be shared across processes. If
563032b8ab6SVijay Mahadevan  *  pcomm is NULL, then this is a serial vector. A unique tag name
564032b8ab6SVijay Mahadevan  *  will be returned in tag_name in either case.
565032b8ab6SVijay Mahadevan  *
566032b8ab6SVijay Mahadevan  *  The tag names have the format _PETSC_VEC_N where N is some integer.
567032b8ab6SVijay Mahadevan  *
568032b8ab6SVijay Mahadevan  *  NOTE: The tag_name is allocated in this routine; The user needs to free
569032b8ab6SVijay Mahadevan  *        the character array.
570032b8ab6SVijay Mahadevan  */
571*f6829af0SVijay Mahadevan PetscErrorCode DMVecCreateTagName_Moab_Private(moab::ParallelComm *pcomm,char** tag_name)
572032b8ab6SVijay Mahadevan {
573032b8ab6SVijay Mahadevan   moab::ErrorCode mberr;
574032b8ab6SVijay Mahadevan   PetscErrorCode  ierr;
575032b8ab6SVijay Mahadevan   PetscInt        n,global_n;
576032b8ab6SVijay Mahadevan   moab::Tag indexTag;
577032b8ab6SVijay Mahadevan 
578032b8ab6SVijay Mahadevan   PetscFunctionBegin;
579032b8ab6SVijay Mahadevan   const char*       PVEC_PREFIX      = "__PETSC_VEC_";
580032b8ab6SVijay Mahadevan   ierr = PetscMalloc(PETSC_MAX_PATH_LEN, tag_name);CHKERRQ(ierr);
581032b8ab6SVijay Mahadevan 
582efd17f3eSVijay Mahadevan   /* Check to see if there are any PETSc vectors defined */
583032b8ab6SVijay Mahadevan   moab::Interface  *mbiface = pcomm->get_moab();
584032b8ab6SVijay Mahadevan   moab::EntityHandle rootset = mbiface->get_root_set();
585032b8ab6SVijay Mahadevan 
586efd17f3eSVijay Mahadevan   /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
587032b8ab6SVijay Mahadevan   mberr = mbiface->tag_get_handle("__PETSC_VECS__",1,moab::MB_TYPE_INTEGER,indexTag,
588032b8ab6SVijay Mahadevan                                   moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT,0);MBERRNM(mberr);
589032b8ab6SVijay Mahadevan   mberr = mbiface->tag_get_data(indexTag, &rootset, 1, &n);
590032b8ab6SVijay Mahadevan   if (mberr == moab::MB_TAG_NOT_FOUND) n=0;  /* this is the first temporary vector */
591032b8ab6SVijay Mahadevan   else MBERRNM(mberr);
592032b8ab6SVijay Mahadevan 
593032b8ab6SVijay Mahadevan   /* increment the new value of n */
594032b8ab6SVijay Mahadevan   ++n;
595032b8ab6SVijay Mahadevan 
596efd17f3eSVijay Mahadevan   /* Make sure that n is consistent across all processes */
597032b8ab6SVijay Mahadevan   ierr = MPI_Allreduce(&n,&global_n,1,MPI_INT,MPI_MAX,pcomm->comm());CHKERRQ(ierr);
598032b8ab6SVijay Mahadevan 
599efd17f3eSVijay Mahadevan   /* Set the new name accordingly and return */
600032b8ab6SVijay Mahadevan   ierr = PetscSNPrintf(*tag_name, PETSC_MAX_PATH_LEN-1, "%s_%D", PVEC_PREFIX, global_n);CHKERRQ(ierr);
601032b8ab6SVijay Mahadevan   mberr = mbiface->tag_set_data(indexTag, &rootset, 1, (const void*)&global_n);MBERRNM(mberr);
602032b8ab6SVijay Mahadevan   PetscFunctionReturn(0);
603032b8ab6SVijay Mahadevan }
604032b8ab6SVijay Mahadevan 
605032b8ab6SVijay Mahadevan 
606032b8ab6SVijay Mahadevan #undef __FUNCT__
607*f6829af0SVijay Mahadevan #define __FUNCT__ "DMCreateGlobalVector_Moab"
608*f6829af0SVijay Mahadevan PetscErrorCode DMCreateGlobalVector_Moab(DM dm,Vec *gvec)
609*f6829af0SVijay Mahadevan {
610*f6829af0SVijay Mahadevan   PetscErrorCode  ierr;
611*f6829af0SVijay Mahadevan   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
612*f6829af0SVijay Mahadevan 
613*f6829af0SVijay Mahadevan   PetscFunctionBegin;
614*f6829af0SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
615*f6829af0SVijay Mahadevan   PetscValidPointer(gvec,2);
616*f6829af0SVijay Mahadevan   ierr = DMCreateVector_Moab_Private(dm,PETSC_NULL,dmmoab->vowned,PETSC_TRUE,PETSC_TRUE,gvec);CHKERRQ(ierr);
617*f6829af0SVijay Mahadevan   PetscFunctionReturn(0);
618*f6829af0SVijay Mahadevan }
619*f6829af0SVijay Mahadevan 
620*f6829af0SVijay Mahadevan 
621*f6829af0SVijay Mahadevan #undef __FUNCT__
622*f6829af0SVijay Mahadevan #define __FUNCT__ "DMCreateLocalVector_Moab"
623*f6829af0SVijay Mahadevan PetscErrorCode DMCreateLocalVector_Moab(DM dm,Vec *lvec)
624*f6829af0SVijay Mahadevan {
625*f6829af0SVijay Mahadevan   PetscErrorCode  ierr;
626*f6829af0SVijay Mahadevan   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
627*f6829af0SVijay Mahadevan 
628*f6829af0SVijay Mahadevan   PetscFunctionBegin;
629*f6829af0SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
630*f6829af0SVijay Mahadevan   PetscValidPointer(lvec,2);
631*f6829af0SVijay Mahadevan   ierr = DMCreateVector_Moab_Private(dm,PETSC_NULL,dmmoab->vlocal,PETSC_FALSE,PETSC_TRUE,lvec);CHKERRQ(ierr);
632*f6829af0SVijay Mahadevan   PetscFunctionReturn(0);
633*f6829af0SVijay Mahadevan }
634*f6829af0SVijay Mahadevan 
635*f6829af0SVijay Mahadevan 
636*f6829af0SVijay Mahadevan #undef __FUNCT__
637*f6829af0SVijay Mahadevan #define __FUNCT__ "DMVecDuplicate_Moab"
638*f6829af0SVijay Mahadevan PetscErrorCode DMVecDuplicate_Moab(Vec x,Vec *y)
639*f6829af0SVijay Mahadevan {
640*f6829af0SVijay Mahadevan   PetscErrorCode ierr;
641*f6829af0SVijay Mahadevan   DM             dm;
642*f6829af0SVijay Mahadevan   PetscContainer  moabdata;
643*f6829af0SVijay Mahadevan   Vec_MOAB        *vmoab;
644*f6829af0SVijay Mahadevan 
645*f6829af0SVijay Mahadevan   PetscFunctionBegin;
646*f6829af0SVijay Mahadevan   PetscValidHeaderSpecific(x,VEC_CLASSID,1);
647*f6829af0SVijay Mahadevan   PetscValidPointer(y,2);
648*f6829af0SVijay Mahadevan 
649*f6829af0SVijay Mahadevan   /* Get the Vec_MOAB struct for the original vector */
650*f6829af0SVijay Mahadevan   ierr = PetscObjectQuery((PetscObject)x,"MOABData", (PetscObject*) &moabdata);CHKERRQ(ierr);
651*f6829af0SVijay Mahadevan   ierr = PetscContainerGetPointer(moabdata, (void**)&vmoab);CHKERRQ(ierr);
652*f6829af0SVijay Mahadevan 
653*f6829af0SVijay Mahadevan   ierr = VecGetDM(x, &dm);CHKERRQ(ierr);
654*f6829af0SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
655*f6829af0SVijay Mahadevan 
656*f6829af0SVijay Mahadevan   ierr = DMCreateVector_Moab_Private(dm,PETSC_NULL,vmoab->tag_range,vmoab->is_global_vec,PETSC_TRUE,y);CHKERRQ(ierr);
657*f6829af0SVijay Mahadevan   ierr = VecSetDM(*y, dm);CHKERRQ(ierr);
658*f6829af0SVijay Mahadevan   PetscFunctionReturn(0);
659*f6829af0SVijay Mahadevan }
660*f6829af0SVijay Mahadevan 
661*f6829af0SVijay Mahadevan 
662*f6829af0SVijay Mahadevan #undef __FUNCT__
663*f6829af0SVijay Mahadevan #define __FUNCT__ "DMVecUserDestroy_Moab"
664*f6829af0SVijay Mahadevan PetscErrorCode DMVecUserDestroy_Moab(void *user)
665032b8ab6SVijay Mahadevan {
666efd17f3eSVijay Mahadevan   Vec_MOAB        *vmoab=(Vec_MOAB*)user;
667032b8ab6SVijay Mahadevan   PetscErrorCode  ierr;
668032b8ab6SVijay Mahadevan   moab::ErrorCode merr;
669032b8ab6SVijay Mahadevan 
670032b8ab6SVijay Mahadevan   PetscFunctionBegin;
671032b8ab6SVijay Mahadevan   if(vmoab->new_tag && vmoab->tag) {
672efd17f3eSVijay Mahadevan     /* Tag was created via a call to VecDuplicate, delete the underlying tag in MOAB */
673032b8ab6SVijay Mahadevan     merr = vmoab->mbiface->tag_delete(vmoab->tag);MBERRNM(merr);
674032b8ab6SVijay Mahadevan   }
675032b8ab6SVijay Mahadevan   delete vmoab->tag_range;
676032b8ab6SVijay Mahadevan   vmoab->tag = PETSC_NULL;
677032b8ab6SVijay Mahadevan   vmoab->mbiface = PETSC_NULL;
678032b8ab6SVijay Mahadevan   vmoab->pcomm = PETSC_NULL;
679032b8ab6SVijay Mahadevan   ierr = PetscFree(vmoab);CHKERRQ(ierr);
680032b8ab6SVijay Mahadevan   PetscFunctionReturn(0);
681032b8ab6SVijay Mahadevan }
682032b8ab6SVijay Mahadevan 
683efd17f3eSVijay Mahadevan 
684efd17f3eSVijay Mahadevan #undef __FUNCT__
685db66d124SVijay Mahadevan #define __FUNCT__ "DMGlobalToLocalBegin_Moab"
686db66d124SVijay Mahadevan PetscErrorCode  DMGlobalToLocalBegin_Moab(DM dm,Vec g,InsertMode mode,Vec l)
687db66d124SVijay Mahadevan {
688db66d124SVijay Mahadevan   PetscErrorCode    ierr;
689db66d124SVijay Mahadevan   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
690db66d124SVijay Mahadevan 
691db66d124SVijay Mahadevan   PetscFunctionBegin;
6926d9eb265SVijay Mahadevan   ierr = VecScatterBegin(dmmoab->ltog_sendrecv,g,l,mode,SCATTER_REVERSE);CHKERRQ(ierr);
693db66d124SVijay Mahadevan   PetscFunctionReturn(0);
694db66d124SVijay Mahadevan }
695db66d124SVijay Mahadevan 
696db66d124SVijay Mahadevan 
697db66d124SVijay Mahadevan #undef __FUNCT__
698db66d124SVijay Mahadevan #define __FUNCT__ "DMGlobalToLocalEnd_Moab"
699db66d124SVijay Mahadevan PetscErrorCode  DMGlobalToLocalEnd_Moab(DM dm,Vec g,InsertMode mode,Vec l)
700db66d124SVijay Mahadevan {
701db66d124SVijay Mahadevan   PetscErrorCode    ierr;
702db66d124SVijay Mahadevan   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
703db66d124SVijay Mahadevan 
704db66d124SVijay Mahadevan   PetscFunctionBegin;
7056d9eb265SVijay Mahadevan   ierr = VecScatterEnd(dmmoab->ltog_sendrecv,g,l,mode,SCATTER_REVERSE);CHKERRQ(ierr);
706db66d124SVijay Mahadevan   PetscFunctionReturn(0);
707db66d124SVijay Mahadevan }
708db66d124SVijay Mahadevan 
709db66d124SVijay Mahadevan 
710db66d124SVijay Mahadevan #undef __FUNCT__
711db66d124SVijay Mahadevan #define __FUNCT__ "DMLocalToGlobalBegin_Moab"
712db66d124SVijay Mahadevan PetscErrorCode  DMLocalToGlobalBegin_Moab(DM dm,Vec l,InsertMode mode,Vec g)
713db66d124SVijay Mahadevan {
714db66d124SVijay Mahadevan   PetscErrorCode    ierr;
715db66d124SVijay Mahadevan   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
716db66d124SVijay Mahadevan 
717db66d124SVijay Mahadevan   PetscFunctionBegin;
7186d9eb265SVijay Mahadevan   ierr = VecScatterBegin(dmmoab->ltog_sendrecv,l,g,mode,SCATTER_FORWARD);CHKERRQ(ierr);
719db66d124SVijay Mahadevan   PetscFunctionReturn(0);
720db66d124SVijay Mahadevan }
721db66d124SVijay Mahadevan 
722db66d124SVijay Mahadevan 
723db66d124SVijay Mahadevan #undef __FUNCT__
724db66d124SVijay Mahadevan #define __FUNCT__ "DMLocalToGlobalEnd_Moab"
725db66d124SVijay Mahadevan PetscErrorCode  DMLocalToGlobalEnd_Moab(DM dm,Vec l,InsertMode mode,Vec g)
726db66d124SVijay Mahadevan {
727db66d124SVijay Mahadevan   PetscErrorCode    ierr;
728db66d124SVijay Mahadevan   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
729db66d124SVijay Mahadevan 
730db66d124SVijay Mahadevan   PetscFunctionBegin;
7316d9eb265SVijay Mahadevan   ierr = VecScatterEnd(dmmoab->ltog_sendrecv,l,g,mode,SCATTER_FORWARD);CHKERRQ(ierr);
732db66d124SVijay Mahadevan   PetscFunctionReturn(0);
733db66d124SVijay Mahadevan }
734db66d124SVijay Mahadevan 
735