Actual source code: commonmpvec.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: #include <../src/vec/vec/impls/mpi/pvecimpl.h>   /*I  "petscvec.h"   I*/

  6: /*
  7:   This is used in VecGhostGetLocalForm and VecGhostRestoreLocalForm to ensure
  8:   that the state is updated if either vector has changed since the last time
  9:   one of these functions was called.  It could apply to any PetscObject, but
 10:   VecGhost is quite different from other objects in that two separate vectors
 11:   look at the same memory.

 13:   In principle, we could only propagate state to the local vector on
 14:   GetLocalForm and to the global vector on RestoreLocalForm, but this version is
 15:   more conservative (i.e. robust against misuse) and simpler.

 17:   Note that this function is correct and changes nothing if both arguments are the
 18:   same, which is the case in serial.
 19: */
 20: static PetscErrorCode VecGhostStateSync_Private(Vec g,Vec l)
 21: {
 22:   PetscErrorCode   ierr;
 23:   PetscObjectState gstate,lstate;

 26:   PetscObjectStateGet((PetscObject)g,&gstate);
 27:   PetscObjectStateGet((PetscObject)l,&lstate);
 28:   PetscObjectStateSet((PetscObject)g,PetscMax(gstate,lstate));
 29:   PetscObjectStateSet((PetscObject)l,PetscMax(gstate,lstate));
 30:   return(0);
 31: }


 36: /*@
 37:     VecGhostGetLocalForm - Obtains the local ghosted representation of
 38:     a parallel vector (obtained with VecCreateGhost(), VecCreateGhostWithArray()
 39:     or VecCreateSeq()). Returns NULL if the Vec is not ghosted.

 41:     Logically Collective

 43:     Input Parameter:
 44: .   g - the global vector

 46:     Output Parameter:
 47: .   l - the local (ghosted) representation, NULL if g is not ghosted

 49:     Notes:
 50:     This routine does not actually update the ghost values, but rather it
 51:     returns a sequential vector that includes the locations for the ghost
 52:     values and their current values. The returned vector and the original
 53:     vector passed in share the same array that contains the actual vector data.

 55:     To update the ghost values from the locations on the other processes one must call
 56:     VecGhostUpdateBegin() and VecGhostUpdateEnd() before accessing the ghost values. Thus normal
 57:     usage is
 58: $     VecGhostUpdateBegin(x,INSERT_VALUES,SCATTER_FORWARD);
 59: $     VecGhostUpdateEnd(x,INSERT_VALUES,SCATTER_FORWARD);
 60: $     VecGhostGetLocalForm(x,&xlocal);
 61: $     VecGetArray(xlocal,&xvalues);
 62: $        // access the non-ghost values in locations xvalues[0:n-1] and ghost values in locations xvalues[n:n+nghost];
 63: $     VecRestoreArray(xlocal,&xvalues);
 64: $     VecGhostRestoreLocalForm(x,&xlocal);

 66:     One should call VecGhostRestoreLocalForm() or VecDestroy() once one is
 67:     finished using the object.

 69:     Level: advanced

 71:    Concepts: vectors^ghost point access

 73: .seealso: VecCreateGhost(), VecGhostRestoreLocalForm(), VecCreateGhostWithArray()

 75: @*/
 76: PetscErrorCode  VecGhostGetLocalForm(Vec g,Vec *l)
 77: {
 79:   PetscBool      isseq,ismpi;


 85:   PetscObjectTypeCompare((PetscObject)g,VECSEQ,&isseq);
 86:   PetscObjectTypeCompare((PetscObject)g,VECMPI,&ismpi);
 87:   if (ismpi) {
 88:     Vec_MPI *v = (Vec_MPI*)g->data;
 89:     *l = v->localrep;
 90:   } else if (isseq) {
 91:     *l = g;
 92:   } else {
 93:     *l = NULL;
 94:   }
 95:   if (*l) {
 96:     VecGhostStateSync_Private(g,*l);
 97:     PetscObjectReference((PetscObject)*l);
 98:   }
 99:   return(0);
100: }

104: /*@
105:     VecGhostIsLocalForm - Checks if a given vector is the local form of a global vector

107:     Not Collective

109:     Input Parameter:
110: +   g - the global vector
111: -   l - the local vector

113:     Output Parameter:
114: .   flg - PETSC_TRUE if local vector is local form

116:     Level: advanced

118:    Concepts: vectors^ghost point access

120: .seealso: VecCreateGhost(), VecGhostRestoreLocalForm(), VecCreateGhostWithArray(), VecGhostGetLocalForm()

122: @*/
123: PetscErrorCode VecGhostIsLocalForm(Vec g,Vec l,PetscBool *flg)
124: {
126:   PetscBool      isseq,ismpi;


132:   *flg = PETSC_FALSE;
133:   PetscObjectTypeCompare((PetscObject)g,VECSEQ,&isseq);
134:   PetscObjectTypeCompare((PetscObject)g,VECMPI,&ismpi);
135:   if (ismpi) {
136:     Vec_MPI *v = (Vec_MPI*)g->data;
137:     if (l == v->localrep) *flg = PETSC_TRUE;
138:   } else if (isseq) {
139:     if (l == g) *flg = PETSC_TRUE;
140:   } else SETERRQ(PetscObjectComm((PetscObject)g),PETSC_ERR_ARG_WRONG,"Global vector is not ghosted");
141:   return(0);
142: }

