[petsc-users] Testing DMPlex new functionalities

Matthew Knepley knepley at gmail.com
Thu Jan 22 16:12:02 CST 2015


On Thu, Jan 22, 2015 at 12:08 PM, Luc Berger-Vergiat <lb2653 at columbia.edu>
wrote:

> Hi all,
> I'm trying to test some functionalities in DMPlex, specifically I want to
> do a 2d fem simulation using a mesh generated by GMSH.
> At this point I am trying to understand how Dirichlet BC can be applied.
>
> Can you explain to me what arguments I have to pass to:
> DMPlexAddBoundary()?
>

  ierr = DMPlexAddBoundary(dm, PETSC_TRUE,"wall","marker",(void
(*)())top_displacement,1,&id,user)

PETSC_TRUE means that these unknowns will be eliminated from the global
system

"wall" is the name of the BC

"marker" is the name of the DMLabel that you want to use to identify the
mesh points

1 is the number of label IDs that identify constrained mesh points

id is the label ID which identifies constrained mesh points

top_displacement is the function used to set boundary values


> I am attaching a little bit of c code and my mesh.geo and mesh.msh so that
> you understand the context.
> The BCs I want to apply are displacement of 0.2 in the y direction if
> y==1.0
> 0.0 in the x direction if y==0.0


I might do this:

void top_displacement(const PetscReal x[], PetscScalar *u, void *ctx)
{
  *u = x[1] < PETSC_MACHINE_EPSILON ? 0.0 : 0.2;
}

  Thanks,

     Matt


>
> --
> Best,
> Luc
>
>


-- 
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/20150122/c5d3db08/attachment.html>


More information about the petsc-users mailing list