Hello Jean,

Thanks for your reply.

In my case, since I am physically moving the mesh ( as I was discussing here 
<https://groups.google.com/forum/#!topic/dealii/8I21UqIWAmU>), even after 
calling 
fe_values.reinit(cell), I am getting the gradients on the geometry before 
mesh movement which should not be the case.



Also, note that once I move the mesh, this happens only if I use the same 
fe_values object. But if I create a new FEValues object (lets call it 
fe_values_new), then if I use fe_values_new.get_function_gradients(), I get 
the gradients on the updated geometry. Code given below for reference.


Code 1: (Does not work)
// Everything is in one big function

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.reinit(cell);

fe_values.get_function_gradients();



Code 2:(It works)
// Everything is in one big function

FEValues<dim> fe_values  = new FEValues<dim>(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)

delete fe_values;

fe_values = new FEValues<dim> (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>
);

fe_values.reinit(cell);

fe_values.get_function_gradients();






Does that mean the instant FEValues object is created, it takes that as a 
reference configuration and do all the gradients on this reference 
configuration ?

I also read the posts related to this issue here 
<https://groups.google.com/forum/#!searchin/dealii/mesh$20fevalues/dealii/KvqA-fQAYgM/yQdzNZsRSOUJ>
 and here 
<https://groups.google.com/forum/#!searchin/dealii/mesh$20fevalues/dealii/B2v7qSLrnck/0e44LC4bCQAJ>
.
 

 


On Tuesday, May 31, 2016 at 1:41:10 AM UTC-4, Jean-Paul Pelteret wrote:
>
> 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