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.
##
#  CMake script for the step-1 tutorial program:
##

# Set the name of the project and target:
SET(TARGET mesh_movement)

# Declare all source files the target consists of:
SET(TARGET_SRC
  ${TARGET}.cc
  # You can specify additional files here!
  )

# Usually, you will not need to modify anything beyond this point...

CMAKE_MINIMUM_REQUIRED(VERSION 2.8.8)

FIND_PACKAGE(deal.II 8.4 QUIET
  HINTS ${deal.II_DIR} ${DEAL_II_DIR} ../ ../../ $ENV{DEAL_II_DIR}
  )
IF(NOT ${deal.II_FOUND})
  MESSAGE(FATAL_ERROR "\n"
    "*** Could not locate deal.II. ***\n\n"
    "You may want to either pass a flag -DDEAL_II_DIR=/path/to/deal.II to 
cmake\n"
    "or set an environment variable \"DEAL_II_DIR\" that contains this path."
    )
ENDIF()

DEAL_II_INITIALIZE_CACHED_VARIABLES()
PROJECT(${TARGET})
DEAL_II_INVOKE_AUTOPILOT()
#include <deal.II/dofs/dof_handler.h>
#include <deal.II/fe/fe_q.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/grid_out.h>
#include <deal.II/grid/grid_tools.h>
#include <deal.II/grid/tria.h>
#include <deal.II/grid/tria_accessor.h>
#include <deal.II/grid/tria_iterator.h>

#include <iostream>
#include <fstream>

using namespace dealii;

template<int dim>
double
integral_volume (const DoFHandler<dim> &dof_handler,
                 FEValues<dim> &fe_values)
{
	double vol = 0.0;
	for (typename DoFHandler<dim>::active_cell_iterator
			 cell  = dof_handler.begin_active();
			 cell != dof_handler.end();
			 ++cell)
	{
		fe_values.reinit(cell);
		for (unsigned int qp_cell=0; qp_cell<fe_values.get_quadrature().size(); ++qp_cell)
			vol += fe_values.JxW(qp_cell);
	}
	return vol;
}

int main()
{
	const int dim = 2;
	Triangulation<dim> tria;
	const int degree = 1;
  FE_Q<dim> fe(degree);
	QGauss<dim> qf_cell(degree+1);
	DoFHandler<dim> dof_handler(tria);

  // Single cell of unit length
	GridGenerator::hyper_cube(tria, 0.0, 1.0);
	tria.refine_global(1);
	dof_handler.distribute_dofs(fe);
	
	FEValues<dim> fe_values (fe, qf_cell, 
		                       update_JxW_values);
	
	// Compute volume of the grid by integration
	{
		std::cout 
		  << "Volume before mesh movement: " 
			<< "  Using GridTools: " << GridTools::volume(tria)
			<< "  Calculated: " << integral_volume(dof_handler, fe_values) 
			<< std::endl;
		std::ofstream output_file ("mesh-initial.eps");
		GridOut().write_eps(tria,output_file);
	}
		
		
	
	// Move one of the grid vertices
	for (typename Triangulation<dim>::active_cell_iterator
			 cell  = tria.begin_active();
			 cell != tria.end();
			 ++cell)
	{
		// Decrease the x-length of the geometry by 25%
		// By moving only the right-most vertices
		// Note: This distorts only the right two
		// cells
		for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
		  if (std::abs(cell->vertex(v)[0] - 1.0) < 1e-9)
			  cell->vertex(v)[0] *= 0.75;
	}
	
	// Now we measure the volume again. 
	// Since the actual geometry has been updates
	// (i.e. the "reference" configuration has changed)
	// the computed volume will change as well
	{
		std::cout 
		  << "Volume after mesh movement: " 
			<< "  Using GridTools: " << GridTools::volume(tria)
			<< "  Calculated: " << integral_volume(dof_handler, fe_values) 
			<< std::endl;
		std::ofstream output_file ("mesh-final.eps");
		GridOut().write_eps(tria,output_file);
	}

	return 0;
}

Reply via email to