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 78 79 80 81 82 83 84
|
<!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/VecGetArray.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>VecGetArray</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.4.2 2013-07-02</b></div>
<A NAME="VecGetArray"><H1>VecGetArray</H1></A>
Returns a pointer to a contiguous array that contains this processor's portion of the vector data. For the standard PETSc vectors, <A HREF="../Vec/VecGetArray.html#VecGetArray">VecGetArray</A>() returns a pointer to the local data array and does not use any copies. If the underlying vector data is not stored in a contiquous array this routine will copy the data to a contiquous array and return a pointer to that. You MUST call <A HREF="../Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</A>() when you no longer need access to the array.
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscvec.h"
PetscErrorCode VecGetArray(Vec x,PetscScalar **a)
</PRE>
Logically Collective on <A HREF="../Vec/Vec.html#Vec">Vec</A>
<P>
<H3><FONT COLOR="#CC3333">Input Parameter</FONT></H3>
<DT><B>x </B> -the vector
<br>
<P>
<H3><FONT COLOR="#CC3333">Output Parameter</FONT></H3>
<DT><B>a </B> -location to put pointer to the array
<br>
<P>
<H3><FONT COLOR="#CC3333">Fortran Note</FONT></H3>
This routine is used differently from Fortran 77
<pre>
<A HREF="../Vec/Vec.html#Vec">Vec</A> x
</pre>
<pre>
<A HREF="../Sys/PetscScalar.html#PetscScalar">PetscScalar</A> x_array(1)
</pre>
<pre>
<A HREF="../Sys/PetscOffset.html#PetscOffset">PetscOffset</A> i_x
</pre>
<pre>
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> ierr
</pre>
<pre>
call <A HREF="../Vec/VecGetArray.html#VecGetArray">VecGetArray</A>(x,x_array,i_x,ierr)
</pre>
<pre>
</pre>
<pre>
Access first local entry in vector with
</pre>
<pre>
value = x_array(i_x + 1)
</pre>
<pre>
</pre>
<pre>
...... other code
</pre>
<pre>
call <A HREF="../Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</A>(x,x_array,i_x,ierr)
</pre>
For Fortran 90 see <A HREF="../Vec/VecGetArrayF90.html#VecGetArrayF90">VecGetArrayF90</A>()
<P>
See the Fortran chapter of the users manual and
petsc/src/snes/examples/tutorials/ex5f.F for details.
<P>
<P>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../Vec/VecRestoreArray.html#VecRestoreArray">VecRestoreArray</A>(), <A HREF="../Vec/VecGetArrayRead.html#VecGetArrayRead">VecGetArrayRead</A>(), <A HREF="../Vec/VecGetArrays.html#VecGetArrays">VecGetArrays</A>(), <A HREF="../Vec/VecGetArrayF90.html#VecGetArrayF90">VecGetArrayF90</A>(), <A HREF="../Vec/VecPlaceArray.html#VecPlaceArray">VecPlaceArray</A>(), <A HREF="../Vec/VecGetArray2d.html#VecGetArray2d">VecGetArray2d</A>()
<BR><P><B><P><B><FONT COLOR="#CC3333">Level:</FONT></B>beginner
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/vec/vec/interface/rvector.c.html#VecGetArray">src/vec/vec/interface/rvector.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>
<A HREF="../../../src/vec/vec/examples/tutorials/ex3.c.html">src/vec/vec/examples/tutorials/ex3.c.html</A><BR>
<A HREF="../../../src/vec/vec/examples/tutorials/ex6.c.html">src/vec/vec/examples/tutorials/ex6.c.html</A><BR>
<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/ex18.c.html">src/vec/vec/examples/tutorials/ex18.c.html</A><BR>
<A HREF="../../../src/vec/vec/examples/tutorials/ex19.c.html">src/vec/vec/examples/tutorials/ex19.c.html</A><BR>
<A HREF="../../../src/vec/vec/examples/tutorials/ex4f.F.html">src/vec/vec/examples/tutorials/ex4f.F.html</A><BR>
<A HREF="../../../src/mat/examples/tutorials/ex4.c.html">src/mat/examples/tutorials/ex4.c.html</A><BR>
<A HREF="../../../src/mat/examples/tutorials/ex12.c.html">src/mat/examples/tutorials/ex12.c.html</A><BR>
<A HREF="../../../src/dm/examples/tutorials/ex2.c.html">src/dm/examples/tutorials/ex2.c.html</A><BR>
<A HREF="../../../src/dm/examples/tutorials/ex5.c.html">src/dm/examples/tutorials/ex5.c.html</A><BR>
</BODY></HTML>
|