On 4/6/22 05:46, Léonhard YU wrote:
The number of degree of freedom of 8-node cell is 8 according to the
reference book, but in deal.II it is 8*3.
I do not understand that. I guess that each column in matrix S is for
one DoF, but I am not sure.
I will point you again at the vector-valued documentation module that
has many examples. In your case, you have x-, y-, and z-displacements at
each vertex of the cell, so you have 8*3 degrees of freedom on each cell.
The cell stiffness and mass matrix are:
for (unsigned int i = 0; i < dofs_per_cell; ++i)
for (unsigned int j = 0; j < dofs_per_cell; ++j)
for (unsigned int q_point = 0; q_point < n_q_points;
++q_point)
{
Tensor<1,6> eps_u_i;
This suggests that you are using Voigt notation. This is common in the
engineering literature, but not typically the perspective we take in
(the more mathematically oriented) deal.II. Rather, we tend to use
notation that emphasizes that eps_u_i is a rank-2 symmetric tensor, and
so use SymmetricTensor<2,3> instead.
It is not inherently wrong to do as you do but you will find more
examples in the deal.II tutorials if you use tensor notation.
cell_stiffness_matrix(i, j) *= eps_u_i * stress_strain_D
* eps_u_j * fe_values.JxW(q_point);
I am certain that you want to use
> cell_stiffness_matrix(i, j) += eps_u_i * stress_strain_D
> * eps_u_j * fe_values.JxW(q_point);
instead (note the += instead of *=).
if (fe.system_to_component_index(i).first ==
fe.system_to_component_index(j).first)
{
cell_mass_matrix(i, j) +=
rho *
fe_values.shape_value(i, q_point) *
fe_values.shape_value(j, q_point) *
fe_values.JxW(q_point);
}
}
This looks correct.
Now it reports convergence only after 1 step, the eigenvalues are all
drop in the interval of spurious eigenvalues.
Yes. That's because the mistake above means that your stiffness matrix
entries are probably all zero, with the exception of those you set to
something else by hand.
Best
W.
--
------------------------------------------------------------------------
Wolfgang Bangerth email: [email protected]
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 [email protected].
To view this discussion on the web visit
https://groups.google.com/d/msgid/dealii/4ee2593a-e91c-14a8-e2a7-2e12c5b8eaff%40colostate.edu.