> Unfortunately, visualizing is not enough to my purposes. I need the 
> tangential component of the gradient at the mesh faces to perform a 
> post-processing phase where I compute some relevant quantities to be used 
> later in solving a different problem.
>
In general, you can create a quadrature rule that has support points at 
face centers on the reference cell [0,1]^{dim}, i.e.
  std::vector<Point<dim-1> > reference_face_centers(1);
  reference_face_center[0] = Point<1> (.5) //2d
  reference_face_center[0] = Point<1> (.5, .5) //3d
  //fill reference_face_centers
  Quadrature<dim>(reference_face_centers);
Then, create at FEFaceValues object using that Quadrature and a vector to 
store the gradient:
  FEValues<dim> fe_dummy_values(fe, reference_face_centers, 
update_gradients | update_normal_vectors);
  std::vector<Tensor<1,dim> > gradient (1);
and loop over all cells and faces as shown in step-7[1]. Finally, you can 
query the location of the quadrature point and the gradient via
  fe_face_values.reinit (cell, face_number);
  fe_face_values.get_function_gradients(solution, gradient)
  std::cout << "The gradient at " << fe_face_values.quadrature_point(0) << 
" is " << gradient[0] << std::endl;
  std::cout << "The normal vector at " << 
fe_face_values.quadrature_point(0) << " is " << 
fe_face_values.normal_vector(0) << std::endl;

Alternatively, you can compute the projection of your discrete solution 
into the FE_Nedelec<dim>(0) space by solving for u satisfying
(v,u)=(v, \nabla solution) for all v in FE_Nedelec<dim>(0)
For doing that, you would need to define an additional FE_Nedelec object, 
create a mass matrix and evaluate the gradient of your original solution
at the quadrature points you are using for building the mass matrix via 
fe_values.get_function_gradients(...).

Which way the more appropriate one is depends on how exactly you want to 
perform postprocessing.

Best,
Daniel

[1] https://www.dealii.org/8.4.0/doxygen/deal.II/step_7.html

-- 
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