I've directed the developers towards this thread, but I'll make a comment 
about it.

On Wednesday, June 1, 2016 at 4:17:22 PM UTC+2, RAJAT ARORA wrote:
>
> Hi Jean,
>
> I am just thinking about it because there may be possibility of bug there 
> as well if reinit(cell) and new fe_values object have something in common 
> like catching data or anything.
> Even I am not sure about the internal working of the library, so I am just 
> making a guess.
>
> Also, in the bug report, we haven't mentioned that the problem does not 
> show up if we use a new fe_values object instead of the old one.
>
> I think we should make a comment about this as well. What do you say ? 
>
>
>
>
> On Wednesday, June 1, 2016 at 9:00:45 AM UTC-4, Jean-Paul Pelteret wrote:
>>
>> You're welcome. 
>>
>> As to why creating a new FEValues object solves the problem - perhaps it 
>> has something to do with the way that FEValues caches data. When you have 
>> just a single FEValues, it references "state" information, but when you 
>> create a new one the equivalent data is completely reevaluated. However, I 
>> can't be certain as I'm not familiar with the inner working of these 
>> classes.
>>
>> On Wednesday, June 1, 2016 at 2:20:35 PM UTC+2, RAJAT ARORA wrote:
>>>
>>> Hello Jean,
>>>
>>> Thanks for your help. It would have been difficult for me to find out 
>>> without your support.
>>>
>>> Also, just to clarify, if the bug is indeed with Dof_handler then, why 
>>> making a new fe_values object (after mesh movement) in my code gave 
>>>  correct gradients ?
>>>
>>>
>>> On Wednesday, June 1, 2016 at 2:10:45 AM UTC-4, Jean-Paul Pelteret wrote:
>>>>
>>>> Dear Rajat,
>>>>
>>>> Thanks for posting an example along the lines of your implementation. 
>>>> Well, from the output you've shown there, I'm not surprised that there's 
>>>> some problem. If I'd thrown in some asserts then it would have been more 
>>>> clear, but the two volume calculations "Using GridTools" and "Calculated" 
>>>> should always be equal. In your test case, this is not true after the mesh 
>>>> movement. Obviously if this is not the case, then the geometry change is 
>>>> not being propagated to the DoFHandler, and therefore any subsequent 
>>>> gradient computations will not work correctly.
>>>>
>>>> Anyway, after some sleuthing I've found that this must in fact be a 
>>>> subtle bug in deal.II. If looks like you need to refine the grid at least 
>>>> once in order for the movement of cell vertices to be registered in the 
>>>> DoFHandler (see attached file, lines 103 and 104). With no refinement the 
>>>> asserts fail, but with at least one level of refinement they pass. I can 
>>>> confirm that the same is the case in your example.
>>>>
>>>> I've opened up a bug report on your behalf here 
>>>> <https://github.com/dealii/dealii/issues/2661>. So, in summary, you'll 
>>>> just have to wait for a bit until the problem is fixed, or use at least 
>>>> one 
>>>> level of global refinement :-)
>>>> J-P
>>>>
>>>> On Wednesday, June 1, 2016 at 2:04:19 AM UTC+2, RAJAT ARORA wrote:
>>>>>
>>>>> Hello Jean,
>>>>>
>>>>> First of all, thanks a lot for taking time to write to me and 
>>>>> explaining this issue.
>>>>> I am also very keen on solving this issue. 
>>>>>
>>>>> I looked at the code you sent and indeed it looks like what you say is 
>>>>> correct. 
>>>>> However, If I calculate the gradients, it is giving wrong results and 
>>>>> giving gradients on the mesh before movement.
>>>>>
>>>>> I took the gradients before and after moving the mesh and it is giving 
>>>>> me the same values which has been happening with me.
>>>>>
>>>>> I am also attaching my file and CMakelists.txt for you to review. 
>>>>>
>>>>> I am running in debug mode but the issue is still unresolved.
>>>>>
>>>>> I am taking the gradient of "f" which I initialized to the material 
>>>>> coordinates.
>>>>>
>>>>> Is this a bug in the library or the version that I am using (8.3.2) ?
>>>>>
>>>>>
>>>>> gradient at qpoint 0 before mesh movement:: 1 0 0 1
>>>>> Volume before mesh movement:   Using GridTools: 1  Calculated: 1
>>>>>
>>>>> gradient at qpoint 0 after mesh movement :: 1 0 0 1
>>>>> Volume after mesh movement:   Using GridTools: 0.95  Calculated: 1
>>>>>
>>>>> gradient at qpoint 0 after mesh movement and with new fe_values object 
>>>>> :: 1.05263 0 0 1
>>>>>
>>>>>
>>>>>  
>>>>>
>>>>>
>>>>> On Tuesday, May 31, 2016 at 3:19:13 PM UTC-4, Jean-Paul Pelteret wrote:
>>>>>>
>>>>>> Hi Rajat,
>>>>>>
>>>>>> Without your code, its difficult to say what's going wrong here. It 
>>>>>> should (and does) work as I've explained. In summary, FEValues gets to 
>>>>>> know 
>>>>>> about any changes you've made to the underlying because its reinit() 
>>>>>> call 
>>>>>> takes an *active* cell which always knows its geometry (you 
>>>>>> manipulate this directly when you move the mesh).
>>>>>>
>>>>>> I've written a small test case to prove to myself that your "code 1" 
>>>>>> should work as expected, and indeed this is the case. I've attached it 
>>>>>> for 
>>>>>> you to investigate, and this is the expected result:
>>>>>> Volume before mesh movement:   Using GridTools: 1  Calculated: 1
>>>>>> Volume after mesh movement:   Using GridTools: 0.75  Calculated: 0.75
>>>>>>
>>>>>> One thing that *might* be the problem is that you might be iterating 
>>>>>> over inactive cells with "dof_handler.begin()" instead of 
>>>>>> "dof_handler.begin_active()". If you only have a single cell and 
>>>>>> accidentally do this, then even in debug mode you do get a result but it 
>>>>>> is 
>>>>>> incorrect:
>>>>>> Volume before mesh movement:   Using GridTools: 1  Calculated: 1
>>>>>> Volume after mesh movement:   Using GridTools: 0.75  Calculated: 1
>>>>>>
>>>>>> However, once you perform some refinement then this error would be 
>>>>>> picked up in debug mode (are you running your tests in debug mode?).
>>>>>>
>>>>>> Does that help? If you solve the issue, then I'd be interested to 
>>>>>> know what the problem was.
>>>>>> J-P
>>>>>>
>>>>>> On Tuesday, May 31, 2016 at 4:17:56 PM UTC+2, RAJAT ARORA wrote:
>>>>>>>
>>>>>>> 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