Dear Simon,
Thank you for your reply. Your suggestion definitely works. I used
functions that output one component when I was interpolating scalar
potentials. It works well.
However, the FE_Nedelec is different. It implements edge elements. They are
vector-based. That is, functions are represented by a superposition of
vector-valued shape functions:
\vec{A} = \sum u_i vec{N}_i .
Therefore, the output of the "value" method in “class CustomFunction :
public Function” must be vector-valued. Three components in a
three-dimensional space. Otherwise, there is no point in interpolation.
This kind of vectorial approximations is a bread-and-butter topic in
magnetics. See, for example, equation (7) in:
https://ieeexplore.ieee.org/document/497322
In short, the source, i.e the current vector potential, must be projected
on the space spanned by the vector-valued shape functions. Otherwise, the
simulation is numerically unstable. The last, from my experience, is
definitely true.
Best,
John
On Fri, Apr 16, 2021 at 5:28 PM [email protected] <[email protected]>
wrote:
> Hi,
> I don't have any experience with FE_Nedelec, but the error message states
> that the function you are trying to interpolate has 1 component, while the
> element has dim components.
>
> If you implemented your own Function, you need to make sure it has
> dim-components by calling the constructor of Function with dim as
> inargument:
>
> template <int dim>
> class CustomFunction : public Function<dim>
> {
> public:
> CustomFunction()
> : Function<dim>(dim) // The function has dim components
> {}
>
> double
> value(const Point<dim> & point,
> const unsigned int component = 0) const override
> {
> // Return value will depend on component here
> return 0;
> }
> };
>
>
> Best,
> Simon
>
> On Friday, April 16, 2021 at 5:07:01 PM UTC+2 [email protected] wrote:
>
>> Hello,
>>
>> It seems I am unable to find a function similar to
>> VectorTools::interpolate for FE_Nedelec finite elements. The existing
>> implementation of this functions gives the following run-time error:
>>
>> *An error occurred in line <556> of file
>> </home/john/dealii-9.1.1/include/deal.II*
>>
>> */numerics/vector_tools.templates.h> in function*
>>
>> * void dealii::VectorTools::interpolate(const dealii::Mapping<dim,
>> spacedim>&,*
>>
>> * const DoFHandlerType<dim, spacedim>&, const dealii::Function<spacedim,
>> typename*
>>
>> * VectorType::value_type>&, VectorType&, const dealii::ComponentMask&)
>> [with int *
>>
>> *dim = 3; int spacedim = 3; VectorType = dealii::Vector<double>;
>> DoFHandlerType =*
>>
>> * dealii::DoFHandler; typename VectorType::value_type = double]*
>>
>> *The violated condition was: *
>>
>> * dof_handler.get_fe().n_components() == function.n_components*
>>
>> *Additional information: *
>>
>> * Dimension 3 not equal to 1.*
>>
>>
>> Looks like VectorTools::interpolate does not like the vector-valued
>> finite elements. Could you please point me in the right direction?
>>
>>
>> Best,
>>
>> John
>>
> --
> 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/wG3j4yD1rpw/unsubscribe.
> To unsubscribe from this group and all its topics, send an email to
> [email protected].
> To view this discussion on the web visit
> https://groups.google.com/d/msgid/dealii/6d051c9e-8840-4e5d-8021-f8007a1dc543n%40googlegroups.com
> <https://groups.google.com/d/msgid/dealii/6d051c9e-8840-4e5d-8021-f8007a1dc543n%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/CAH3L9AVtN%3Dtaqi2YiTjfoyfhxo-%2B9XnNzpBqNXSHqk5wX4iUtQ%40mail.gmail.com.