Jane,

> I continued to find out why I wasn't getting the correct applied Dirichlet 
> values on the boundary for a code very similar to step-20, where the 
> Dirichlet 
> condition is applied weakly using
> 
> for (unsigned int face_no=0;
> face_no<GeometryInfo<dim>::faces_per_cell;
> ++face_no)
> if (cell->at_boundary(face_no))
> {
> fe_face_values.reinit 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classFEFaceValues.html#af6e079ca7429d54433343d50bd334c3c>
>  
> (cell, face_no);
> pressure_boundary_values
> .value_list (fe_face_values.get_quadrature_points 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classFEValuesBase.html#a5f8732ebe2d3c6746f6de26a79cb1e45>(),
> boundary_values);
> for (unsigned int q=0; q<n_face_q_points; ++q)
> for (unsigned int i=0; i<dofs_per_cell; ++i)
> local_rhs(i) += -(fe_face_values[velocities].value (i, q) *
> fe_face_values.normal_vector 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classFEValuesBase.html#a130eea0fa89263d93b20521addc830c7>(q)
>  
> *
> boundary_values[q] *
> fe_face_values.JxW 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classFEValuesBase.html#ad097580a2f71878695096cc73b271b9d>(q));
> }
> 
> 
> I then looked at step-20 - I used the exact code but solved directly instead, 
> giving me the same results as in the tutorial (for the errors etc).
> 
> Even in step-20, the boundary values aren't correct for most of the tests, or 
> rather, have a lot of error in itself.
> For example, at the point (1,1) which is on the boundary, p = -1.1 with the 
> given test problem.
> 
> These are the values I obtained having run the code:
> at lowest degree (0), refinement level (RL) 3, p=-0.941992, which is rather 
> far off the -1.1 value for the Dirichlet condition applied. even at RL6, 
> p=1.07976
> I tried the next degree up (1), at RL3, p= -1.09984. eventually at RL6 for 
> this degree, we get -1.1.
> Similarly, at degree 2, at RL3. p is still not accurate at p = -1.10004

Individual values are usually not particularly meaningful. What matters is the 
convergence towards the correct value. Can you produce a table for each order 
in which you show the value as a function of refinement level?

I will add that it's not clear to me that the finite element method applied to 
this problem guarantees that the error *at a single point* converges to zero. 
We know that the error converges to zero when measured in integral quantities 
such as the L2 norm, but that mean not imply pointwise convergence.


> How else can we imposed the Dirichlet condition (without having to use very 
> high refinement levels and degrees) on the boundary for this problem?

I don't know of other ways. The question I would ask first is whether what you 
have *converges*. That's all you can really hope for. The second question is 
whether the error is *small*. If the method converges, then the error *will* 
be small if only the mesh is fine enough; but "fine enough" may be smaller 
than you can computationally afford, and in that case one can start to ask 
about alternatives (such as using the Q2-Q1 element used in step-43 instead of 
the RT element you have here), but the first step is to unambiguously 
establish convergence.

Best
  WB

-- 
------------------------------------------------------------------------
Wolfgang Bangerth          email:                 bange...@colostate.edu
                            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 dealii+unsubscr...@googlegroups.com.
For more options, visit https://groups.google.com/d/optout.

Reply via email to