Konstantin,

What is the correct way to set smooth coefficients in space for high p-order
of elements simulation in deal.ii?

E.g. in Step-6 nonconstant coefficient is a step function  alpha=20 for R<0.5
and alpha=1 otherwise (0<R<1).  How to set it to be a smooth function, e.g.
some polynomial function of order 'n' like alpha=20*(1-R^n)+1? So it should be
enough to use few high-order elements for the whole model to converge to the
smooth solution.

I`ve checked a number of tutorials
(e.g. 
http://www.dealii.org/developer/doxygen/deal.II/step_4.html#Righthandsideandboundaryvalues
, equation data section in step-27, and few others - all of them use a
point-wise approximation of an analytic function  (so it is tuned for
h-refinements), which is practically the same as a step function for a given
mesh.

No, that's not true. When you say "pointwise", what I think you really mean is that we only evaluate the coefficient at individual (quadrature) points. But that's good enough, because when you do quadrature -- say for an integrand on a cell K,

  A_{ij}^K  =  \sum_q  a(x_q) \nabla phi_i(x_q) \nabla phi_j(x_q) JxW(q)

then that is equivalent to computing the integral of a polynomial that is equal to

  a(x_q) \nabla phi_i(x_q) \nabla phi_j(x_q)

i.e., that *interpolates* these point values. This polynomial is of course smooth. In particular, this corresponds to a *smooth* interpolation of the coefficient a(x_q), even though you only evaluate it at quadrature points.

What kind of polynomial you use to interpolate the integrand depends on how many quadrature points you use. If you need a very accurate interpolation of the coefficient, just use a higher order Gauss formula.


So in some sense, you've got it exactly the wrong way around: If, in step-6, we use a discontinuous coefficient, then the finite element method using quadrature really sees only a polynomial approximation of that discontinuous coefficient on every mesh. (Though that polynomial approximation of course tends to a discontinuous function under mesh refinement.)

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