On 03/14/2018 01:39 PM, Jane Lee wrote:
Hi All,

Just a quick question on convergence rates/error analysis.

I have a solver which does primarily what they do in step-20. The model I am trying to solve takes the pressure solution from the equations in step-20 and uses its gradient on the right hand side of a new equation.

I use
fe_values_p[pressure].get_function_gradients (solution_darcy, grad_p_values);

on my test problem, where i just do v=grad p (v is a vector), to solve after the pressure is found from step-20 so that with the mass matrix as the system matrix, on the rhs i have: local_rhs(i) += (grad_pf_values[q_point][component_i] * fe_values_new.shape_value(i,q_point)*fe_values_new.JxW(q_point);

I have noticed that convergence rates are very low, at p=1 order, convergence rate is between 2 and 3, when it should be 4. Further, i noticed that increasing the degree of the test problem also did not change the error at all.

So you are saying that

  || v - v_h ||  =  O(h)

or put differently,

  || grad p - v_h ||  =  O(h)

?

I am using
fe (FE_RaviartThomas <https://www.dealii.org/8.4.0/doxygen/deal.II/classFE__RaviartThomas.html><dim>(degree1), 1, FE_DGQ <https://www.dealii.org/8.4.0/doxygen/deal.II/classFE__DGQ.html><dim>(degree1), 1)
as normal as in step-20 for the u and pressure variable,
and for my test problem, for the solution v, i am using
FE_Q<dim>(degree2)

What are the values for degree1 and degree2 here? If degree1 is low, then choosing degree2 high will not fix anything: If grad p_h is already wrong, then projecting it to a high order space isn't going to bring back any accuracy previously lost.

Have you evaluated the error in p independently? What convergence rate do you get for it? What rate do you get for grad p? You can't expect the rate for v to be any better than the rate for p_h.


So my main question is, how much error is get_function_gradients actually accounting for?

None. It evaluates a finite element field at quadrature points exactly.

Best
 W.

--
------------------------------------------------------------------------
Wolfgang Bangerth          email:                 [email protected]
                           www: http://www.math.colostate.edu/~bangerth/

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

Reply via email to