Actual source code: commonmpvec.c

petsc-3.3-p7 2013-05-11
  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: {
 23:   PetscInt       gstate,lstate;

 26:   PetscObjectStateQuery((PetscObject)g,&gstate);
 27:   PetscObjectStateQuery((PetscObject)l,&lstate);
 28:   PetscObjectSetState((PetscObject)g,PetscMax(gstate,lstate));
 29:   PetscObjectSetState((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 PETSC_NULL if the Vec is not ghosted.

 41:     Not Collective

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

 46:     Output Parameter:
 47: .   l - the local (ghosted) representation, PETSC_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 = PETSC_NULL;
 94:   }
 95:   if (*l) {
 96:     VecGhostStateSync_Private(g,*l);
 97:     PetscObjectReference((PetscObject)*l);
 98:   }
 99:   return(0);
100: }

104: /*@
105:     VecGhostRestoreLocalForm - Restores the local ghosted representation of 
106:     a parallel vector obtained with VecGhostGetLocalForm().

108:     Not Collective

110:     Input Parameter:
111: +   g - the global vector
112: -   l - the local (ghosted) representation

114:     Notes:
115:     This routine does not actually update the ghost values, but rather it
116:     returns a sequential vector that includes the locations for the ghost values
117:     and their current values.

119:     Level: advanced

121: .seealso: VecCreateGhost(), VecGhostGetLocalForm(), VecCreateGhostWithArray()
122: @*/
123: PetscErrorCode  VecGhostRestoreLocalForm(Vec g,Vec *l)
124: {

128:   if (*l) {
129:     VecGhostStateSync_Private(g,*l);
130:     PetscObjectDereference((PetscObject)*l);
131:   }
132:   return(0);
133: }

137: /*@
138:    VecGhostUpdateBegin - Begins the vector scatter to update the vector from
139:    local representation to global or global representation to local.

141:    Neighbor-wise Collective on Vec

143:    Input Parameters:
144: +  g - the vector (obtained with VecCreateGhost() or VecDuplicate())
145: .  insertmode - one of ADD_VALUES or INSERT_VALUES
146: -  scattermode - one of SCATTER_FORWARD or SCATTER_REVERSE

148:    Notes:
149:    Use the following to update the ghost regions with correct values from the owning process
150: .vb
151:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
152:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
153: .ve

155:    Use the following to accumulate the ghost region values onto the owning processors
156: .vb
157:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
158:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
159: .ve

161:    To accumulate the ghost region values onto the owning processors and then update
162:    the ghost regions correctly, call the later followed by the former, i.e.,
163: .vb
164:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
165:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
166:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
167:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
168: .ve

170:    Level: advanced

172: .seealso: VecCreateGhost(), VecGhostUpdateEnd(), VecGhostGetLocalForm(),
173:           VecGhostRestoreLocalForm(),VecCreateGhostWithArray()

175: @*/
176: PetscErrorCode  VecGhostUpdateBegin(Vec g,InsertMode insertmode,ScatterMode scattermode)
177: {
178:   Vec_MPI        *v;
180:   PetscBool      ismpi,isseq;

184:   PetscObjectTypeCompare((PetscObject)g,VECMPI,&ismpi);
185:   PetscObjectTypeCompare((PetscObject)g,VECSEQ,&isseq);
186:   if (ismpi) {
187:     v  = (Vec_MPI*)g->data;
188:     if (!v->localrep) SETERRQ(((PetscObject)g)->comm,PETSC_ERR_ARG_WRONG,"Vector is not ghosted");
189:     if (!v->localupdate) return(0);
190:     if (scattermode == SCATTER_REVERSE) {
191:       VecScatterBegin(v->localupdate,v->localrep,g,insertmode,scattermode);
192:     } else {
193:       VecScatterBegin(v->localupdate,g,v->localrep,insertmode,scattermode);
194:     }
195:   } else if (isseq) {
196:     /* Do nothing */
197:   } else SETERRQ(((PetscObject)g)->comm,PETSC_ERR_ARG_WRONG,"Vector is not ghosted");
198:   return(0);
199: }

203: /*@
204:    VecGhostUpdateEnd - End the vector scatter to update the vector from
205:    local representation to global or global representation to local.

207:    Neighbor-wise Collective on Vec

209:    Input Parameters:
210: +  g - the vector (obtained with VecCreateGhost() or VecDuplicate())
211: .  insertmode - one of ADD_VALUES or INSERT_VALUES
212: -  scattermode - one of SCATTER_FORWARD or SCATTER_REVERSE

214:    Notes:

216:    Use the following to update the ghost regions with correct values from the owning process
217: .vb
218:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
219:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
220: .ve

222:    Use the following to accumulate the ghost region values onto the owning processors
223: .vb
224:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
225:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
226: .ve

228:    To accumulate the ghost region values onto the owning processors and then update
229:    the ghost regions correctly, call the later followed by the former, i.e.,
230: .vb
231:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
232:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
233:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
234:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
235: .ve

237:    Level: advanced

239: .seealso: VecCreateGhost(), VecGhostUpdateBegin(), VecGhostGetLocalForm(),
240:           VecGhostRestoreLocalForm(),VecCreateGhostWithArray()

242: @*/
243: PetscErrorCode  VecGhostUpdateEnd(Vec g,InsertMode insertmode,ScatterMode scattermode)
244: {
245:   Vec_MPI        *v;
247:   PetscBool      ismpi;

251:   PetscObjectTypeCompare((PetscObject)g,VECMPI,&ismpi);
252:   if (ismpi) {
253:     v  = (Vec_MPI*)g->data;
254:     if (!v->localrep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Vector is not ghosted");
255:     if (!v->localupdate) return(0);
256:     if (scattermode == SCATTER_REVERSE) {
257:       VecScatterEnd(v->localupdate,v->localrep,g,insertmode,scattermode);
258:     } else {
259:       VecScatterEnd(v->localupdate,g,v->localrep,insertmode,scattermode);
260:     }
261:   }
262:   return(0);
263: }

267: PetscErrorCode VecSetOption_MPI(Vec v,VecOption op,PetscBool  flag)
268: {
270:   if (op == VEC_IGNORE_OFF_PROC_ENTRIES) {
271:     v->stash.donotstash = flag;
272:   } else if (op == VEC_IGNORE_NEGATIVE_INDICES) {
273:     v->stash.ignorenegidx = flag;
274:   }
275:   return(0);
276: }


281: PetscErrorCode VecResetArray_MPI(Vec vin)
282: {
283:   Vec_MPI        *v = (Vec_MPI *)vin->data;

287:   v->array         = v->unplacedarray;
288:   v->unplacedarray = 0;
289:   if (v->localrep) {
290:     VecResetArray(v->localrep);
291:   }
292:   return(0);
293: }