146: /*@
147:     VecGhostRestoreLocalForm - Restores the local ghosted representation of
148:     a parallel vector obtained with VecGhostGetLocalForm().

150:     Logically Collective

152:     Input Parameter:
153: +   g - the global vector
154: -   l - the local (ghosted) representation

156:     Notes:
157:     This routine does not actually update the ghost values, but rather it
158:     returns a sequential vector that includes the locations for the ghost values
159:     and their current values.

161:     Level: advanced

163: .seealso: VecCreateGhost(), VecGhostGetLocalForm(), VecCreateGhostWithArray()
164: @*/
165: PetscErrorCode  VecGhostRestoreLocalForm(Vec g,Vec *l)
166: {

170:   if (*l) {
171:     VecGhostStateSync_Private(g,*l);
172:     PetscObjectDereference((PetscObject)*l);
173:   }
174:   return(0);
175: }

179: /*@
180:    VecGhostUpdateBegin - Begins the vector scatter to update the vector from
181:    local representation to global or global representation to local.

183:    Neighbor-wise Collective on Vec

185:    Input Parameters:
186: +  g - the vector (obtained with VecCreateGhost() or VecDuplicate())
187: .  insertmode - one of ADD_VALUES or INSERT_VALUES
188: -  scattermode - one of SCATTER_FORWARD or SCATTER_REVERSE

190:    Notes:
191:    Use the following to update the ghost regions with correct values from the owning process
192: .vb
193:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
194:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
195: .ve

197:    Use the following to accumulate the ghost region values onto the owning processors
198: .vb
199:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
200:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
201: .ve

203:    To accumulate the ghost region values onto the owning processors and then update
204:    the ghost regions correctly, call the later followed by the former, i.e.,
205: .vb
206:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
207:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
208:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
209:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
210: .ve

212:    Level: advanced

214: .seealso: VecCreateGhost(), VecGhostUpdateEnd(), VecGhostGetLocalForm(),
215:           VecGhostRestoreLocalForm(),VecCreateGhostWithArray()

217: @*/
218: PetscErrorCode  VecGhostUpdateBegin(Vec g,InsertMode insertmode,ScatterMode scattermode)
219: {
220:   Vec_MPI        *v;
222:   PetscBool      ismpi,isseq;

226:   PetscObjectTypeCompare((PetscObject)g,VECMPI,&ismpi);
227:   PetscObjectTypeCompare((PetscObject)g,VECSEQ,&isseq);
228:   if (ismpi) {
229:     v = (Vec_MPI*)g->data;
230:     if (!v->localrep) SETERRQ(PetscObjectComm((PetscObject)g),PETSC_ERR_ARG_WRONG,"Vector is not ghosted");
231:     if (!v->localupdate) return(0);
232:     if (scattermode == SCATTER_REVERSE) {
233:       VecScatterBegin(v->localupdate,v->localrep,g,insertmode,scattermode);
234:     } else {
235:       VecScatterBegin(v->localupdate,g,v->localrep,insertmode,scattermode);
236:     }
237:   } else if (isseq) {
238:     /* Do nothing */
239:   } else SETERRQ(PetscObjectComm((PetscObject)g),PETSC_ERR_ARG_WRONG,"Vector is not ghosted");
240:   return(0);
241: }

245: /*@
246:    VecGhostUpdateEnd - End the vector scatter to update the vector from
247:    local representation to global or global representation to local.

249:    Neighbor-wise Collective on Vec

251:    Input Parameters:
252: +  g - the vector (obtained with VecCreateGhost() or VecDuplicate())
253: .  insertmode - one of ADD_VALUES or INSERT_VALUES
254: -  scattermode - one of SCATTER_FORWARD or SCATTER_REVERSE

256:    Notes:

258:    Use the following to update the ghost regions with correct values from the owning process
259: .vb
260:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
261:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
262: .ve

264:    Use the following to accumulate the ghost region values onto the owning processors
265: .vb
266:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
267:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
268: .ve

270:    To accumulate the ghost region values onto the owning processors and then update
271:    the ghost regions correctly, call the later followed by the former, i.e.,
272: .vb
273:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
274:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
275:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
276:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
277: .ve

279:    Level: advanced

281: .seealso: VecCreateGhost(), VecGhostUpdateBegin(), VecGhostGetLocalForm(),
282:           VecGhostRestoreLocalForm(),VecCreateGhostWithArray()

284: @*/
285: PetscErrorCode  VecGhostUpdateEnd(Vec g,InsertMode insertmode,ScatterMode scattermode)
286: {
287:   Vec_MPI        *v;
289:   PetscBool      ismpi;

293:   PetscObjectTypeCompare((PetscObject)g,VECMPI,&ismpi);
294:   if (ismpi) {
295:     v = (Vec_MPI*)g->data;
296:     if (!v->localrep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Vector is not ghosted");
297:     if (!v->localupdate) return(0);
298:     if (scattermode == SCATTER_REVERSE) {
299:       VecScatterEnd(v->localupdate,v->localrep,g,insertmode,scattermode);
300:     } else {
301:       VecScatterEnd(v->localupdate,g,v->localrep,insertmode,scattermode);
302:     }
303:   }
304:   return(0);
305: }

309: PetscErrorCode VecSetOption_MPI(Vec v,VecOption op,PetscBool flag)
310: {
312:   if (op == VEC_IGNORE_OFF_PROC_ENTRIES)      v->stash.donotstash   = flag;
313:   else if (op == VEC_IGNORE_NEGATIVE_INDICES) v->stash.ignorenegidx = flag;
314:   return(0);
315: }


320: PetscErrorCode VecResetArray_MPI(Vec vin)
321: {
322:   Vec_MPI        *v = (Vec_MPI*)vin->data;

326:   v->array         = v->unplacedarray;
327:   v->unplacedarray = 0;
328:   if (v->localrep) {
329:     VecResetArray(v->localrep);
330:   }
331:   return(0);
332: }