Newton's method: compute residual in Step 15 and 57

52 views
Skip to first unread message

Lixing Zhu

unread,
Apr 10, 2022, 2:00:14 PM4/10/22
to deal.II User Group
Hi all,

I am confused with Newton's method's residual computation (i.e., L2 norm of the system residual vector) in Step-15 and Step-57.

In Step-15, to eliminate the residual component corresponding to the Dirichlet boundary condition, DOFTools::extract_boundary_dofs() is used to explicitly make those components zero.

However, in Step-57, there seems no explicit treatment to eliminate the residual components on the Dirichlet boundary in the computation of residuals in each iteration within Newton's method. Is this because of the zero_constraints?

Regards,
Lixing

SebG

unread,
Apr 11, 2022, 3:37:14 AM4/11/22
to deal.II User Group
Dear Lixing,

I think you are right. In Newton's method, the entries of the residual vector which correspond to degrees of freedom that are constrained through a Dirichlet boundary condition should vanish. In step-57, this is achieved using the member variable AffineConstraints<double> zero_constraints and the constraints are enforced during the assembly process by calling AffineConstraints<double>::distribute_local_to_global. Since inhomogeneous boundary conditions are applied in step-57, a switch between the inhomogeneous and homogeneous constraints is necessary. The inhomogeneous constraints are actually only used in the very first assembly of the linear system (system matrix and right-hand side vector). In all subsequent assemblies of the residual vector and the linear system, the homogeneous constraints are used. This is achieved by

const AffineConstraints<double> &constraints_used =
    initial_step ? nonzero_constraints : zero_constraints;

In step-15, Dirichlet boundary conditions and hanging node constraints are not treated using a joint AffineConstraints<double> object. The hanging node constraints are enforced through an explicit condensation after assembling the system matrix and right-hand side vector by means of an AffineConstraints<double> object. The Dirichlet boundary conditions are then enforced by calling  MatrixTools::apply_boundary_values and by explicitly setting them to zero in the residual vector. In step-15, only homogeneous constraints are applied because the solution automatically satisfies them due to the call to the method set_boundary_values.

Both approaches are equivalent, but should not be confused or mixed. I prefer the way of step-57, but this is a personal flavor. I hope my comments are helpful.

Best wishes,
Sebastian

Lixing Zhu

unread,
Apr 11, 2022, 11:40:30 AM4/11/22
to deal.II User Group
Dear Sebastian,

Many thanks for your detailed explanation.

Best,
Lixing
Reply all
Reply to author
Forward
0 new messages