If you physically move the mesh vertices, like you were discussing here 
<https://groups.google.com/forum/#!topic/dealii/8I21UqIWAmU>, then when you 
call fe_values.reinit(cell), the mapping between the current cell geometry 
(as defined by the triangulation) and isoparametric cell are accounted for 
and used when computing gradient quantities. 

For example, compare the update_quadrature_point_history functions in 
step-18 
<https://dealii.org/8.4.0/doxygen/deal.II/step_18.html#TopLevelupdate_quadrature_point_history>
 
and step-44 
<https://dealii.org/8.4.0/doxygen/deal.II/step_44.html#Solidupdate_qph_incremental>.
 
In step-44 the mesh is not physically moved, so "get_function_gradients" 
applied to the total displacement vector will return elements of 
\frac{\partial \mathbf{u}}{\partial \mathbf{X}} (i.e. the material 
gradient), while in step-18 the mesh is moved and "get_function_gradients" 
returns elements of \frac{\partial \mathbf{u}}{\partial \mathbf{x}} (i.e. 
the spatial gradient). However, you will note that in step-18, only the 
incremental displacement gradient is computed as the stresses are 
systematically incremented using this result.

On Tuesday, May 31, 2016 at 2:54:32 AM UTC+2, RAJAT ARORA wrote:
>
> Hello all,
>
> I have a question about FEValues<dim> object that we use in the 
> computation of element stiffness and other places.
>
>
> If I declare fe_values object and get then move the mesh and then try to 
> call fe_values.get_function_gradients()
> will this give me gradients on the deformed mesh or the mesh that was 
> chosen to begin with ? (Code structure below)
>
>
> FEValues<dim> fe_values (fe, quadrature_formula,
>  update_values 
> <https://dealii.org/8.4.1/doxygen/deal.II/group__feaccess.html#ggaa94b67d2fdcc390690c523f28019e52fa4057ca2f127aa619c65886a9d3ad4aea>
>  
> | update_gradients 
> <https://dealii.org/8.4.1/doxygen/deal.II/group__feaccess.html#ggaa94b67d2fdcc390690c523f28019e52facbcc430975fa6af05f75ca786dc6fe20>
>  
> |
>  update_quadrature_points 
> <https://dealii.org/8.4.1/doxygen/deal.II/group__feaccess.html#ggaa94b67d2fdcc390690c523f28019e52fad5c9ff886b9615349a5d04a6c782df4a>
>  
> | update_JxW_values 
> <https://dealii.org/8.4.1/doxygen/deal.II/group__feaccess.html#ggaa94b67d2fdcc390690c523f28019e52fa714204722e9eeb43aadbd0d5ddc48c85>
> );
>
> // do some stuff
>
> // move mesh like step 17, adding displacements to cell->vertex(i)
>
> fe_values.get_function_gradients();
>
>
>
>
>
> If it gives gradients on the initial mesh, then how does it know those 
> coordinates once you change the coordinates of the 
> vertices of the triangulation ?
>
>  
>

-- 
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.
For more options, visit https://groups.google.com/d/optout.

Reply via email to