From: David Wells Date: Wed, 15 Apr 2020 19:13:02 +0000 (-0400) Subject: step-8: Use AffineConstraints::distribute_local_to_global(). X-Git-Tag: v9.2.0-rc1~211^2~6 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=fd5b35f7d77a3671d78b49912a3315f41df9b540;p=dealii.git step-8: Use AffineConstraints::distribute_local_to_global(). --- diff --git a/examples/step-8/step-8.cc b/examples/step-8/step-8.cc index 32dabe4e87..4b658b76c0 100644 --- a/examples/step-8/step-8.cc +++ b/examples/step-8/step-8.cc @@ -98,7 +98,7 @@ namespace Step8 FESystem fe; - AffineConstraints hanging_node_constraints; + AffineConstraints constraints; SparsityPattern sparsity_pattern; SparseMatrix system_matrix; @@ -230,16 +230,19 @@ namespace Step8 void ElasticProblem::setup_system() { dof_handler.distribute_dofs(fe); - hanging_node_constraints.clear(); - DoFTools::make_hanging_node_constraints(dof_handler, - hanging_node_constraints); - hanging_node_constraints.close(); + constraints.clear(); + DoFTools::make_hanging_node_constraints(dof_handler, constraints); + VectorTools::interpolate_boundary_values(dof_handler, + 0, + Functions::ZeroFunction(dim), + constraints); + constraints.close(); DynamicSparsityPattern dsp(dof_handler.n_dofs(), dof_handler.n_dofs()); DoFTools::make_sparsity_pattern(dof_handler, dsp, - hanging_node_constraints, - /*keep_constrained_dofs = */ true); + constraints, + /*keep_constrained_dofs = */ false); sparsity_pattern.copy_from(dsp); system_matrix.reinit(sparsity_pattern); @@ -409,43 +412,11 @@ namespace Step8 // The transfer from local degrees of freedom into the global matrix // and right hand side vector does not depend on the equation under // consideration, and is thus the same as in all previous - // examples. The same holds for the elimination of hanging nodes from - // the matrix and right hand side, once we are done with assembling - // the entire linear system: + // examples. cell->get_dof_indices(local_dof_indices); - for (unsigned int i = 0; i < dofs_per_cell; ++i) - { - for (unsigned int j = 0; j < dofs_per_cell; ++j) - system_matrix.add(local_dof_indices[i], - local_dof_indices[j], - cell_matrix(i, j)); - - system_rhs(local_dof_indices[i]) += cell_rhs(i); - } + constraints.distribute_local_to_global( + cell_matrix, cell_rhs, local_dof_indices, system_matrix, system_rhs); } - - hanging_node_constraints.condense(system_matrix); - hanging_node_constraints.condense(system_rhs); - - // The interpolation of the boundary values needs a small modification: - // since the solution function is vector-valued, so need to be the - // boundary values. The Functions::ZeroFunction constructor - // accepts a parameter that tells it that it shall represent a vector - // valued, constant zero function with that many components. By default, - // this parameter is equal to one, in which case the - // Functions::ZeroFunction object would represent a scalar - // function. Since the solution vector has dim components, we - // need to pass dim as number of components to the zero - // function as well. - std::map boundary_values; - VectorTools::interpolate_boundary_values(dof_handler, - 0, - Functions::ZeroFunction(dim), - boundary_values); - MatrixTools::apply_boundary_values(boundary_values, - system_matrix, - solution, - system_rhs); } @@ -467,7 +438,7 @@ namespace Step8 cg.solve(system_matrix, solution, system_rhs, preconditioner); - hanging_node_constraints.distribute(solution); + constraints.distribute(solution); }