Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Simplify some code that uses SolutionTransfer. #15438

Merged
merged 1 commit into from
Jun 29, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
15 changes: 7 additions & 8 deletions examples/step-31/step-31.cc
Original file line number Diff line number Diff line change
Expand Up @@ -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<TrilinosWrappers::MPI::Vector> x_temperature(2);
x_temperature[0] = temperature_solution;
x_temperature[1] = old_temperature_solution;
const std::vector<TrilinosWrappers::MPI::Vector> x_temperature = {
temperature_solution, old_temperature_solution};
TrilinosWrappers::MPI::BlockVector x_stokes = stokes_solution;

SolutionTransfer<dim, TrilinosWrappers::MPI::Vector> temperature_trans(
Expand All @@ -1971,13 +1970,13 @@ namespace Step31
triangulation.execute_coarsening_and_refinement();
setup_dofs();

std::vector<TrilinosWrappers::MPI::Vector> tmp(2);
tmp[0].reinit(temperature_solution);
tmp[1].reinit(temperature_solution);
std::vector<TrilinosWrappers::MPI::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.
Expand Down
28 changes: 12 additions & 16 deletions examples/step-32/step-32.cc
Original file line number Diff line number Diff line change
Expand Up @@ -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<const TrilinosWrappers::MPI::Vector *> x_temperature(2);
x_temperature[0] = &temperature_solution;
x_temperature[1] = &old_temperature_solution;
std::vector<const TrilinosWrappers::MPI::BlockVector *> x_stokes(2);
x_stokes[0] = &stokes_solution;
x_stokes[1] = &old_stokes_solution;
const std::vector<const TrilinosWrappers::MPI::Vector *> x_temperature = {
&temperature_solution, &old_temperature_solution};
const std::vector<const TrilinosWrappers::MPI::BlockVector *> x_stokes = {
&stokes_solution, &old_stokes_solution};

triangulation.prepare_coarsening_and_refinement();

Expand All @@ -3319,27 +3317,25 @@ namespace Step32
TrilinosWrappers::MPI::Vector distributed_temp1(temperature_rhs);
TrilinosWrappers::MPI::Vector distributed_temp2(temperature_rhs);

std::vector<TrilinosWrappers::MPI::Vector *> tmp(2);
tmp[0] = &(distributed_temp1);
tmp[1] = &(distributed_temp2);
std::vector<TrilinosWrappers::MPI::Vector *> tmp = {&distributed_temp1,
&distributed_temp2};
temperature_trans.interpolate(tmp);

// enforce constraints to make the interpolated solution conforming on
// the new mesh:
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<TrilinosWrappers::MPI::BlockVector *> stokes_tmp(2);
stokes_tmp[0] = &(distributed_stokes);
stokes_tmp[1] = &(old_distributed_stokes);
std::vector<TrilinosWrappers::MPI::BlockVector *> stokes_tmp = {
&distributed_stokes, &old_distributed_stokes};

stokes_trans.interpolate(stokes_tmp);

Expand All @@ -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);
}
}
}
Expand Down
27 changes: 7 additions & 20 deletions examples/step-33/step-33.cc
Original file line number Diff line number Diff line change
Expand Up @@ -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<Vector<double>> transfer_in;
std::vector<Vector<double>> transfer_out;

transfer_in.push_back(old_solution);
transfer_in.push_back(predictor);
const std::vector<Vector<double>> transfer_in = {old_solution, predictor};

triangulation.prepare_coarsening_and_refinement();

Expand All @@ -2280,26 +2276,17 @@ namespace Step33
dof_handler.clear();
dof_handler.distribute_dofs(fe);

{
Vector<double> new_old_solution(1);
Vector<double> 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<Vector<double>> transfer_out = {
Vector<double>(dof_handler.n_dofs()),
Vector<double>(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());
}

Expand Down
3 changes: 2 additions & 1 deletion include/deal.II/numerics/solution_transfer.h
Original file line number Diff line number Diff line change
Expand Up @@ -91,7 +91,8 @@ DEAL_II_NAMESPACE_OPEN
* std::vector<Vector<double> > solutions(n_vectors, Vector<double> (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.
*
* <li> If the grid has cells that will be coarsened, then use @p
* SolutionTransfer as follows:
Expand Down