[Dune] Position of dofs for finite elements

Timo Koch timo.koch at iws.uni-stuttgart.de
Fri Jun 19 13:22:12 CEST 2020


-- 
_________________________________________________

Timo Koch                                      phone: +49 711 685 64676
IWS, Universität Stuttgart                  fax:   +49 711 685 60430
Pfaffenwaldring 61         email: timo.koch at iws.uni-stuttgart.de
D-70569 Stuttgart             url: www.iws.uni-stuttgart.de/en/lh2/
_________________________________________________

> On 19. Jun 2020, at 12:38, Dedner, Andreas <A.S.Dedner at warwick.ac.uk> wrote:
> 
> Hi.
> I might have misunderstood you but the way I see it the term `dof position` doesn't actually make any sense in the concept of `dune-localfunctions`. The interface there is too generic . As an example I could define a continuous Lagrange space using moments on the edges instead of point values. The only thing available then is the information to which subentity a dof is connected. The issue that the local keys are often 'the wrong way around' has been discussed a few times... But your approach of using the local interpolation to get some point coordinates would of course fail with that type of space.
> Since dune-localfunctions is generic we didn't add a concept of `Lagrangepoints`.

Hi Andreas,

thanks a lot for the explanation, that already makes sense a lot more sense now. That means that our interface `problem.dirichlet(element, localPos)` is not very useful, if we want it to be more general.
We’ll try to revisit the approach.

Timo


> 
> In this generic setting the only way to get dof values is to use the local interpolation as you do but for the Dirichlet function. I.e. you interpolate your localized dirichlet data on to the element (we don't have a method to only perform the interpolation to a subentity). That gives you a local dof vector and then you extract the dofs you want to constrain using the local dof key.
> 
> As I said, I'm not sure if this helps because as I said I might have misunderstood your problem (not being too familiar with dune-functions doesn't help of course...)
> 
> Best
> Andreas
> From: Dune <dune-bounces at lists.dune-project.org> on behalf of Timo Koch <timo.koch at iws.uni-stuttgart.de>
> Sent: 19 June 2020 11:13
> To: dune at lists.dune-project.org <dune at lists.dune-project.org>
> Subject: [Dune] Position of dofs for finite elements
>  
> Dear Dune users & developers,
> 
> we are using dune-functions and the core modules in a simple FEM code. We want to strongly enforce Dirichlet boundary conditions.
> The boundary conditions are specified locally for each element as a function of the element and the local position.
> 
> For globally defined analytical functions, we could use a global interpolate to get the coefficient / Dirichlet values. But how do we do this for “piece-wise”-defined boundary conditions?
> We came up with something like this (restriction to simple Lagrange basis, nothing composite and assuming all dofs are Dirichlet for now):
> 
> Dune::BlockVector<Dune::FieldVector<double, numEq>> dirichletValues;
> auto getDirichletValues = [&] (auto localIndex, const auto& localView, const auto& intersection)
> {
>     const auto& element = localView.element();
> 
>     // get the local coordinates of the dofs of the element
>     std::vector<double> coords;
>     const auto& interp = localView.tree().finiteElement().localInterpolation();
>     interp.interpolate([&] (const LocalCoordinate& x) { return x[0]; }, coords);
>     localDofCoords.resize(coords.size());
>     for (unsigned int i = 0; i < coords.size(); ++i)
>         localDofCoords[i][0] = coords[i];
> 
>     interp.interpolate([&] (const LocalCoordinate& x) { return x[1]; }, coords);
>     for (unsigned int i = 0; i < coords.size(); ++i)
>         localDofCoords[i][1] = coords[i];
> 
>     const auto index = localView.index(localIndex);
>     for (unsigned int eqIdx = 0; eqIdx < numEq; ++eqIdx)
>         dirichletValues[index] = problem.dirichlet(element, localDofCoords[localIndex]);
> };
> 
> Dune::Functions::forEachBoundaryDOF(lagrangeBasis, getDirichletValues);
> 
> which is far from nice. In this way the coordinates are recomputed every time although it would be enough to compute them once per element type. Also we could figure out how to get the full position instead of only the x,y,z-components separately.
> There seems to be no obvious way to get the local position of a DOF within an element, although this sounds like a very basic thing to do.
> With the `localKey` it is possible to obtain the subEntity and its center position via the reference element. But this doesn’t (necessarily) correspond to the DOF position of course.
> 
> Is there some obvious way that we are missing on how to solve this problem without writing an excessive amount of code?
> Do you have any suggestion or comments on how “piece-wise” defined strongly enforced Dirichlet boundary conditions could be realized?
> 
> Thanks a lot for your help
> Timo
> 
> -- 
> _________________________________________________
> 
> Timo Koch                                      phone: +49 711 685 64676
> IWS, Universität Stuttgart                  fax:   +49 711 685 60430
> Pfaffenwaldring 61         email: timo.koch at iws.uni-stuttgart.de <mailto:timo.koch at iws.uni-stuttgart.de>
> D-70569 Stuttgart             url: www.iws.uni-stuttgart.de/en/lh2/ <http://www.iws.uni-stuttgart.de/en/lh2/>
> _________________________________________________

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <https://lists.dune-project.org/pipermail/dune/attachments/20200619/c1369625/attachment.htm>


More information about the Dune mailing list