1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecGhostGetLocalForm.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>VecGhostGetLocalForm</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.7.5 2017-01-01</b></div>
<div id="bugreport" align=right><a href="mailto:petsc-maint@mcs.anl.gov?subject=Typo or Error in Documentation &body=Please describe the typo or error in the documentation: petsc-3.7.5 v3.7.5 docs/manualpages/Vec/VecGhostGetLocalForm.html "><small>Report Typos and Errors</small></a></div>
<A NAME="VecGhostGetLocalForm"><H1>VecGhostGetLocalForm</H1></A>
Obtains the local ghosted representation of a parallel vector (obtained with <A HREF="../Vec/VecCreateGhost.html#VecCreateGhost">VecCreateGhost</A>(), <A HREF="../Vec/VecCreateGhostWithArray.html#VecCreateGhostWithArray">VecCreateGhostWithArray</A>() or <A HREF="../Vec/VecCreateSeq.html#VecCreateSeq">VecCreateSeq</A>()). Returns NULL if the <A HREF="../Vec/Vec.html#Vec">Vec</A> is not ghosted.
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscvec.h"
PetscErrorCode VecGhostGetLocalForm(Vec g,Vec *l)
</PRE>
Logically Collective
<P>
<H3><FONT COLOR="#CC3333">Input Parameter</FONT></H3>
<DT><B>g </B> -the global vector
<br>
<P>
<H3><FONT COLOR="#CC3333">Output Parameter</FONT></H3>
<DT><B>l </B> -the local (ghosted) representation, NULL if g is not ghosted
<br>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
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.
<P>
To update the ghost values from the locations on the other processes one must call
<A HREF="../Vec/VecGhostUpdateBegin.html#VecGhostUpdateBegin">VecGhostUpdateBegin</A>() and <A HREF="../Vec/VecGhostUpdateEnd.html#VecGhostUpdateEnd">VecGhostUpdateEnd</A>() before accessing the ghost values. Thus normal
usage is
<pre>
<A HREF="../Vec/VecGhostUpdateBegin.html#VecGhostUpdateBegin">VecGhostUpdateBegin</A>(x,<A HREF="../Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</A>,<A HREF="../Vec/SCATTER_FORWARD.html#SCATTER_FORWARD">SCATTER_FORWARD</A>);
</pre>
<pre>
<A HREF="../Vec/VecGhostUpdateEnd.html#VecGhostUpdateEnd">VecGhostUpdateEnd</A>(x,<A HREF="../Sys/INSERT_VALUES.html#INSERT_VALUES">INSERT_VALUES</A>,<A HREF="../Vec/SCATTER_FORWARD.html#SCATTER_FORWARD">SCATTER_FORWARD</A>);
</pre>
<pre>
<A HREF="../Vec/VecGhostGetLocalForm.html#VecGhostGetLocalForm">VecGhostGetLocalForm</A>(x,&xlocal);
</pre>
<pre>
<A HREF="../Vec/VecGetArray.html#VecGetArray">VecGetArray</A>(xlocal,&xvalues);
</pre>
<pre>
// access the non-ghost values in locations xvalues[0:n-1] and ghost values in locations xvalues[n:n+nghost];
</pre>
<pre>
<A HREF="../Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</A>(xlocal,&xvalues);
</pre>
<pre>
<A HREF="../Vec/VecGhostRestoreLocalForm.html#VecGhostRestoreLocalForm">VecGhostRestoreLocalForm</A>(x,&xlocal);
</pre>
<P>
One should call <A HREF="../Vec/VecGhostRestoreLocalForm.html#VecGhostRestoreLocalForm">VecGhostRestoreLocalForm</A>() or <A HREF="../Vec/VecDestroy.html#VecDestroy">VecDestroy</A>() once one is
finished using the object.
<P>
<P>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../Vec/VecCreateGhost.html#VecCreateGhost">VecCreateGhost</A>(), <A HREF="../Vec/VecGhostRestoreLocalForm.html#VecGhostRestoreLocalForm">VecGhostRestoreLocalForm</A>(), <A HREF="../Vec/VecCreateGhostWithArray.html#VecCreateGhostWithArray">VecCreateGhostWithArray</A>()
<BR>
<P>
<P><B><P><B><FONT COLOR="#CC3333">Level:</FONT></B>advanced
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/vec/vec/impls/mpi/commonmpvec.c.html#VecGhostGetLocalForm">src/vec/vec/impls/mpi/commonmpvec.c</A>
<BR><A HREF="./index.html">Index of all Vec routines</A>
<BR><A HREF="../../index.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
<P><H3><FONT COLOR="#CC3333">Examples</FONT></H3>
<A HREF="../../../src/vec/vec/examples/tutorials/ex9.c.html">src/vec/vec/examples/tutorials/ex9.c.html</A><BR>
<A HREF="../../../src/vec/vec/examples/tutorials/ex9f.F.html">src/vec/vec/examples/tutorials/ex9f.F.html</A><BR>
<A HREF="../../../src/vec/vec/examples/tutorials/ex14f.F.html">src/vec/vec/examples/tutorials/ex14f.F.html</A><BR>
</BODY></HTML>
|