Hi.

There are several point_value functions. Right now, you are calling one that 
only works for scalar elements:

https://www.dealii.org/current/doxygen/deal.II/namespaceVectorTools.html#a7be5c7eed52308898dfaad91c4cff204

You need to call the more general function:

https://www.dealii.org/current/doxygen/deal.II/namespaceVectorTools.html#acd358e9b110ccbf4a7f76796d206b9c7

by passing a vector with each components value as the last argument to the 
function:

Vector<double> value(n_solution_components);
VectorTools::point_value(dof_handler,
                         solution,
                         point,
                         value);

Best,
Simon

On 20/07/2023 12:54, 'Jost Arndt' via deal.II User Group wrote:
Dear everyone,

I have trouble understanding in how to evaluate a solution at a given point 
correctly.
Described in Step-3 for example the function point_value is used to evaluate 
the solution at (0.3, 0.3), and I want to do something quite similar.
However, my solution is vector valued (i.e. I use FESystem ), and I get the 
following error message

 >>An error occurred in line <181> of file 
<./include/deal.II/numerics/vector_tools_point_value.templates.h> in function
     typename VectorType::value_type dealii::VectorTools::point_value(const dealii::Mapping<dim, 
spacedim>&, const dealii::DoFHandler<dim, spacedim>&, const VectorType&, const 
dealii::Point<spacedim>&) [with int dim =
2; VectorType = dealii::Vector<double>; int spacedim = 2; typename 
VectorType::value_type = double]
The violated condition was:
     dof.get_fe(0).n_components() == 1
Additional information:
     Finite element is not scalar as is necessary for this function


However in the documentation it says "Evaluate a possibly vector-valued finite element 
function defined by the given DoFHandler 
<https://dealii.org/developer/doxygen/deal.II/classDoFHandler.html> and nodal vector 
fe_function at the given point point". So I am curious how to evaluate my solution at a 
few fixed points? Did I miss something?

Performance-wise it does not have to be perfect at all, as I want to evaluate 
only ~400 points on a few hundered solutions (time-steps).

Best,

Jost

--
The deal.II project is located at http://www.dealii.org/ 
<http://www.dealii.org/>
For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en 
<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] 
<mailto:[email protected]>.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/7e2a3e13-aef1-4a0e-a357-8dafc3bdf958n%40googlegroups.com
 
<https://groups.google.com/d/msgid/dealii/7e2a3e13-aef1-4a0e-a357-8dafc3bdf958n%40googlegroups.com?utm_medium=email&utm_source=footer>.

--
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].
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/229072fe-6f09-1470-35b6-ee34ddc308fd%40gmail.com.

Reply via email to