Actual source code: commonmpvec.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2:  #include <../src/vec/vec/impls/mpi/pvecimpl.h>

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

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

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

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


 32: /*@
 33:     VecGhostGetLocalForm - Obtains the local ghosted representation of
 34:     a parallel vector (obtained with VecCreateGhost(), VecCreateGhostWithArray()
 35:     or VecCreateSeq()). Returns NULL if the Vec is not ghosted.

 37:     Logically Collective

 39:     Input Parameter:
 40: .   g - the global vector

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

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

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

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

 65:     Level: advanced

 67:    Concepts: vectors^ghost point access

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

 71: @*/
 72: PetscErrorCode  VecGhostGetLocalForm(Vec g,Vec *l)
 73: {
 75:   PetscBool      isseq,ismpi;


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

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

101:     Not Collective

103:     Input Parameter:
104: +   g - the global vector
105: -   l - the local vector

107:     Output Parameter:
108: .   flg - PETSC_TRUE if local vector is local form

110:     Level: advanced

112:    Concepts: vectors^ghost point access

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

116: @*/
117: PetscErrorCode VecGhostIsLocalForm(Vec g,Vec l,PetscBool *flg)
118: {
120:   PetscBool      isseq,ismpi;


126:   *flg = PETSC_FALSE;
127:   PetscObjectTypeCompare((PetscObject)g,VECSEQ,&isseq);
128:   PetscObjectTypeCompare((PetscObject)g,VECMPI,&ismpi);
129:   if (ismpi) {
130:     Vec_MPI *v = (Vec_MPI*)g->data;
131:     if (l == v->localrep) *flg = PETSC_TRUE;
132:   } else if (isseq) {
133:     if (l == g) *flg = PETSC_TRUE;
134:   } else SETERRQ(PetscObjectComm((PetscObject)g),PETSC_ERR_ARG_WRONG,"Global vector is not ghosted");
135:   return(0);
136: }

138: /*@
139:     VecGhostRestoreLocalForm - Restores the local ghosted representation of
140:     a parallel vector obtained with VecGhostGetLocalForm().

142:     Logically Collective

144:     Input Parameter:
145: +   g - the global vector
146: -   l - the local (ghosted) representation

148:     Notes:
149:     This routine does not actually update the ghost values, but rather it
150:     returns a sequential vector that includes the locations for the ghost values
151:     and their current values.

153:     Level: advanced

155: .seealso: VecCreateGhost(), VecGhostGetLocalForm(), VecCreateGhostWithArray()
156: @*/
157: PetscErrorCode  VecGhostRestoreLocalForm(Vec g,Vec *l)
158: {

162:   if (*l) {
163:     VecGhostStateSync_Private(g,*l);
164:     PetscObjectDereference((PetscObject)*l);
165:   }
166:   return(0);
167: }

169: /*@
170:    VecGhostUpdateBegin - Begins the vector scatter to update the vector from
171:    local representation to global or global representation to local.

173:    Neighbor-wise Collective on Vec

175:    Input Parameters:
176: +  g - the vector (obtained with VecCreateGhost() or VecDuplicate())
177: .  insertmode - one of ADD_VALUES or INSERT_VALUES
178: -  scattermode - one of SCATTER_FORWARD or SCATTER_REVERSE

180:    Notes:
181:    Use the following to update the ghost regions with correct values from the owning process
182: .vb
183:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
184:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
185: .ve

187:    Use the following to accumulate the ghost region values onto the owning processors
188: .vb
189:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
190:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
191: .ve

193:    To accumulate the ghost region values onto the owning processors and then update
194:    the ghost regions correctly, call the later followed by the former, i.e.,
195: .vb
196:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
197:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
198:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
199:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
200: .ve

202:    Level: advanced

204: .seealso: VecCreateGhost(), VecGhostUpdateEnd(), VecGhostGetLocalForm(),
205:           VecGhostRestoreLocalForm(),VecCreateGhostWithArray()

207: @*/
208: PetscErrorCode  VecGhostUpdateBegin(Vec g,InsertMode insertmode,ScatterMode scattermode)
209: {
210:   Vec_MPI        *v;
212:   PetscBool      ismpi,isseq;

216:   PetscObjectTypeCompare((PetscObject)g,VECMPI,&ismpi);
217:   PetscObjectTypeCompare((PetscObject)g,VECSEQ,&isseq);
218:   if (ismpi) {
219:     v = (Vec_MPI*)g->data;
220:     if (!v->localrep) SETERRQ(PetscObjectComm((PetscObject)g),PETSC_ERR_ARG_WRONG,"Vector is not ghosted");
221:     if (!v->localupdate) return(0);
222:     if (scattermode == SCATTER_REVERSE) {
223:       VecScatterBegin(v->localupdate,v->localrep,g,insertmode,scattermode);
224:     } else {
225:       VecScatterBegin(v->localupdate,g,v->localrep,insertmode,scattermode);
226:     }
227:   } else if (isseq) {
228:     /* Do nothing */
229:   } else SETERRQ(PetscObjectComm((PetscObject)g),PETSC_ERR_ARG_WRONG,"Vector is not ghosted");
230:   return(0);
231: }

233: /*@
234:    VecGhostUpdateEnd - End the vector scatter to update the vector from
235:    local representation to global or global representation to local.

237:    Neighbor-wise Collective on Vec

239:    Input Parameters:
240: +  g - the vector (obtained with VecCreateGhost() or VecDuplicate())
241: .  insertmode - one of ADD_VALUES or INSERT_VALUES
242: -  scattermode - one of SCATTER_FORWARD or SCATTER_REVERSE

244:    Notes:

246:    Use the following to update the ghost regions with correct values from the owning process
247: .vb
248:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
249:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
250: .ve

252:    Use the following to accumulate the ghost region values onto the owning processors
253: .vb
254:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
255:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
256: .ve

258:    To accumulate the ghost region values onto the owning processors and then update
259:    the ghost regions correctly, call the later followed by the former, i.e.,
260: .vb
261:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
262:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
263:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
264:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
265: .ve

267:    Level: advanced

269: .seealso: VecCreateGhost(), VecGhostUpdateBegin(), VecGhostGetLocalForm(),
270:           VecGhostRestoreLocalForm(),VecCreateGhostWithArray()

272: @*/
273: PetscErrorCode  VecGhostUpdateEnd(Vec g,InsertMode insertmode,ScatterMode scattermode)
274: {
275:   Vec_MPI        *v;
277:   PetscBool      ismpi;

281:   PetscObjectTypeCompare((PetscObject)g,VECMPI,&ismpi);
282:   if (ismpi) {
283:     v = (Vec_MPI*)g->data;
284:     if (!v->localrep) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Vector is not ghosted");
285:     if (!v->localupdate) return(0);
286:     if (scattermode == SCATTER_REVERSE) {
287:       VecScatterEnd(v->localupdate,v->localrep,g,insertmode,scattermode);
288:     } else {
289:       VecScatterEnd(v->localupdate,g,v->localrep,insertmode,scattermode);
290:     }
291:   }
292:   return(0);
293: }