← Back to team overview

dolfin team mailing list archive

Re: displacement condition on individual nodes

 

Quoting Ola Skavhaug <skavhaug@xxxxxxxxx>:

> On Mon, Aug 10, 2009 at 5:46 AM, Catherine Micek<mice0012@xxxxxxx> wrote:
> >
> > On Aug 9, 2009, at 3:24 AM, Kristian Oelgaard wrote:
> >
> >> Quoting Peter Brune <prbrune@xxxxxxxxx>:
> >>
> >>> As a quick solution, you could, knowing the coordinates of a node,
> >>> constrain
> >>> everything within some epsilon smaller than the mesh feature size
> >>> of the
> >>> node (such as DOLFIN_EPS).
> >>
> >> Yes, and then make sure to use the "pointwise" method when creating
> >> the
> >> DirichletBC. See comment in dolfin/fem/DirichletBC.h.
> >> You can also look at demo/pde/dg/advection-diffusion/main.cpp for
> >> an example on
> >> how to use the "geometric" approach. The majority of DOLFIN demos
> >> uses the
> >> "topological" approach which is the default.
> >
> > This seems simple enough, but I am having a little difficulty with
> > it.  I am coding with the python interface, so I looked at the python
> > version of the demo in demo/pde/dg/advection-diffusion/.  In there,
> > the Dirichlet boundary is specified as the line x=1 on the boundary,
> > i.e.
> >
> > class DirichletBoundary(SubDomain):
> >     def inside(self, x, on_boundary):
> >        return (abs(x[0] - 1.0) < DOLFIN_EPS and on_boundary
> >
> > and I'm unclear as to why it's necessary to specify "pointwise" as
> > opposed to "topological" boundary conditions.  So this condition
> > would fix all the nodes on the line x=1; it doesn't appear to be
> > picking out points.  Is it necessary here because we're using DG
> > elements?

Yes, you can use "geometric" for DG elements where you would normally use
"topological" for other elements. In this case it's not necessary to use
"pointwise" because we fix all nodes on the line x=1. Note that "topological"
is faster than "geometric" which is faster than "pointwise", how big the
difference is in practice I don't know.

> > The other question is related to syntax.  Say I want to fix the point
> > (1,0).  When I set the boundary as
> >
> > class DirichletBoundary(SubDomain):
> >     def inside(self, x, on_boundary):
> >        return ( (abs(x[0] - 1.0) and abs(x[1] - 0.0) < DOLFIN_EPS and
> > on_boundary)
>
> This is not correct, since abs(x[0]-1.0]) almost always evaluates to True.
> Try
> return ( (abs(x[0] - 1.0) < DOLFIN_EPS and abs(x[1] - 0.0) <
> DOLFIN_EPS and on_boundary)
>
> Ola
>
> > and set the Dirichlet boundary condition with the option set to
> > "geometric," I am not getting the correct solution.  Is this not the
> > correct way to set the boundary conditions?

No, you want to set the dofs at a point, so you should use "pointwise". Also
note that the variable 'on_boundary' is ALWAYS False when using the "pointwise"
method, so you should modify Ola's suggestion to:

class DirichletBoundary(SubDomain):
  def inside(self, x, on_boundary):
    return abs(x[0] - 1.0) < DOLFIN_EPS and abs(x[1] - 0.0) < DOLFIN_EPS

Kristian

> > Again, any suggestions would be appreciated -- thanks!
> >
> > Katy
> >
> >
> >>
> >> Kristian
> >>
> >>> - Peter
> >>>
> >>> On Sat, Aug 8, 2009 at 11:27 PM, Catherine Micek
> >>> <mice0012@xxxxxxx> wrote:
> >>>
> >>>> Hi,
> >>>> I'm working with a group that's using Fenics to study the
> >>>> equations of
> >>>> linear elasticity (with a pure displacement formulation), and I
> >>>> have a
> >>>> question about boundary conditions.  Our goal is to study the
> >>>> pure traction
> >>>> problem, and this means we have to be careful in formulating the
> >>>> boundary
> >>>> conditions.  In order to prevent a singularity in our equations,
> >>>> we need to
> >>>> constrain certain nodes so as to prevent translations and rigid
> >>>> rotations.
> >>>>  There are many examples in the fenics demos about how to
> >>>> prescribe mixed
> >>>> conditions, but I have only seen examples where the Dirichlet
> >>>> condition
> >>>> occurs along an edge (as opposed to occurring at a few individual
> >>>> nodes).
> >>>>  So my question is, given a mesh, how does one prescribe Dirichlet
> >>>> conditions for individual nodes?  Are there any demos that
> >>>> address this?
> >>>>
> >>>> Any help would be greatly appreciated -- thanks!
> >>>>
> >>>> Katy
> >>>>
> >>>> ---
> >>>> Catherine (Katy) A. Micek
> >>>> Graduate Assistant, U of MN Mathematics Department
> >>>> http://www.math.umn.edu/~mice0012 <http://www.math.umn.edu/%
> >>>> 7Emice0012>
> >>>>
> >>>>
> >>>>
> >>>>
> >>>>
> >>>> _______________________________________________
> >>>> DOLFIN-dev mailing list
> >>>> DOLFIN-dev@xxxxxxxxxx
> >>>> http://www.fenics.org/mailman/listinfo/dolfin-dev
> >>>>
> >>>>
> >>>
> >>
> >>
> >
> > _______________________________________________
> > DOLFIN-dev mailing list
> > DOLFIN-dev@xxxxxxxxxx
> > http://www.fenics.org/mailman/listinfo/dolfin-dev
> >
>
>
>
> --
> Ola Skavhaug
>




Follow ups

References