Actual source code: commonmpvec.c

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

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

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

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

 21:   PetscFunctionBegin;
 22:   PetscCall(PetscObjectStateGet((PetscObject)g, &gstate));
 23:   PetscCall(PetscObjectStateGet((PetscObject)l, &lstate));
 24:   PetscCall(PetscObjectStateSet((PetscObject)g, PetscMax(gstate, lstate)));
 25:   PetscCall(PetscObjectStateSet((PetscObject)l, PetscMax(gstate, lstate)));
 26:   PetscFunctionReturn(PETSC_SUCCESS);
 27: }

 29: /*@
 30:   VecGhostGetLocalForm - Obtains the local ghosted representation of
 31:   a parallel vector (obtained with `VecCreateGhost()`, `VecCreateGhostWithArray()` or `VecCreateSeq()`).

 33:   Logically Collective

 35:   Input Parameter:
 36: . g - the global vector

 38:   Output Parameter:
 39: . l - the local (ghosted) representation,`NULL` if `g` is not ghosted

 41:   Level: advanced

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

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

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

 65: .seealso: [](ch_vectors), `Vec`, `VecType`, `VecCreateGhost()`, `VecGhostRestoreLocalForm()`, `VecCreateGhostWithArray()`
 66: @*/
 67: PetscErrorCode VecGhostGetLocalForm(Vec g, Vec *l)
 68: {
 69:   PetscBool isseq, ismpi;

 71:   PetscFunctionBegin;
 73:   PetscAssertPointer(l, 2);

 75:   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECSEQ, &isseq));
 76:   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECMPI, &ismpi));
 77:   if (ismpi) {
 78:     Vec_MPI *v = (Vec_MPI *)g->data;
 79:     *l         = v->localrep;
 80:   } else if (isseq) {
 81:     *l = g;
 82:   } else {
 83:     *l = NULL;
 84:   }
 85:   if (*l) {
 86:     PetscCall(VecGhostStateSync_Private(g, *l));
 87:     PetscCall(PetscObjectReference((PetscObject)*l));
 88:   }
 89:   PetscFunctionReturn(PETSC_SUCCESS);
 90: }

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

 95:   Not Collective

 97:   Input Parameters:
 98: + g - the global vector
 99: - l - the local vector

101:   Output Parameter:
102: . flg - `PETSC_TRUE` if `l` is the local form

104:   Level: advanced

106: .seealso: [](ch_vectors), `Vec`, `VecType`, `VecCreateGhost()`, `VecGhostRestoreLocalForm()`, `VecCreateGhostWithArray()`, `VecGhostGetLocalForm()`
107: @*/
108: PetscErrorCode VecGhostIsLocalForm(Vec g, Vec l, PetscBool *flg)
109: {
110:   PetscBool isseq, ismpi;

112:   PetscFunctionBegin;

116:   *flg = PETSC_FALSE;
117:   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECSEQ, &isseq));
118:   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECMPI, &ismpi));
119:   if (ismpi) {
120:     Vec_MPI *v = (Vec_MPI *)g->data;
121:     if (l == v->localrep) *flg = PETSC_TRUE;
122:   } else if (isseq) {
123:     if (l == g) *flg = PETSC_TRUE;
124:   } else SETERRQ(PetscObjectComm((PetscObject)g), PETSC_ERR_ARG_WRONG, "Global vector is not ghosted");
125:   PetscFunctionReturn(PETSC_SUCCESS);
126: }

128: /*@
129:   VecGhostRestoreLocalForm - Restores the local ghosted representation of
130:   a parallel vector obtained with `VecGhostGetLocalForm()`.

132:   Logically Collective

134:   Input Parameters:
135: + g - the global vector
136: - l - the local (ghosted) representation

138:   Level: advanced

140:   Note:
141:   This routine does not actually update the ghost values, but rather it
142:   returns a sequential vector that includes the locations for the ghost values
143:   and their current values.

145: .seealso: [](ch_vectors), `Vec`, `VecType`, `VecCreateGhost()`, `VecGhostGetLocalForm()`, `VecCreateGhostWithArray()`
146: @*/
147: PetscErrorCode VecGhostRestoreLocalForm(Vec g, Vec *l)
148: {
149:   PetscFunctionBegin;
150:   if (*l) {
151:     PetscCall(VecGhostStateSync_Private(g, *l));
152:     PetscCall(PetscObjectDereference((PetscObject)*l));
153:   }
154:   PetscFunctionReturn(PETSC_SUCCESS);
155: }

