Dear Nan,

basically, you can see this in both ways. Only two of the four shape functions in the cell have support on one particular face. Thus, you end up with summing up the contributions of two shape functions and two times zero.

However, FEFaceValues::shape_value expects the index of the cell shape functions and not of the shape functions on the face. Thus, this code is correct, indeed.

Best Regards,
Markus



On 10.08.2012 16:58, Zhang, Nan wrote:
Hi there,

I am studying how to code up the Neumann boundary condition with step-7. I am confused by following segment of code in step-7
**********************************************************************************************************
for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
        if (cell->face(face)->at_boundary()
            &&
            (cell->face(face)->boundary_indicator() == 1))
          {
                fe_face_values.reinit (cell, face);
for (unsigned int q_point=0; q_point<n_face_q_points; ++q_point)
              {
                const double neumann_value
= (exact_solution.gradient (fe_face_values.quadrature_point(q_point)) *
                     fe_face_values.normal_vector(q_point));

                for (unsigned int i=0; i<dofs_per_cell; ++i)
                  cell_rhs(i) += (neumann_value *
fe_face_values.shape_value(i,q_point) *
fe_face_values.JxW(q_point));
              }
          }
*******************************************************************************************************
If I understand it correctly, for the 2D problem, a particular cell, you detect one of four faces of this cell with indicator 1. Then you do the integration on this face. You need two values on two n_face_q_points (bilinear shape function) using fe_face_values. But what is fe_face_values.shape_value(i,q_point) where i=1->4? Because you are using shape function on a face, the i in the fe_face_values.shape_value(i,q_point) is the dofs_per_face and can only be 1 and 2. However, you use dofs_per_cell that is 1 to 4. I am really confused. Is the shape function for the Neumann boundary, in the weak form f(K,i), a shape function on face or on cell?

--
Dr. Nan Zhang
Postdoctoral Research Associate
Department of Geological Science
Brown University
324 Brook Street
Box 1846
Providence, RI 02912
HomePage: http://www.geo.brown.edu/People/Postdocs/Zhang/NanZhangHome.html


_______________________________________________
dealii mailing list http://poisson.dealii.org/mailman/listinfo/dealii

Reply via email to