projectin of gradients of the discrete solution to nodes - two different discontinuous approaches

66 views
Skip to first unread message

Simon

unread,
May 22, 2021, 3:22:09 AM5/22/21
to deal.II User Group

Dear all,

I am projecting gradients of my displacement solution, for instance stresses, which are only available to me at the quadrature points, to the nodes in order to get fields. (However I project many more variables, so I transformed my problem in projecting many scalar variables, but anyway this is not important regarding my question).

-I am doing the projection once by solving a global minimization problem, i.e. minimizing the integral over the whole domain of the squared difference diff = "values at qps - values at nodes". This leads to the typical mass matrix (phi_i, phi_j) which is also provided in the VectorTools namespace.
I solve this min.-problem once with a continuous FE_Q und once with a discontinuous FE_DGQ, i.e. I call the very same function body, I just give different DoFHandlers as input.

-A second approach for the same discontinuous FE_DGQ is intruduced in step-18 (extensions) by using compute_projection_from_quadrature_points_matrix(), i.e. a local method. (My poly_degree is constructed in a way that I have always as many dofs_per_cell as qps.)

My question is if the output, i.e. the nodal values, of both discontinuous approaches
(1. given a FE_DGQ to the global min.-problem vs. 2. the local compute_projection_from_quadrature_points_matrix() approach)
should be the same, except of numerical errors?

I am aware that I solve a linear system in the first case whereas in the second case not, but my idea was the following: The first approach minimizes the squared difference, but if I have as many dofs_per_cell as qps then of course the squared difference can be minimized to zero for each qp, i.e. the qp values can directly be transferred to the nodes. And if my understanding is correct, this is exactly what the second approach does. So in the first approach I do not do it "directly" but due to its definition and the number of dofs it should do the same as second.
Of course I compared this approaches in my program. However depending on the mesh size there are deviations up to percent, with finer meshes this difference reduces.
Can this deviation be argumented away with the standard argument "this is the numerics..." or is there is a mathematical difference and both approaches do something different?

Thanks for the input!

Best
Simon

Wolfgang Bangerth

unread,
May 24, 2021, 11:54:05 AM5/24/21
to dea...@googlegroups.com
On 5/22/21 1:22 AM, Simon wrote:
>
> I am aware that I solve a linear system in the first case whereas in the
> second case not, but my idea was the following: The first approach minimizes
> the squared difference, but if I have as many dofs_per_cell as qps then of
> course the squared difference can be minimized to zero for each qp, i.e. the
> qp values can directly be transferred to the nodes. And if my understanding is
> correct, this is exactly what the second approach does. So in the first
> approach I do not do it "directly" but due to its definition and the number of
> dofs it should do the same as second.
> Of course I compared this approaches in my program. However depending on the
> mesh size there are deviations up to percent, with finer meshes this
> difference reduces.
> Can this deviation be argumented away with the standard argument "this is the
> numerics..." or is there is a mathematical difference and both approaches do
> something different?

For discontinuous elements, the mass matrix you compute is block diagonal and
the projection can be computed cell-by-cell -- so it shouldn't make a
difference whether you solve the global problem or do it one cell at a time.
If you get different results, it would be useful to investigate how exactly
they are different.

For continuous elements, the projection is not exact in general, and whatever
you compute locally on one cell has to be reconciled with what you compute on
neighboring cells. The global project is one way to do this reconciliation.

Best
W.


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

Simon Wiesheier

unread,
May 25, 2021, 3:41:45 AM5/25/21
to dea...@googlegroups.com
Thank you!

I failed to see that I used the FE_DGPMonomial<dim> instead of FE_DGQ<dim>. With the latter it makes as you said acutally no difference.

However I am wondering why FE_DGPMonomial<dim> gives my quite good results as well. As I said, for "fine meshes" the results of the local and global approach are nearly the same and for coarse meshes, i.e. just one or two times globally refined meshes, the biggest deviation also reaches not more than 5%.
I didn´t have to use FE_DGPMonomial so far. I´ve seen that for degree=1 this element has only 3 dofs per cell for a scalar problem in 2D, i.e I can not extrapolate the 4 qp-values to the nodes. But let´s assume I use this element for my local approach.
-What is the relations of the resulting 3 nodel values to the 4 qp-values? Is it still something like minimizing the squared difference or is the idea a different one?
-And why is there a difference at all between the local and global approach? Should the global approach not result in a block-diagonal matrix as well?

I couldn´t figure out the purpose of the disontinuous legendre-polynomials by having a look at the documentation.

-Simon

--
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+un...@googlegroups.com.
To view this discussion on the web visit https://groups.google.com/d/msgid/dealii/e91a3cb7-d193-b5e4-850f-900ca6a5dea3%40colostate.edu.

Wolfgang Bangerth

unread,
May 25, 2021, 3:30:33 PM5/25/21
to dea...@googlegroups.com
On 5/25/21 1:41 AM, Simon Wiesheier wrote:
>
> However I am wondering why FE_DGPMonomial<dim> gives my quite good results as
> well. As I said, for "fine meshes" the results of the local and global
> approach are nearly the same and for coarse meshes, i.e. just one or two times
> globally refined meshes, the biggest deviation also reaches not more than 5%.
> I didn´t have to use FE_DGPMonomial so far. I´ve seen that for degree=1 this
> element has only 3 dofs per cell for a scalar problem in 2D, i.e I can not
> extrapolate the 4 qp-values to the nodes. But let´s assume I use this element
> for my local approach.
> -What is the relations of the resulting 3 nodel values to the 4 qp-values? Is
> it still something like minimizing the squared difference or is the idea a
> different one?

Yes, that's what the *projection* does: It minimizes the difference in some
kind of norm.


> -And why is there a difference at all between the local and global approach?
> Should the global approach not result in a block-diagonal matrix as well?

The projection is local (=the matrix is block diagonal) if the space is
discontinuous. That is, you should get the same result whether you do the
local or the global projection for both FE_DGQ and FE_DGPMonomial, along with
all other discontinuous elements.

If you get different results between local and global projection for
FE_DGPMonomial, it would be interesting to investigate why.
Reply all
Reply to author
Forward
0 new messages