[petsc-users] DMPLex ghost indexing

Matthew Knepley knepley at gmail.com
Wed Jul 22 17:10:40 CDT 2015


On Wed, Jul 22, 2015 at 3:53 PM, Adrian Croucher <a.croucher at auckland.ac.nz>
wrote:

> hi
>
> I have a distributed DMPlex representing a finite volume mesh, and first I
> want to check if I understand the way the indexing of ghost cells works,
> based on my experiments.
>
> As I understand it, there are two types of ghost cells:
>
> 1) those created by DMPlexDistribute() (which I call with overlap = 1),
> representing the cells over the edge of the processor partition boundary
> 2) those created by DMPlexConstructGhostCells(), which I call to add ghost
> cells on the physical boundary of the mesh (using a DMLabel), for applying
> boundary conditions
>
> The type 1 'partition' ghost cells only appear in local vectors, not
> global, and are added on the end after the interior cells.
>

The are only in local vectors, but they can be anywhere in the list of
cells. They are marked in the ghost label. It might be that they are now at
the end.


> The type 2 'boundary' ghost cells appear in both local and global vectors.
> In a local vector they are added on the end after the partition ghost
> cells. In a global vector they appear straight after the interior cells.
>

Yes, they are required to be after all interior/partition cells.


> Is this correct? If so, I have a question:
>
> If I'm computing fluxes between my finite volume cells, I need a local
> vector to get access to both types of ghost cells.
>

Yes.


> But if I'm doing another calculation (computing fluid properties in the
> cells, in this case) that does not need access to the partition ghosts, but
> needs to be carried out for all interior cells *and* the boundary ghost
> cells, I figure it would be preferable to do this with a global vector (to
> save needless parallel communication scattering a global vector into a
> local one).
>

I would use whichever one makes more send in the flow of the algorithm.


> However if I use DMPlexGetHeightStratum() to get the start and end cells,
> the indices it returns apply to a local vector with the partition ghosts
> included.
>

Yes.


> Is it safe just to subtract off the number of partition ghost cells from
> the end cell index returned by DMPlexGetHeightStratum(), if I want to loop
> over all interior and boundary ghosts cells in a global vector?
>

Nope. If you really do not want partition ghost cells, you would need

    ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
    ierr = DMLabelGetValue(ghostLabel, cell, &ghost);CHKERRQ(ierr);
    if (ghost >= 0) continue;

I only use this to avoid ghost faces right now.

  Thanks,

     Matt


> Cheers, Adrian
>
> --
> Dr Adrian Croucher
> Department of Engineering Science
> University of Auckland
> New Zealand
> tel 64-9-373-7599 ext 84611
>
>


-- 
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/20150722/251f1c26/attachment.html>


More information about the petsc-users mailing list