From: Wolfgang Bangerth Date: Thu, 22 Jun 2023 17:05:52 +0000 (-0600) Subject: Simplify some code that uses SolutionTransfer. X-Git-Tag: relicensing~849^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=761c11eaa6d65d553f4fad69961a00f922553a3c;p=dealii.git Simplify some code that uses SolutionTransfer. --- diff --git a/examples/step-31/step-31.cc b/examples/step-31/step-31.cc index a9cc1faa1b..7222de350a 100644 --- a/examples/step-31/step-31.cc +++ b/examples/step-31/step-31.cc @@ -1943,9 +1943,8 @@ namespace Step31 // and temperature DoFHandler objects, by attaching them to the old dof // handlers. With this at place, we can prepare the triangulation and the // data vectors for refinement (in this order). - std::vector x_temperature(2); - x_temperature[0] = temperature_solution; - x_temperature[1] = old_temperature_solution; + const std::vector x_temperature = { + temperature_solution, old_temperature_solution}; TrilinosWrappers::MPI::BlockVector x_stokes = stokes_solution; SolutionTransfer temperature_trans( @@ -1971,13 +1970,13 @@ namespace Step31 triangulation.execute_coarsening_and_refinement(); setup_dofs(); - std::vector tmp(2); - tmp[0].reinit(temperature_solution); - tmp[1].reinit(temperature_solution); + std::vector tmp = { + TrilinosWrappers::MPI::Vector(temperature_solution), + TrilinosWrappers::MPI::Vector(temperature_solution)}; temperature_trans.interpolate(x_temperature, tmp); - temperature_solution = tmp[0]; - old_temperature_solution = tmp[1]; + temperature_solution = std::move(tmp[0]); + old_temperature_solution = std::move(tmp[1]); // After the solution has been transferred we then enforce the constraints // on the transferred solution. diff --git a/examples/step-32/step-32.cc b/examples/step-32/step-32.cc index 575f6d4a22..e6f3c5af77 100644 --- a/examples/step-32/step-32.cc +++ b/examples/step-32/step-32.cc @@ -3294,12 +3294,10 @@ namespace Step32 // remainder of the function further down below is then concerned with // setting up the data structures again after mesh refinement and // restoring the solution vectors on the new mesh. - std::vector x_temperature(2); - x_temperature[0] = &temperature_solution; - x_temperature[1] = &old_temperature_solution; - std::vector x_stokes(2); - x_stokes[0] = &stokes_solution; - x_stokes[1] = &old_stokes_solution; + const std::vector x_temperature = { + &temperature_solution, &old_temperature_solution}; + const std::vector x_stokes = { + &stokes_solution, &old_stokes_solution}; triangulation.prepare_coarsening_and_refinement(); @@ -3319,9 +3317,8 @@ namespace Step32 TrilinosWrappers::MPI::Vector distributed_temp1(temperature_rhs); TrilinosWrappers::MPI::Vector distributed_temp2(temperature_rhs); - std::vector tmp(2); - tmp[0] = &(distributed_temp1); - tmp[1] = &(distributed_temp2); + std::vector tmp = {&distributed_temp1, + &distributed_temp2}; temperature_trans.interpolate(tmp); // enforce constraints to make the interpolated solution conforming on @@ -3329,17 +3326,16 @@ namespace Step32 temperature_constraints.distribute(distributed_temp1); temperature_constraints.distribute(distributed_temp2); - temperature_solution = distributed_temp1; - old_temperature_solution = distributed_temp2; + temperature_solution = std::move(distributed_temp1); + old_temperature_solution = std::move(distributed_temp2); } { TrilinosWrappers::MPI::BlockVector distributed_stokes(stokes_rhs); TrilinosWrappers::MPI::BlockVector old_distributed_stokes(stokes_rhs); - std::vector stokes_tmp(2); - stokes_tmp[0] = &(distributed_stokes); - stokes_tmp[1] = &(old_distributed_stokes); + std::vector stokes_tmp = { + &distributed_stokes, &old_distributed_stokes}; stokes_trans.interpolate(stokes_tmp); @@ -3348,8 +3344,8 @@ namespace Step32 stokes_constraints.distribute(distributed_stokes); stokes_constraints.distribute(old_distributed_stokes); - stokes_solution = distributed_stokes; - old_stokes_solution = old_distributed_stokes; + stokes_solution = std::move(distributed_stokes); + old_stokes_solution = std::move(old_distributed_stokes); } } } diff --git a/examples/step-33/step-33.cc b/examples/step-33/step-33.cc index 0e5b75d23f..b54c95f895 100644 --- a/examples/step-33/step-33.cc +++ b/examples/step-33/step-33.cc @@ -2264,11 +2264,7 @@ namespace Step33 // examples, so we won't comment much on the following code. The last // three lines simply re-set the sizes of some other vectors to the now // correct size: - std::vector> transfer_in; - std::vector> transfer_out; - - transfer_in.push_back(old_solution); - transfer_in.push_back(predictor); + const std::vector> transfer_in = {old_solution, predictor}; triangulation.prepare_coarsening_and_refinement(); @@ -2280,26 +2276,17 @@ namespace Step33 dof_handler.clear(); dof_handler.distribute_dofs(fe); - { - Vector new_old_solution(1); - Vector new_predictor(1); - - transfer_out.push_back(new_old_solution); - transfer_out.push_back(new_predictor); - transfer_out[0].reinit(dof_handler.n_dofs()); - transfer_out[1].reinit(dof_handler.n_dofs()); - } - + std::vector> transfer_out = { + Vector(dof_handler.n_dofs()), + Vector(dof_handler.n_dofs())}; soltrans.interpolate(transfer_in, transfer_out); - old_solution.reinit(transfer_out[0].size()); - old_solution = transfer_out[0]; - - predictor.reinit(transfer_out[1].size()); - predictor = transfer_out[1]; + old_solution = std::move(transfer_out[0]); + predictor = std::move(transfer_out[1]); current_solution.reinit(dof_handler.n_dofs()); current_solution = old_solution; + right_hand_side.reinit(dof_handler.n_dofs()); } diff --git a/include/deal.II/numerics/solution_transfer.h b/include/deal.II/numerics/solution_transfer.h index 4edd1632bf..b527113c9e 100644 --- a/include/deal.II/numerics/solution_transfer.h +++ b/include/deal.II/numerics/solution_transfer.h @@ -91,7 +91,8 @@ DEAL_II_NAMESPACE_OPEN * std::vector > solutions(n_vectors, Vector (n)); * soltrans.refine_interpolate(solutions_old, solutions); * @endcode - * This is used in several of the tutorial programs, for example step-31. + * This is used in several of the tutorial programs, for example step-31 + * and step-33. * *
  • If the grid has cells that will be coarsened, then use @p * SolutionTransfer as follows: