Re: [deal.II] Query regarding DoFTools::dof_indices_with_subdomain_association()

2019-10-09 Thread Daniel Arndt
Vachan,

[...]
>
However, the documentation of the former function states
>
> Note that this function is of questionable use for DoFHandler objects
>> built on parallel::distributed::Triangulation since in that case ownership
>> of individual degrees of freedom by MPI processes is controlled by the DoF
>> handler object, not based on some geometric algorithm in conjunction with
>> subdomain id. In particular, the degrees of freedom identified by the
>> functions in this namespace as associated with a subdomain are not the same
>> the DoFHandler class identifies as those it owns.
>
>
> What is the meaning of the second line? I suppose this has something to do
> with the FiniteElement type of the dof handler. Does this mean to say that
> irrespective of whether the FiniteElement has dofs on cell faces or not,
> this function approves a dof if it is geometrically residing on a face?
> This behaviour is actually what I want because although FE_DGQ element
> doesn't attach any dofs to cell faces, I want to get dof indices of the
> dofs (geometrically) lying on the face.
>

DoFTools::dof_indices_with_subdomain_association returns the degrees of
freedom of all the cells that have the given subdomain id. For
parallel::distributed::Triangulation objects the subdomain id is the the
MPI rank and this is the only valid input.
In this case, the function simply returns all the degrees of freedom on
locally owned cells. If a finite element defines degrees of freedoms that
are associated with a face, these degrees of freedom are part of all the
IndexSet objects in case the face is located on the interface between two
subdomains. For a parallel::distributed::Triangulation object, these are
the degrees of freedom on the faces to ghosted cells.

FE_DGQ doesn't define any face degrees of freedom. All its degrees of
freedom are defined to be inside of a cell. You could still ask for the
unit support points and figure out which of them are geometrically on a
given face.
However, note that as soon as gradients are involved all the degrees of
freedom contribute to values on faces.

Best,
Daniel

-- 
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 dealii+unsubscr...@googlegroups.com.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/CAOYDWbLMDhoL4frLmb50yuLxOgoJnMMwmHh74Yp1uhkv%2B7Z_3Q%40mail.gmail.com.


[deal.II] Query regarding DoFTools::dof_indices_with_subdomain_association()

2019-10-09 Thread vachan potluri
Hello all,

I am writing an MPI parallel DG code for linear advection equation, just 
for understanding the parallel programming paradigm of deal.II.

To evaluate numerical flux at a face connecting two subdomains, I only 
require the solution (from a different mpi process) at dofs which lie on 
the face; the dofs which lie inside the neighboring cell (of the neighbor 
subdomain) are not required. Consequently, I think the ghost indices for my 
solution vector would be those provided by the function 
DoFTools::dof_indices_with_subdomain_association() 
,
 
rather than DoFTools::extract_locally_relevant_dofs() 
.
 
However, the documentation of the former function states

Note that this function is of questionable use for DoFHandler objects built 
> on parallel::distributed::Triangulation since in that case ownership of 
> individual degrees of freedom by MPI processes is controlled by the DoF 
> handler object, not based on some geometric algorithm in conjunction with 
> subdomain id. In particular, the degrees of freedom identified by the 
> functions in this namespace as associated with a subdomain are not the same 
> the DoFHandler class identifies as those it owns.


What is the meaning of the second line? I suppose this has something to do 
with the FiniteElement type of the dof handler. Does this mean to say that 
irrespective of whether the FiniteElement has dofs on cell faces or not, 
this function approves a dof if it is geometrically residing on a face? 
This behaviour is actually what I want because although FE_DGQ element 
doesn't attach any dofs to cell faces, I want to get dof indices of the 
dofs (geometrically) lying on the face.

Thanks!
Vachan

-- 
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 dealii+unsubscr...@googlegroups.com.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/f1a506a4-f652-43fd-9b02-c948f1eca49d%40googlegroups.com.