xref: /petsc/src/vec/vec/impls/mpi/commonmpvec.c (revision 0e30de3c61834e263b2aa8210e130ed660b4065b)
1 #include <../src/vec/vec/impls/mpi/pvecimpl.h> /*I  "petscvec.h"   I*/
2 
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.
9 
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.
13 
14   Note that this function is correct and changes nothing if both arguments are the
15   same, which is the case in serial.
16 */
VecGhostStateSync_Private(Vec g,Vec l)17 static PetscErrorCode VecGhostStateSync_Private(Vec g, Vec l)
18 {
19   PetscObjectState gstate, lstate;
20 
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 }
28 
29 /*@
30   VecGhostGetLocalForm - Obtains the local ghosted representation of
31   a parallel vector (obtained with `VecCreateGhost()`, `VecCreateGhostWithArray()` or `VecCreateSeq()`).
32 
33   Logically Collective
34 
35   Input Parameter:
36 . g - the global vector
37 
38   Output Parameter:
39 . l - the local (ghosted) representation,`NULL` if `g` is not ghosted
40 
41   Level: advanced
42 
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.
48 
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      VecGetArrayRead(xlocal, &xvalues);
57         // access the non-ghost values in locations xvalues[0:n-1] and ghost values in locations xvalues[n:n+nghost];
58      VecRestoreArrayRead(xlocal, &xvalues);
59      VecGhostRestoreLocalForm(x, &xlocal);
60 .ve
61 
62   One must call `VecGhostRestoreLocalForm()` once finished using the object.
63 
64 .seealso: [](ch_vectors), `VecGhostUpdateBegin()`, `VecGhostUpdateEnd()`, `Vec`, `VecType`, `VecCreateGhost()`, `VecGhostRestoreLocalForm()`, `VecCreateGhostWithArray()`
65 @*/
VecGhostGetLocalForm(Vec g,Vec * l)66 PetscErrorCode VecGhostGetLocalForm(Vec g, Vec *l)
67 {
68   PetscBool isseq, ismpi;
69 
70   PetscFunctionBegin;
71   PetscValidHeaderSpecific(g, VEC_CLASSID, 1);
72   PetscAssertPointer(l, 2);
73 
74   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECSEQ, &isseq));
75   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECMPI, &ismpi));
76   if (ismpi) {
77     Vec_MPI *v = (Vec_MPI *)g->data;
78     *l         = v->localrep;
79   } else if (isseq) {
80     *l = g;
81   } else {
82     *l = NULL;
83   }
84   if (*l) PetscCall(VecGhostStateSync_Private(g, *l));
85   PetscFunctionReturn(PETSC_SUCCESS);
86 }
87 
88 /*@
89   VecGhostIsLocalForm - Checks if a given vector is the local form of a global vector
90 
91   Not Collective
92 
93   Input Parameters:
94 + g - the global vector
95 - l - the local vector
96 
97   Output Parameter:
98 . flg - `PETSC_TRUE` if `l` is the local form
99 
100   Level: advanced
101 
102 .seealso: [](ch_vectors), `Vec`, `VecType`, `VecCreateGhost()`, `VecGhostRestoreLocalForm()`, `VecCreateGhostWithArray()`, `VecGhostGetLocalForm()`
103 @*/
VecGhostIsLocalForm(Vec g,Vec l,PetscBool * flg)104 PetscErrorCode VecGhostIsLocalForm(Vec g, Vec l, PetscBool *flg)
105 {
106   PetscBool isseq, ismpi;
107 
108   PetscFunctionBegin;
109   PetscValidHeaderSpecific(g, VEC_CLASSID, 1);
110   PetscValidHeaderSpecific(l, VEC_CLASSID, 2);
111 
112   *flg = PETSC_FALSE;
113   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECSEQ, &isseq));
114   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECMPI, &ismpi));
115   if (ismpi) {
116     Vec_MPI *v = (Vec_MPI *)g->data;
117     if (l == v->localrep) *flg = PETSC_TRUE;
118   } else if (isseq) {
119     if (l == g) *flg = PETSC_TRUE;
120   } else SETERRQ(PetscObjectComm((PetscObject)g), PETSC_ERR_ARG_WRONG, "Global vector is not ghosted");
121   PetscFunctionReturn(PETSC_SUCCESS);
122 }
123 
124 /*@
125   VecGhostRestoreLocalForm - Restores the local ghosted representation of
126   a parallel vector obtained with `VecGhostGetLocalForm()`.
127 
128   Logically Collective
129 
130   Input Parameters:
131 + g - the global vector
132 - l - the local (ghosted) representation
133 
134   Level: advanced
135 
136 .seealso: [](ch_vectors), `VecGhostUpdateBegin()`, `VecGhostUpdateEnd()`, `Vec`, `VecType`, `VecCreateGhost()`, `VecGhostGetLocalForm()`, `VecCreateGhostWithArray()`
137 @*/
VecGhostRestoreLocalForm(Vec g,Vec * l)138 PetscErrorCode VecGhostRestoreLocalForm(Vec g, Vec *l)
139 {
140   PetscFunctionBegin;
141   if (*l) PetscCall(VecGhostStateSync_Private(g, *l));
142   PetscFunctionReturn(PETSC_SUCCESS);
143 }
144 
145 /*@
146   VecGhostUpdateBegin - Begins the vector scatter to update the vector from
147   local representation to global or global representation to local.
148 
149   Neighbor-wise Collective
150 
151   Input Parameters:
152 + g           - the vector (obtained with `VecCreateGhost()` or `VecDuplicate()`)
153 . insertmode  - one of `ADD_VALUES`, `MAX_VALUES`, `MIN_VALUES` or `INSERT_VALUES`
154 - scattermode - one of `SCATTER_FORWARD` (update ghosts) or `SCATTER_REVERSE` (update local values from ghosts)
155 
156   Level: advanced
157 
158   Notes:
159   Use the following to update the ghost regions with correct values from the owning process
160 .vb
161        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
162        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
163 .ve
164 
165   Use the following to accumulate the ghost region values onto the owning processors
166 .vb
167        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
168        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
169 .ve
170 
171   To accumulate the ghost region values onto the owning processors and then update
172   the ghost regions correctly, call the latter followed by the former, i.e.,
173 .vb
174        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
175        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
176        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
177        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
178 .ve
179 
180 .seealso: [](ch_vectors), `Vec`, `VecType`, `VecCreateGhost()`, `VecGhostUpdateEnd()`, `VecGhostGetLocalForm()`,
181           `VecGhostRestoreLocalForm()`, `VecCreateGhostWithArray()`
182 @*/
VecGhostUpdateBegin(Vec g,InsertMode insertmode,ScatterMode scattermode)183 PetscErrorCode VecGhostUpdateBegin(Vec g, InsertMode insertmode, ScatterMode scattermode)
184 {
185   Vec_MPI  *v;
186   PetscBool ismpi, isseq;
187 
188   PetscFunctionBegin;
189   PetscValidHeaderSpecific(g, VEC_CLASSID, 1);
190   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECMPI, &ismpi));
191   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECSEQ, &isseq));
192   if (ismpi) {
193     v = (Vec_MPI *)g->data;
194     PetscCheck(v->localrep, PetscObjectComm((PetscObject)g), PETSC_ERR_ARG_WRONG, "Vector is not ghosted");
195     if (!v->localupdate) PetscFunctionReturn(PETSC_SUCCESS);
196     if (scattermode == SCATTER_REVERSE) {
197       PetscCall(VecScatterBegin(v->localupdate, v->localrep, g, insertmode, scattermode));
198     } else {
199       PetscCall(VecScatterBegin(v->localupdate, g, v->localrep, insertmode, scattermode));
200     }
201   } else if (isseq) {
202     /* Do nothing */
203   } else SETERRQ(PetscObjectComm((PetscObject)g), PETSC_ERR_ARG_WRONG, "Vector is not ghosted");
204   PetscFunctionReturn(PETSC_SUCCESS);
205 }
206 
207 /*@
208   VecGhostUpdateEnd - End the vector scatter to update the vector from
209   local representation to global or global representation to local.
210 
211   Neighbor-wise Collective
212 
213   Input Parameters:
214 + g           - the vector (obtained with `VecCreateGhost()` or `VecDuplicate()`)
215 . insertmode  - one of `ADD_VALUES`, `MAX_VALUES`, `MIN_VALUES` or `INSERT_VALUES`
216 - scattermode - one of `SCATTER_FORWARD` (update ghosts) or `SCATTER_REVERSE` (update local values from ghosts)
217 
218   Level: advanced
219 
220   Notes:
221   Use the following to update the ghost regions with correct values from the owning process
222 .vb
223        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
224        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
225 .ve
226 
227   Use the following to accumulate the ghost region values onto the owning processors
228 .vb
229        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
230        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
231 .ve
232 
233   To accumulate the ghost region values onto the owning processors and then update
234   the ghost regions correctly, call the later followed by the former, i.e.,
235 .vb
236        VecGhostUpdateBegin(v,ADD_VALUES,SCATTER_REVERSE);
237        VecGhostUpdateEnd(v,ADD_VALUES,SCATTER_REVERSE);
238        VecGhostUpdateBegin(v,INSERT_VALUES,SCATTER_FORWARD);
239        VecGhostUpdateEnd(v,INSERT_VALUES,SCATTER_FORWARD);
240 .ve
241 
242 .seealso: [](ch_vectors), `Vec`, `VecType`, `VecCreateGhost()`, `VecGhostUpdateBegin()`, `VecGhostGetLocalForm()`,
243           `VecGhostRestoreLocalForm()`, `VecCreateGhostWithArray()`
244 @*/
VecGhostUpdateEnd(Vec g,InsertMode insertmode,ScatterMode scattermode)245 PetscErrorCode VecGhostUpdateEnd(Vec g, InsertMode insertmode, ScatterMode scattermode)
246 {
247   Vec_MPI  *v;
248   PetscBool ismpi;
249 
250   PetscFunctionBegin;
251   PetscValidHeaderSpecific(g, VEC_CLASSID, 1);
252   PetscCall(PetscObjectTypeCompare((PetscObject)g, VECMPI, &ismpi));
253   if (ismpi) {
254     v = (Vec_MPI *)g->data;
255     PetscCheck(v->localrep, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Vector is not ghosted");
256     if (!v->localupdate) PetscFunctionReturn(PETSC_SUCCESS);
257     if (scattermode == SCATTER_REVERSE) {
258       PetscCall(VecScatterEnd(v->localupdate, v->localrep, g, insertmode, scattermode));
259     } else {
260       PetscCall(VecScatterEnd(v->localupdate, g, v->localrep, insertmode, scattermode));
261     }
262   }
263   PetscFunctionReturn(PETSC_SUCCESS);
264 }
265