petsc-3.14.6 2021-03-30
VecGhostGetLocalForm
Obtains the local ghosted representation of a parallel vector (obtained with VecCreateGhost(), VecCreateGhostWithArray() or VecCreateSeq()). Returns NULL if the Vec is not ghosted.
Synopsis
#include "petscvec.h"
PetscErrorCode VecGhostGetLocalForm(Vec g,Vec *l)
Logically Collective
Input Parameter
Output Parameter
| l | - the local (ghosted) representation, NULL if g is not ghosted
|
Notes
This routine does not actually update the ghost values, but rather it
returns a sequential vector that includes the locations for the ghost
values and their current values. The returned vector and the original
vector passed in share the same array that contains the actual vector data.
To update the ghost values from the locations on the other processes one must call
VecGhostUpdateBegin() and VecGhostUpdateEnd() before accessing the ghost values. Thus normal
usage is
VecGhostUpdateBegin(x,INSERT_VALUES,SCATTER_FORWARD);
VecGhostUpdateEnd(x,INSERT_VALUES,SCATTER_FORWARD);
VecGhostGetLocalForm(x,&xlocal);
VecGetArray(xlocal,&xvalues);
// access the non-ghost values in locations xvalues[0:n-1] and ghost values in locations xvalues[n:n+nghost];
VecRestoreArray(xlocal,&xvalues);
VecGhostRestoreLocalForm(x,&xlocal);
One should call VecGhostRestoreLocalForm() or VecDestroy() once one is
finished using the object.
See Also
VecCreateGhost(), VecGhostRestoreLocalForm(), VecCreateGhostWithArray()
Level
advanced
Location
src/vec/vec/impls/mpi/commonmpvec.c
Examples
src/vec/vec/tutorials/ex9.c.html
src/vec/vec/tutorials/ex9f.F90.html
src/vec/vec/tutorials/ex14f.F90.html
Index of all Vec routines
Table of Contents for all manual pages
Index of all manual pages