Using DoFTools::map_dofs_to_support_points for two-field with FE_DGPMonomial

33 views
Skip to first unread message

Masoud Ahmadi

unread,
Feb 7, 2023, 1:41:17 PM2/7/23
to deal.II User Group

Dear dealii community,


I am trying to use "DoFTools::map_dofs_to_support_points" function for a two-field problem. The components of the first (vectorial) field are approximated using "FE_Q" and the second (scalar) field with "FE_DGPMonomial". The "component_mask" parameter is used to select only the dofs associated with the first field.


However I get the following error message (debug mode):

    void dealii::DoFTools::internal::(anonymous namespace)::map_dofs_to_support_points(const hp::MappingCollection<dim, spacedim> &, const DoFHandler<dim, spacedim> &, std::map<types::global_dof_index, Point<spacedim>> &, const dealii::ComponentMask &) [dim = 2, spacedim = 2]

The violated condition was:

    (fe_collection[fe_index].n_dofs_per_cell() == 0) || (fe_collection[fe_index].has_support_points())

Additional information:

    You are trying to access the support points of a finite element that

    either has no support points at all, or for which the corresponding

    tables have not been implemented.


From inspecting the method responsible for throwing the error is appears that all fields in the collection are checked via the Assert. Is the method applicable to problems where some fields have support points and others don't but are masked? The mask appears to only be appropriate to mask fields with support points.


Please find attached a minimal code for a two-field formulation that demonstrates the problem.

test_mdtsp.zip

Wolfgang Bangerth

unread,
Feb 7, 2023, 4:11:15 PM2/7/23
to dea...@googlegroups.com

Masoud,

> I am trying to use "DoFTools::map_dofs_to_support_points" function for a
> two-field problem. The components of the first (vectorial) field are
> approximated using "FE_Q" and the second (scalar) field with
> "FE_DGPMonomial". The "component_mask" parameter is used to select only
> the dofs associated with the first field.
>
>
> However I get the following error message (debug mode):
>
> /void dealii::DoFTools::internal::(anonymous
> namespace)::map_dofs_to_support_points(const hp::MappingCollection<dim,
> spacedim> &, const DoFHandler<dim, spacedim> &,
> std::map<types::global_dof_index, Point<spacedim>> &, const
> dealii::ComponentMask &) [dim = 2, spacedim = 2]/
>
> /The violated condition was: /
>
> /(fe_collection[fe_index].n_dofs_per_cell() == 0) ||
> (fe_collection[fe_index].has_support_points())/
>
> /Additional information: /
>
> /You are trying to access the support points of a finite element that/
>
> /either has no support points at all, or for which the corresponding/
>
> /tables have not been implemented./
>
>
> From inspecting the method responsible for throwing the error is
> appears that all fields in the collection are checked via the Assert. Is
> the method applicable to problems where some fields have support points
> and others don't but are masked? The mask appears to only be appropriate
> to mask fields with support points.

You are correct: The DoFTools::map_dofs_to_support_points first asks the
finite element object for its support points, and only then applies the
mask. It doesn't know how to break down the element to its base/
component elements if only some of the base elements have support
points, and the error message you see is the result.

I don't have an alternative suggestion either. You'd have to implement
this kind of functionality (i.e., looking through the base elements of a
FESystem) if you really need the map from DoF indices to support points.

Best
W.

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

Masoud Ahmadi

unread,
Feb 11, 2023, 6:09:01 AM2/11/23
to deal.II User Group
Thanks for the reply.
My goal of using this function is simply to get the DoF of some specific nodes which I know their coordinates; so, are there any alternative options that I won't have to use this function?

Bests
MA

Wolfgang Bangerth

unread,
Feb 13, 2023, 12:15:35 PM2/13/23
to dea...@googlegroups.com
On 2/11/23 04:09, Masoud Ahmadi wrote:
>
> My goal of using this function is simply to get the DoF of some specific nodes
> which I know their coordinates; so, are there any alternative options that I
> won't have to use this function?

In order to determine the DoF indices for a given vertex, you can use
cell->vertex_dof_index(v,i)
where v=the number of the vertex of the cell, i=the i'th DoF on that vertex.

It is a bit more complicated if you also have DoFs on edges or faces.
Reply all
Reply to author
Forward
0 new messages