On 09/28/2016 11:32 PM, Erik Svensson wrote:
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?
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: [email protected]
www: http://www.math.colostate.edu/~bangerth/
--
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 the Google Groups "deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email
to [email protected].
For more options, visit https://groups.google.com/d/optout.