Dear all,
I am trying to use VectorTools::compute_nonzero_normal_flux_constraints
to apply a normal displacement field to part of the boundary. The code
works if I simply take a ConstantFunction to impose as normal flux.
I have tried defining a more general function as a list of points and
values associated with those points (using vector_value_list). This
definition is based on the output of another simulation (with
post-processing of the solution gradient at the boundary). I used the
vertices located along the boundary as list of points. I don't know if it
makes sense.
I have attached the relevant parts of code below. I have two main questions:
- Are there specific points where the function should be defined? (For
instance the face Gauss points along the boundary)
- I get an error of dimensions mismatch (also copied below), could you help
me understand the reason? I am solving the problem in 2D
Thank you for this great resource. I have been using deal ii for more than
3 years.
Giovanna
- - - - - - - -
std::set<types::boundary_id> normal_flux_boundaries;
normal_flux_boundaries.insert (Xnormal_positive_plane);
typename FunctionMap<dim>::type boundary_functions;
BoundaryValues<dim> pointwise_velocity;
boundary_functions[0] = &pointwise_velocity;
VectorTools::compute_nonzero_normal_flux_constraints ( dof_handler,
0, // first_vector_component,
normal_flux_boundaries,
boundary_functions,
constraints,
StaticMappingQ1< dim >::mapping
);
// Below is the definition of the BoundaryValue function
template <int dim>
class BoundaryValues: public Function<dim>
{
public:
BoundaryValues () {};
virtual void vector_value_list (const std::vector<Point<dim> > &points,
std::vector<Vector<double> > &values)
const;
};
template <int dim>
void BoundaryValues<dim>::vector_value_list(const std::vector<Point<dim> >
&points,
std::vector<Vector<double> > &values) const
{
Assert(values.size()==points.size(),
ExcDimensionMismatch(values.size(),points.size()));
// read the list of points and values from a file
std::ifstream input1("Values_Interface.csv");
unsigned int n_Points;
input1 >> n_Points;
AssertThrow (input1, ExcIO());
std::vector<Point<dim> > listNodes(n_Points);
std::vector<Vector<double> > listValues(n_Points);
for (unsigned int node=0; node<n_Points; ++node)
{
double x[dim];
char delim;
double readValue;
AssertThrow (input1, ExcIO());
input1 >> x[0] >> delim >> x[1] >> delim >> readValue;
for (unsigned int d=0; d<dim; ++d)
{
listNodes[node](d) = x[d];
listValues[node](d) = readValue;
}
}
const double tol = 1.e-4 *
GridTools::minimal_cell_diameter(triangulation);
for (unsigned int i=0; i<points.size(); ++i)
{
for (unsigned int j=0; j<n_Points; ++j)
{
if (points[i].distance(listNodes[j]) < tol)
values[i] = listValues[j];
}
}
}
- - - - - - - -
An error occurred in line <76> of file
</Users/gio/development/dealii-8.5.1/include/deal.II/base/function.templates.h>
in function
virtual void dealii::Function<2, double>::vector_value(const Point<dim>
&, Vector<Number> &) const [dim = 2, Number = double]
The violated condition was:
(v.size()) == (this->n_components)
Additional information:
Dimension 2 not equal to 1.
--
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.