[petsc-dev] DMDA_BOUNDARY_GHOSTED

Matthew Knepley knepley at gmail.com
Thu Sep 26 13:13:03 CDT 2013


On Thu, Sep 26, 2013 at 11:10 AM, Mark F. Adams <mfadams at lbl.gov> wrote:

> What Chombo usually does, which is not great but it simplifies things
> greatly, is fill ghost cells with, for instance, the negative of the
> boundary cell for low order homogenous Dirichlet BCs, then run your kernel
> (once), rinse and repeat.  I can do this in my FormFunction routine.  I was
> thinking this would be problem for G-S smoothers.  Chombo calls a setBC
> method in our G-S code to do this but I'm guessing PETSc does not.  Now I'm
> thinking I can just use damped Jacob, using my operator, and that should
> take care of all places where the stencil is used.  This is constant
> coefficient Laplacian and am just trying to get a mock up code to do SR.
>

We have ghost cells in TS ex11, and I think we do them the correct way :)
We have a function you register that fills them up, and we
loop over them explicitly.


> Also, I know this is a dumb question, and I don't really need to know
> this, but I can't understand this, and am very curious how code like this
> can work:
>
>  PetscScalar    **array;
>   ierr = DMDAVecGetArray(da, b, &array);CHKERRQ(ierr);
>   for (j=ys; j<ys+ym; j++) {
>     for (i=xs; i<xs+xm; i++) {
>       array[j][i] = ….
>     }
>   }
>
> This is fantastic.  You seem to be doing arbitrary(ish) indexing into a
> multidimensional, dynamically sized, array.  Just like fortran.  I never
> knew you could do this and can not see how the compiler even implements
> this as you are not even telling it any dimension as far as I can see.  I
> obviously missing something very simple here.
>

This is one of the most successful Barryisms. We construct the
multi-dimensional array by-hand, putting in the correct pointers
to the relevant parts of the 1D array in the Vec.

   Matt


> Mark
>
> On Sep 26, 2013, at 12:36 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
>
> >
> >  For vertex centered finite differences it is easy, one simply does not
> include the boundary nodes in the DMDA and then puts the actual Dirichlet
> values into the ghost slots.
> >
> >  For cell centered since the boundary value is the average of the ghost
> cell and the first cell I don't know how to do this because I don't know
> what one should put into the ghost location. Maybe include the "ghost" cell
> as a real point (don't use boundary_ghosted) and make the equation for that
> last cell force the average of that slot and its neighbor to equal the
> Dirichlet value.
> >
> >
> >   Barry
> >
> > On Sep 26, 2013, at 11:28 AM, Matthew Knepley <knepley at gmail.com> wrote:
> >
> >> On Thu, Sep 26, 2013 at 9:25 AM, Mark F. Adams <mfadams at lbl.gov> wrote:
> >> I am interested in making a cell centered test problem with Dirichlet
> BCs.  I want to go high(ish) order and doing explicit stencils (preferred)
> looks rather messy and I am interest in using ghost cells.  I assume this
> would entail using DMDA_BOUNDARY_GHOSTED.  I do not see any examples or
> documentation on using this other than just creating it.  Any suggestions?
>  I would think that there would have to be callback function for the user
> to set ghost values but I am not seeing anything like this around.  Perhaps
> I am not understand the intent of DMDA_BOUNDARY_GHOSTED.
> >>
> >> My understanding is that it adds a single layer of ghost cells around
> the domain boundary. Then these show up in the gxm, etc. so that
> >> you can put values there if you want.
> >>
> >>  Matt
> >>
> >>
> >> Mark
> >>
> >>
> >>
> >> --
> >> 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
> >
>
>


-- 
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-dev/attachments/20130926/f4253ad3/attachment.html>


More information about the petsc-dev mailing list