157: /*@
158:   VecGhostUpdateBegin - Begins the vector scatter to update the vector from
159:   local representation to global or global representation to local.

161:   Neighbor-wise Collective

163:   Input Parameters:
164: + g           - the vector (obtained with `VecCreateGhost()` or `VecDuplicate()`)
165: . insertmode  - one of `ADD_VALUES`, `MAX_VALUES`, `MIN_VALUES` or `INSERT_VALUES`
166: - scattermode - one of `SCATTER_FORWARD` or `SCATTER_REVERSE`

168:   Level: advanced

170:   Notes:
171:   Use the following to update the ghost regions with correct values from the owning process
172: .vb
173:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
174:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
175: .ve

177:   Use the following to accumulate the ghost region values onto the owning processors
178: .vb
179:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
180:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
181: .ve

183:   To accumulate the ghost region values onto the owning processors and then update
184:   the ghost regions correctly, call the latter followed by the former, i.e.,
185: .vb
186:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
187:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
188:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
189:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
190: .ve

192: .seealso: [](ch_vectors), `Vec`, `VecType`, `VecCreateGhost()`, `VecGhostUpdateEnd()`, `VecGhostGetLocalForm()`,
193:           `VecGhostRestoreLocalForm()`, `VecCreateGhostWithArray()`
194: @*/
195: PetscErrorCode VecGhostUpdateBegin(Vec g, InsertMode insertmode, ScatterMode scattermode)
196: {
197:   Vec_MPI  *v;
198:   PetscBool ismpi, isseq;

200:   PetscFunctionBegin;
202:   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECMPI, &ismpi));
203:   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECSEQ, &isseq));
204:   if (ismpi) {
205:     v = (Vec_MPI *)g->data;
206:     PetscCheck(v->localrep, PetscObjectComm((PetscObject)g), PETSC_ERR_ARG_WRONG, "Vector is not ghosted");
207:     if (!v->localupdate) PetscFunctionReturn(PETSC_SUCCESS);
208:     if (scattermode == SCATTER_REVERSE) {
209:       PetscCall(VecScatterBegin(v->localupdate, v->localrep, g, insertmode, scattermode));
210:     } else {
211:       PetscCall(VecScatterBegin(v->localupdate, g, v->localrep, insertmode, scattermode));
212:     }
213:   } else if (isseq) {
214:     /* Do nothing */
215:   } else SETERRQ(PetscObjectComm((PetscObject)g), PETSC_ERR_ARG_WRONG, "Vector is not ghosted");
216:   PetscFunctionReturn(PETSC_SUCCESS);
217: }

219: /*@
220:   VecGhostUpdateEnd - End the vector scatter to update the vector from
221:   local representation to global or global representation to local.

223:   Neighbor-wise Collective

225:   Input Parameters:
226: + g           - the vector (obtained with `VecCreateGhost()` or `VecDuplicate()`)
227: . insertmode  - one of `ADD_VALUES`, `MAX_VALUES`, `MIN_VALUES` or `INSERT_VALUES`
228: - scattermode - one of `SCATTER_FORWARD` or `SCATTER_REVERSE`

230:   Level: advanced

232:   Notes:
233:   Use the following to update the ghost regions with correct values from the owning process
234: .vb
235:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
236:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
237: .ve

239:   Use the following to accumulate the ghost region values onto the owning processors
240: .vb
241:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
242:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
243: .ve

245:   To accumulate the ghost region values onto the owning processors and then update
246:   the ghost regions correctly, call the later followed by the former, i.e.,
247: .vb
248:        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
249:        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
250:        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
251:        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
252: .ve

254: .seealso: [](ch_vectors), `Vec`, `VecType`, `VecCreateGhost()`, `VecGhostUpdateBegin()`, `VecGhostGetLocalForm()`,
255:           `VecGhostRestoreLocalForm()`, `VecCreateGhostWithArray()`
256: @*/
257: PetscErrorCode VecGhostUpdateEnd(Vec g, InsertMode insertmode, ScatterMode scattermode)
258: {
259:   Vec_MPI  *v;
260:   PetscBool ismpi;

262:   PetscFunctionBegin;
264:   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECMPI, &ismpi));
265:   if (ismpi) {
266:     v = (Vec_MPI *)g->data;
267:     PetscCheck(v->localrep, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Vector is not ghosted");
268:     if (!v->localupdate) PetscFunctionReturn(PETSC_SUCCESS);
269:     if (scattermode == SCATTER_REVERSE) {
270:       PetscCall(VecScatterEnd(v->localupdate, v->localrep, g, insertmode, scattermode));
271:     } else {
272:       PetscCall(VecScatterEnd(v->localupdate, g, v->localrep, insertmode, scattermode));
273:     }
274:   }
275:   PetscFunctionReturn(PETSC_SUCCESS);
276: }