From: Wolfgang Bangerth Date: Fri, 8 Mar 2024 06:20:53 +0000 (-0700) Subject: Rename a few variables. X-Git-Tag: v9.6.0-rc1~495^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=44fbbd9bbb2219a026cbbecc9c933b2d3b6786ea;p=dealii.git Rename a few variables. Just make things look more like we have in most places. --- diff --git a/include/deal.II/numerics/vector_tools_interpolate.templates.h b/include/deal.II/numerics/vector_tools_interpolate.templates.h index 63c2285941..d0b51fc12a 100644 --- a/include/deal.II/numerics/vector_tools_interpolate.templates.h +++ b/include/deal.II/numerics/vector_tools_interpolate.templates.h @@ -862,17 +862,18 @@ namespace VectorTools template DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type) - void interpolate_to_different_mesh(const DoFHandler &dof1, - const VectorType &u1, - const DoFHandler &dof2, - VectorType &u2) + void interpolate_to_different_mesh( + const DoFHandler &dof_handler_1, + const VectorType &u1, + const DoFHandler &dof_handler_2, + VectorType &u2) { - Assert(GridTools::have_same_coarse_mesh(dof1, dof2), + Assert(GridTools::have_same_coarse_mesh(dof_handler_1, dof_handler_2), ExcMessage("The two DoF handlers must represent triangulations that " "have the same coarse meshes")); InterGridMap> intergridmap; - intergridmap.make_mapping(dof1, dof2); + intergridmap.make_mapping(dof_handler_1, dof_handler_2); AffineConstraints dummy; dummy.close(); @@ -885,18 +886,18 @@ namespace VectorTools template DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type) void interpolate_to_different_mesh( - const DoFHandler &dof1, + const DoFHandler &dof_handler_1, const VectorType &u1, - const DoFHandler &dof2, + const DoFHandler &dof_handler_2, const AffineConstraints &constraints, VectorType &u2) { - Assert(GridTools::have_same_coarse_mesh(dof1, dof2), + Assert(GridTools::have_same_coarse_mesh(dof_handler_1, dof_handler_2), ExcMessage("The two DoF handlers must represent triangulations that " "have the same coarse meshes")); InterGridMap> intergridmap; - intergridmap.make_mapping(dof1, dof2); + intergridmap.make_mapping(dof_handler_1, dof_handler_2); interpolate_to_different_mesh(intergridmap, u1, constraints, u2); } @@ -911,15 +912,17 @@ namespace VectorTools const AffineConstraints &constraints, VectorType &u2) { - const DoFHandler &dof1 = intergridmap.get_source_grid(); - const DoFHandler &dof2 = intergridmap.get_destination_grid(); - (void)dof2; + const DoFHandler &dof_handler_1 = + intergridmap.get_source_grid(); + const DoFHandler &dof_handler_2 = + intergridmap.get_destination_grid(); + (void)dof_handler_2; - Assert(dof1.get_fe_collection() == dof2.get_fe_collection(), - ExcMessage( - "The FECollections of both DoFHandler objects must match")); - AssertDimension(u1.size(), dof1.n_dofs()); - AssertDimension(u2.size(), dof2.n_dofs()); + Assert( + dof_handler_1.get_fe_collection() == dof_handler_2.get_fe_collection(), + ExcMessage("The FECollections of both DoFHandler objects must match")); + AssertDimension(u1.size(), dof_handler_1.n_dofs()); + AssertDimension(u2.size(), dof_handler_2.n_dofs()); Vector cache; @@ -932,10 +935,7 @@ namespace VectorTools // Therefore, loop over all cells // (active and inactive) of the source // grid .. - typename DoFHandler::cell_iterator cell1 = dof1.begin(); - const typename DoFHandler::cell_iterator endc1 = dof1.end(); - - for (; cell1 != endc1; ++cell1) + for (const auto &cell1 : dof_handler_1.cell_iterators()) { const typename DoFHandler::cell_iterator cell2 = intergridmap[cell1];