Hello,

I am dealing with a smilar question as Shahin, i.e. computing global
L2-projections and approximations of mass-matrices with a lumped
integration. So I thought it makes sense to continue in this chat.

I produced two fields via L2-projection:
-In the first case I computed the mass matrix as a matrix, i.e.
(phi_i,phi_j).
-In the second case I approximated the mass matrix  with a lumped
integration. Instead of storing a mass matrix I computed a "mass-vector" by
summing up all off-diagonal elements in each row of the corresponding
mass-matrix and adding it do the corresponding diagonal. From the start I
directly used a vector instead of a matrix.

My distribution of the field (stresses) in terms of critical peaks,... is
of course the same.
However I am a little bit surprised that the "lumped nodal values" result
in lower stresses at most places. I actually expected the other way round,
because:
When thinking on a simple bar-element with 2 nodes, a lumped integration
concentrates the total mass of the bar to the nodes so that the inertia
becomes bigger compared with other Quadrature formulas.

So intuitively speaking, what is the effect of the lumped integration? From
a mathematical viewpoint there are no more couplings between the DoFs. But
does this lead to lower nodal values, higher nodal values or is this
context dependent?

Best
Simon


Am Di., 25. Mai 2021 um 11:52 Uhr schrieb Shahin Heydari <
sh1990.heyd...@gmail.com>:

> Thank you so much Sir, appreciate your help.
>
> Regards
> Shahin
>
> On Mon, May 24, 2021 at 5:46 PM Wolfgang Bangerth <bange...@colostate.edu>
> wrote:
>
>> On 5/21/21 1:52 PM, Shahin Heydari wrote:
>> > m_i = sigma (m_ij) = int(phi_i dx)
>>
>> Right. You compute the integral as a sum over the integrals on each cell.
>> So
>> when you do
>>
>>    for (cell = ...)
>>      for (i=...)
>>        {
>>          lump = 0;
>>          for (q=...)
>>            lump += ...
>>
>>          lump_cell_matrix(i,i) = lump;
>>        }
>>
>> then this is correct if you then add lump_cell_matrix(i,i) to your global
>> matrix.
>>
>> You mention in your original question that you are unsure about what to
>> do
>> with Q2 (or higher) elements. This is a difficult question, but there is
>> a
>> good amount of literature on the topic. Just search for "mass lumping
>> higher
>> order elements".
>>
>> Best
>>   W.
>>
>> --
>> ------------------------------------------------------------------------
>> Wolfgang Bangerth          email:                 bange...@colostate.edu
>>                             www: http://www.math.colostate.edu/~bangerth/
>>
>> --
>> 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/1b0908b3-e43c-2d6a-b2a0-09a7a51b4636%40colostate.edu
>> .
>>
> --
> 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/CALUbCFwyACZFLqF2ZEXztJ3H4b2iEhcVGYJ90beAMNPnUVz6aQ%40mail.gmail.com
> <https://groups.google.com/d/msgid/dealii/CALUbCFwyACZFLqF2ZEXztJ3H4b2iEhcVGYJ90beAMNPnUVz6aQ%40mail.gmail.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 dealii+unsubscr...@googlegroups.com.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/CAM50jEtwjTj-BF58cS%3Dn%3DL0E52tFNfbrWYTzspJRnaFRuszfbw%40mail.gmail.com.

Reply via email to