Hi,
I just recently did something similar and as has been pointed out, there 
are as always multiple ways of achieving this.
Since I wanted to not only output the first principal stress, but actually 
use it further for my calculations, I wanted to have a DoF vector 
resembling the field in my mesh. Thus, there are again (at least) 2 
options: 
i) introduce a discontinuous field with corresponding DoFs (would require 
more effort since I am only using continuous FEs)
ii) use a continuous, averaged field, accepting the introduced error 
(constant per element and averaged over neighbors)
If you want to do ii) as well, it boils down to a classical assembly loop,  
and one computes the Cauchy stress tensor (or whatever stress tensor you 
want to get the principal stresses from) as you did in your momentum 
balance equation.
This might look something like the following in integration point q:
{
                            Tensor<2,dim> F_s = identity_tensor + 
gradients_displacement[q];  // displacement gradient
                            Tensor<2,dim> PK2_stress = 
get_PK2_stress<dim>(F_s); // second Piola-Kirchhoff stress tensor, you 
define this function based on your constitutive equation
                            Tensor<2,dim> cauchy_stress = 
symmetrize((1.0/determinant(F_s)) * F_s * PK2_s * transpose(F_s)); // 
Cauchy stress tensor via Piola transform for finite strain theory
                      const std::array<double,dim> eigenvalues_cauchy = 
eigenvalues(cauchy_stress);
                            double max_principal_stress = cauchy_stress[0];
}
Optionally, one could always get a global DoF vector via projecting the 
stress components or the first principal stress (or the quantity you need 
to plot). This would be much much more accurate than what I suggest in ii).
Going through i) or ii), both ways one ends up with (dis-)continuous global 
vectors, possibly averaged over neighbors. These are then easily exported 
via the classical dealii::DataOut<dim> 
data_out; data_out.add_data_vector(...).
If you want to just output the results, maybe try the 
dealii::DataPostprocessor<dim> as mentioned by Wolfgang.

Best regards
Richard


[email protected] schrieb am Donnerstag, 17. November 2022 um 16:21:30 
UTC+1:

> Hello everyone,
>
> It is the first time I am using the library Deal II. I am wondering, is 
> there any way/example to get the principal stress components as an output 
> of the simulation? Thank you very much in advance
>
> Raul
>

-- 
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/313e7875-b17f-4519-bb76-694a8bca16b5n%40googlegroups.com.

Reply via email to