Jane,
> Always setting local_rhs = 0 immediately before the below implementation > would take into account all cases so that would be the best: > > 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/current/doxygen/deal.II/classFEValuesBase.html#a130eea0fa89263d93b20521addc830c7>(q) > > * > boundary_values[q] * > fe_face_values.JxW > <https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#ad097580a2f71878695096cc73b271b9d>(q)); But it's already there, in the marked line (lkine 494 in the current sources): for (const auto &cell : dof_handler.active_cell_iterators()) { fe_values.reinit(cell); local_matrix = 0; local_rhs = 0; // *********************** right_hand_side.value_list(fe_values.get_quadrature_points(), rhs_values); k_inverse.value_list(fe_values.get_quadrature_points(), k_inverse_values); for (unsigned int q = 0; q < n_q_points; ++q) for (unsigned int i = 0; i < dofs_per_cell; ++i) ...cell integration... for (unsigned int face_n = 0; face_n < GeometryInfo<dim>::faces_per_cell; ++face_n) if (cell->at_boundary(face_n)) { fe_face_values.reinit(cell, face_n); pressure_boundary_values.value_list( fe_face_values.get_quadrature_points(), 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) += ... So all that is happening is that we accumulate contributions from all faces of the cell before we put it all into the global matrix. This is the right thing to do. > As for the Stokes case, I am making the analogy for the INHOMOGENEOUS > normal component of normal stress boundary condition. \ > You are right in that the no tangential stress condition requires > nothing, but then i have a nonzero value of the normal component of the > normal stress which needs to be in the weak form. > I can't quite figure out if this is a "dirichlet' condition so I also > set local_rhs = 0 before the boundary implementation, or whether it is > just an extra addition for the cells on that boundary. It is correct that the normal component of the normal stress is nonzero. But it is multiplied by a test function whose normal component is zero (because the normal component of the solution is prescribed). So the product of these two terms is zero, and there is no additional contribution to the weak form. Best W. -- ------------------------------------------------------------------------ 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.