Neumann bc in step-20

76 views
Skip to first unread message

Erik Svensson

unread,
Sep 29, 2016, 1:32:04 AM9/29/16
to deal.II User Group
Hi,

I’m trying to implement homogeneous Neumann bc on part of the boundary in tutorial example step-20.

The Neumann bc should be imposed strongly (as for Dirichlet bc in normal, non-mixed, formulation). Is this correct?

In any case, I’m trying to impose the Neumann bc strongly using the technique from step-22.

Running program I get the error:

---
An error occurred in line <1878> of file </usr/local/dealii/dealii-8.4.1/include/deal.II/numerics/vector_tools.templates.h> in function
    void dealii::VectorTools::{anonymous}::do_interpolate_boundary_values(const M_or_MC<DoFHandlerType:: dimension, DoFHandlerType:: space_dimension>&, const DoFHandlerType&, const typename dealii::FunctionMap<DoFHandlerType:: space_dimension>::type&, std::map<unsigned int, double>&, const dealii::ComponentMask&, dealii::internal::int2type<dim_>) [with DoFHandlerType = dealii::DoFHandler<2>; M_or_MC = dealii::Mapping; int dim_ = 2; typename dealii::FunctionMap<DoFHandlerType:: space_dimension>::type = std::map<unsigned char, const dealii::Function<2, double>*, std::less<unsigned char>, std::allocator<std::pair<const unsigned char, const dealii::Function<2, double>*> > >]
The violated condition was:
    cell->get_fe().is_primitive (i)
The name and call sequence of the exception was:
    ExcMessage ("This function can only deal with requested boundary " "values that correspond to primitive (scalar) base " "elements")
Additional Information:
This function can only deal with requested boundary values that correspond to primitive (scalar) base elements
---

Can someone please help me to interpret the error message, or even better guide me how to proceed implementing the Neumann bc?

Best regards, Erik

Daniel Arndt

unread,
Sep 29, 2016, 5:58:07 AM9/29/16
to deal.II User Group
Erik,


Am Donnerstag, 29. September 2016 07:32:04 UTC+2 schrieb Erik Svensson:
Hi,

I’m trying to implement homogeneous Neumann bc on part of the boundary in tutorial example step-20.

The Neumann bc should be imposed strongly (as for Dirichlet bc in normal, non-mixed, formulation). Is this correct?

In any case, I’m trying to impose the Neumann bc strongly using the technique from step-22.
What exactly are you doing? As far as I see step-22 just uses weak Neumann bc and if theses are homogeneous you normally don't have to implement anything additional.


Running program I get the error:

---
An error occurred in line <1878> of file </usr/local/dealii/dealii-8.4.1/include/deal.II/numerics/vector_tools.templates.h> in function
    void dealii::VectorTools::{anonymous}::do_interpolate_boundary_values(const M_or_MC<DoFHandlerType:: dimension, DoFHandlerType:: space_dimension>&, const DoFHandlerType&, const typename dealii::FunctionMap<DoFHandlerType:: space_dimension>::type&, std::map<unsigned int, double>&, const dealii::ComponentMask&, dealii::internal::int2type<dim_>) [with DoFHandlerType = dealii::DoFHandler<2>; M_or_MC = dealii::Mapping; int dim_ = 2; typename dealii::FunctionMap<DoFHandlerType:: space_dimension>::type = std::map<unsigned char, const dealii::Function<2, double>*, std::less<unsigned char>, std::allocator<std::pair<const unsigned char, const dealii::Function<2, double>*> > >]
The violated condition was:
    cell->get_fe().is_primitive (i)
The name and call sequence of the exception was:
    ExcMessage ("This function can only deal with requested boundary " "values that correspond to primitive (scalar) base " "elements")
Additional Information:
This function can only deal with requested boundary values that correspond to primitive (scalar) base elements
---
For non-primitive FiniteElements (e.g. FE_RaviartThomas) you can't use VectorTools::interpolate_boundary_values. Use VectorTools::project_boundary_values [1] instead.

Best,
Daniel

[1] https://www.dealii.org/8.4.1/doxygen/deal.II/group__constraints.html#ga8981786b2369d1c6db8f5c363e290edc

Wolfgang Bangerth

unread,
Sep 29, 2016, 11:54:11 PM9/29/16
to dea...@googlegroups.com
Erik -- imposing zero Neumann boundary values for the mixed Laplace implies
that the normal velocity is zero on the boundary. This is no problem if you
discretize the velocity variable through the ordinary Lagrange elements, but
in step-20 we use the div-conforming Raviart-Thomas element for which this
function is not equipped to handle the situation because the Raviart-Thomas
element mixes both the x- and y-components of the velocity vector.

I think you should be able to use
VectorTools::project_boundary_values_div_conforming
instead to do what you want to do.

Best
W.


--
------------------------------------------------------------------------
Wolfgang Bangerth email: bang...@colostate.edu
www: http://www.math.colostate.edu/~bangerth/

Erik Svensson

unread,
Sep 30, 2016, 1:39:04 AM9/30/16
to dea...@googlegroups.com
Ok, thanks Wolfgang and Daniel!

I will follow your advise and post updates on this issue.

/Erik



--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en
--- You received this message because you are subscribed to a topic in the Google Groups "deal.II User Group" group.
To unsubscribe from this topic, visit https://groups.google.com/d/topic/dealii/Tvs-V4ESiZE/unsubscribe.
To unsubscribe from this group and all its topics, send an email to dealii+unsubscribe@googlegroups.com.
For more options, visit https://groups.google.com/d/optout.

Erik Svensson

unread,
Oct 5, 2016, 2:14:24 AM10/5/16
to dea...@googlegroups.com
The 'VectorTools::project_boundary_values_div_conforming' worked very well. Thanks a lot Wolfgang!
/Erik

Reply all
Reply to author
Forward
0 new messages