On 4/4/19 4:05 AM, Muhammad Mashhood wrote:
> 
>       for (unsigned int f=0; f<GeometryInfo<dim> :: faces_per_cell; ++f)
>              if (cell->face(f)->at_boundary() )
>              {
>                  fe_face_values.reinit (cell,f);
>                  for (unsigned int i=0; i<dofs_per_cell; ++i)
>                    {
>                      const unsigned int
>                      component_i = fe.system_to_component_index(i).first;
> 
>                      for (unsigned int q_point=0; q_point<n_face_q_points; 
> ++q_point)
>                        cell_rhs(i) += fe_face_values.shape_value(i,q_point) *
>                                       rhs_values[q_point][component_i] *
>                                       fe_face_values.JxW(q_point);
>                                          }
>              }
> 
> 
>   Where the *"rhs_values"* are being extracted from my previously mentioned 
> function "*right_hand_side (fe_values.get_quadrature_points(), rhs_values, 
> cell);"* in a same loop of "*for (; cell!=endc; ++cell)" *and the code for 
> this "*right_hand_side" *function is that which I mentioned previously.
> In my limited knowledge I am going through each face in cell as in 2D there 
> are four edges are present and I want to choose the edge number 3 (which is 
> upper outer edge in my example case).

The way you do things here does not work. You are assembling boundary terms on 
an individual face here. You need to evaluate the rhs_values *for this 
particular face*. In other words, you need to call the function you showed us 
with
   fe_face_values.get_quadrature_points()
(not fe_values.get_quadrature_points()), and not with an iterator to the cell 
but with an iterator to the face.

This way, you do not need the function you showed us guess whether a 
quadrature point is near a particular boundary -- the function already *knows* 
that the quadrature points lie on that boundary, because it is only called for 
such points.

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