[petsc-users] Handling non-solution local vector data

Matthew Knepley knepley at gmail.com
Sat Dec 7 13:34:20 CST 2013


On Sat, Dec 7, 2013 at 1:26 PM, Mark Lohry <mlohry at gmail.com> wrote:

> 1) Is there somewhere on the doc site with a single list of all examples
> and what they do? As a new user I hugely appreciate all the great examples
> and documentation, but it's a little unwieldy to try and find relevant
> examples by navigating the function documents.
>

No, we have not yet setup a nice system for this.


> 2) I'm a little unclear as to the best-practice way to handle
> solution/grid related data that isn't actual solution DOFS or mesh points,
> and doesn't require global communication.
>
> Say I want to work on an (ni X nj ) cell-centered finite volume layout
> with one DOF per cell. I start with
>
> DMDACreate2d(PETSC_COMM_WORLD,
>       DMDA_BOUNDARY_PERIODIC, DMDA_BOUNDARY_NONE, /* periodic in i, no
> conn in j */
>       DMDA_STENCIL_BOX, /* 5-point stencil */
>       ni,nj, /* number of cells */
>       PETSC_DECIDE,PETSC_DECIDE, /* let petsc determine number of procs */
>       1, /* one degree of freedom per cell */
>       1, /* stencil width */
>       NULL,NULL,&da);
>
> Now I can construct the solution vector like so:
>
> struct Field { PetscScalar u; };
> Field **soln;
> Vec global_vec;
> DMDACreateGlobalVector(da,&global_vec);
> DMDAVecGetArray(da,global_vec,&soln); // converts global vector into
> Field**
>
> and access the local elements through soln[j][i]. But say I want to store
> some grid information locally, say the FaceNormal for a cell face:
> struct FaceNormal{ PetscScalar nx,ny; };
> FaceNormal **IFace;
>
> Now it seems that the CreateGlobalVector/CreateLocalVector are just going
> to allocate for 1 DOF because of how the DMDA was defined, so there
> wouldn't be sufficient storage for the 2 variables in FaceNormal. What's an
> appropriate way to do this? Clearly I could just malloc that myself, but
> it'd be nice to keep the same indexing as the petsc vectors.
>

There are at least two answers:

  1) If you do not care about redundant storage, you could create another
2d DMDA using dof = 8, and store all the normals for
       each cell. DMDAs take almost no memory. Your storage would be the
Vec you use to store the normals.

  2) If you care about storage, then create another DMDA, but then create a
PetscSection to assign 2 dof to every face, which you can do with

         DMDACreateSection()

      and use DMSetDefaultSection(). This will make DMCreateGlobalVector()
create a Vec you can use to store all the face normals.

     Matt


> Thanks in advance,
>



-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20131207/2ce90124/attachment.html>


More information about the petsc-users mailing list