23#ifndef INCLUDE_EXADG_FLUID_STRUCTURE_INTERACTION_ACCELERATION_SCHEMES_PARTITIONED_SOLVER_H_
24#define INCLUDE_EXADG_FLUID_STRUCTURE_INTERACTION_ACCELERATION_SCHEMES_PARTITIONED_SOLVER_H_
27#include <exadg/fluid_structure_interaction/acceleration_schemes/linear_algebra.h>
28#include <exadg/fluid_structure_interaction/acceleration_schemes/parameters.h>
29#include <exadg/fluid_structure_interaction/single_field_solvers/fluid.h>
30#include <exadg/fluid_structure_interaction/single_field_solvers/structure.h>
33#include <exadg/utilities/print_solver_results.h>
34#include <exadg/utilities/timer_tree.h>
40template<
int dim,
typename Number>
41class PartitionedSolver
44 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
47 PartitionedSolver(
Parameters const & parameters, MPI_Comm
const & comm);
54 solve(std::function<
void(VectorType &, VectorType
const &,
unsigned int)>
const &
55 apply_dirichlet_neumann_scheme);
58 print_iterations(dealii::ConditionalOStream
const & pcout)
const;
60 std::shared_ptr<TimerTree>
64 get_structure_velocity(VectorType & velocity_structure,
unsigned int const iteration)
const;
68 get_structure_displacement(VectorType & displacement_structure,
69 unsigned int const iteration)
const;
72 check_convergence(VectorType
const & residual)
const;
75 print_solver_info_header(
unsigned int const iteration)
const;
78 print_solver_info_converged(
unsigned int const iteration)
const;
83 dealii::ConditionalOStream pcout;
85 std::shared_ptr<SolverFluid<dim, Number>> fluid;
86 std::shared_ptr<SolverStructure<dim, Number>> structure;
89 std::vector<std::shared_ptr<std::vector<VectorType>>> D_history, R_history, Z_history;
92 std::shared_ptr<TimerTree> timer_tree;
98 std::pair<unsigned int, unsigned long long> partitioned_iterations;
101template<
int dim,
typename Number>
102PartitionedSolver<dim, Number>::PartitionedSolver(
Parameters const & parameters,
103 MPI_Comm
const & comm)
104 : parameters(parameters),
105 pcout(std::cout, dealii::
Utilities::MPI::this_mpi_process(comm) == 0),
106 partitioned_iterations({0, 0})
108 timer_tree = std::make_shared<TimerTree>();
111template<
int dim,
typename Number>
113PartitionedSolver<dim, Number>::setup(std::shared_ptr<SolverFluid<dim, Number>> fluid_,
114 std::shared_ptr<SolverStructure<dim, Number>> structure_)
117 structure = structure_;
120template<
int dim,
typename Number>
122PartitionedSolver<dim, Number>::check_convergence(VectorType
const & residual)
const
124 double const residual_norm = residual.l2_norm();
125 double const ref_norm_abs = std::sqrt(structure->pde_operator->get_number_of_dofs());
126 double const ref_norm_rel = structure->time_integrator->get_velocity_np().l2_norm() *
127 structure->time_integrator->get_time_step_size();
129 bool const converged = (residual_norm < parameters.abs_tol * ref_norm_abs) or
130 (residual_norm < parameters.rel_tol * ref_norm_rel);
135template<
int dim,
typename Number>
137PartitionedSolver<dim, Number>::print_solver_info_header(
unsigned int const iteration)
const
139 if(fluid->time_integrator->print_solver_info())
142 <<
"======================================================================" << std::endl
143 <<
" Partitioned FSI: iteration counter = " << std::left << std::setw(8) << iteration
145 <<
"======================================================================" << std::endl;
149template<
int dim,
typename Number>
151PartitionedSolver<dim, Number>::print_solver_info_converged(
unsigned int const iteration)
const
153 if(fluid->time_integrator->print_solver_info())
156 <<
"Partitioned FSI iteration converged in " << iteration <<
" iterations." << std::endl;
160template<
int dim,
typename Number>
162PartitionedSolver<dim, Number>::print_iterations(dealii::ConditionalOStream
const & pcout)
const
164 std::vector<std::string> names;
165 std::vector<double> iterations_avg;
167 names = {
"Partitioned iterations"};
168 iterations_avg.resize(1);
170 (double)partitioned_iterations.second / std::max(1.0, (
double)partitioned_iterations.first);
172 print_list_of_iterations(pcout, names, iterations_avg);
175template<
int dim,
typename Number>
176std::shared_ptr<TimerTree>
177PartitionedSolver<dim, Number>::get_timings()
const
182template<
int dim,
typename Number>
184PartitionedSolver<dim, Number>::get_structure_velocity(VectorType & velocity_structure,
185 unsigned int iteration)
const
189 if(parameters.initial_guess_coupling_scheme ==
190 InitialGuessCouplingScheme::SolutionExtrapolatedToEndOfTimeStep)
192 structure->time_integrator->extrapolate_velocity_to_np(velocity_structure);
194 else if(parameters.initial_guess_coupling_scheme ==
195 InitialGuessCouplingScheme::ConvergedSolutionOfPreviousTimeStep)
197 velocity_structure = structure->time_integrator->get_velocity_n();
203 "Behavior for this `InitialGuessCouplingScheme` is not defined."));
208 velocity_structure = structure->time_integrator->get_velocity_np();
212template<
int dim,
typename Number>
214PartitionedSolver<dim, Number>::get_structure_displacement(VectorType & displacement_structure,
215 unsigned int const iteration)
const
219 if(this->parameters.initial_guess_coupling_scheme ==
220 InitialGuessCouplingScheme::SolutionExtrapolatedToEndOfTimeStep)
222 structure->time_integrator->extrapolate_displacement_to_np(displacement_structure);
224 else if(this->parameters.initial_guess_coupling_scheme ==
225 InitialGuessCouplingScheme::ConvergedSolutionOfPreviousTimeStep)
227 displacement_structure = structure->time_integrator->get_displacement_n();
233 "Behavior for this `InitialGuessCouplingScheme` is not defined."));
238 displacement_structure = structure->time_integrator->get_displacement_np();
242template<
int dim,
typename Number>
244PartitionedSolver<dim, Number>::solve(
245 std::function<
void(VectorType &, VectorType
const &,
unsigned int)>
const &
246 apply_dirichlet_neumann_scheme)
252 if(parameters.acceleration_method == AccelerationMethod::Aitken)
255 structure->pde_operator->initialize_dof_vector(r_old);
256 structure->pde_operator->initialize_dof_vector(d);
258 bool converged =
false;
260 while(not(converged) and k < parameters.partitioned_iter_max)
262 print_solver_info_header(k);
264 get_structure_displacement(d, k);
266 VectorType d_tilde(d);
267 apply_dirichlet_neumann_scheme(d_tilde, d, k);
270 VectorType r = d_tilde;
272 converged = check_convergence(r);
282 omega = parameters.omega_init;
286 VectorType delta_r = r;
287 delta_r.add(-1.0, r_old);
288 omega *= -(r_old * delta_r) / delta_r.norm_sqr();
294 structure->time_integrator->set_displacement(d);
296 timer_tree->insert({
"Aitken"}, timer.wall_time());
303 else if(parameters.acceleration_method == AccelerationMethod::IQN_ILS)
305 std::shared_ptr<std::vector<VectorType>> D, R;
306 D = std::make_shared<std::vector<VectorType>>();
307 R = std::make_shared<std::vector<VectorType>>();
309 VectorType d, d_tilde, d_tilde_old, r, r_old;
310 structure->pde_operator->initialize_dof_vector(d);
311 structure->pde_operator->initialize_dof_vector(d_tilde);
312 structure->pde_operator->initialize_dof_vector(d_tilde_old);
313 structure->pde_operator->initialize_dof_vector(r);
314 structure->pde_operator->initialize_dof_vector(r_old);
316 unsigned int const q = parameters.reused_time_steps;
317 unsigned int const n = fluid->time_integrator->get_number_of_time_steps();
319 bool converged =
false;
320 while(not(converged) and k < parameters.partitioned_iter_max)
322 print_solver_info_header(k);
324 get_structure_displacement(d, k);
326 apply_dirichlet_neumann_scheme(d_tilde, d, k);
331 converged = check_convergence(r);
339 if(k == 0 and (q == 0 or n == 0))
341 d.add(parameters.omega_init, r);
348 VectorType delta_d_tilde = d_tilde;
349 delta_d_tilde.add(-1.0, d_tilde_old);
350 D->push_back(delta_d_tilde);
352 VectorType delta_r = r;
353 delta_r.add(-1.0, r_old);
354 R->push_back(delta_r);
358 std::vector<VectorType> Q = *R;
359 for(
auto R_q : R_history)
360 for(
auto delta_r : *R_q)
361 Q.push_back(delta_r);
362 std::vector<VectorType> D_all = *D;
363 for(
auto D_q : D_history)
364 for(
auto delta_d : *D_q)
365 D_all.push_back(delta_d);
367 AssertThrow(D_all.size() == Q.size(),
368 dealii::ExcMessage(
"D, Q vectors must have same size."));
370 unsigned int const k_all = Q.size();
375 compute_QR_decomposition(Q, U);
377 std::vector<Number> rhs(k_all, 0.0);
378 for(
unsigned int i = 0; i < k_all; ++i)
379 rhs[i] = -Number(Q[i] * r);
382 std::vector<Number> alpha(k_all, 0.0);
383 backward_substitution(U, alpha, rhs);
387 for(
unsigned int i = 0; i < k_all; ++i)
388 d.add(alpha[i], D_all[i]);
392 d.add(parameters.omega_init, r);
396 d_tilde_old = d_tilde;
399 structure->time_integrator->set_displacement(d);
401 timer_tree->insert({
"IQN-ILS"}, timer.wall_time());
412 D_history.push_back(D);
413 R_history.push_back(R);
414 if(D_history.size() > q)
415 D_history.erase(D_history.begin());
416 if(R_history.size() > q)
417 R_history.erase(R_history.begin());
419 timer_tree->insert({
"IQN-ILS"}, timer.wall_time());
421 else if(parameters.acceleration_method == AccelerationMethod::IQN_IMVLS)
423 std::shared_ptr<std::vector<VectorType>> D, R;
424 D = std::make_shared<std::vector<VectorType>>();
425 R = std::make_shared<std::vector<VectorType>>();
427 std::vector<VectorType> B;
429 VectorType d, d_tilde, d_tilde_old, r, r_old, b, b_old;
430 structure->pde_operator->initialize_dof_vector(d);
431 structure->pde_operator->initialize_dof_vector(d_tilde);
432 structure->pde_operator->initialize_dof_vector(d_tilde_old);
433 structure->pde_operator->initialize_dof_vector(r);
434 structure->pde_operator->initialize_dof_vector(r_old);
435 structure->pde_operator->initialize_dof_vector(b);
436 structure->pde_operator->initialize_dof_vector(b_old);
438 std::shared_ptr<Matrix<Number>> U;
439 std::vector<VectorType> Q;
441 unsigned int const q = parameters.reused_time_steps;
442 unsigned int const n = fluid->time_integrator->get_number_of_time_steps();
444 bool converged =
false;
445 while(not(converged) and k < parameters.partitioned_iter_max)
447 print_solver_info_header(k);
449 get_structure_displacement(d, k);
451 apply_dirichlet_neumann_scheme(d_tilde, d, k);
456 converged = check_convergence(r);
465 inv_jacobian_times_residual(b, D_history, R_history, Z_history, r);
467 if(k == 0 and (q == 0 or n == 0))
469 d.add(parameters.omega_init, r);
479 VectorType delta_d_tilde = d_tilde;
480 delta_d_tilde.add(-1.0, d_tilde_old);
481 D->push_back(delta_d_tilde);
483 VectorType delta_r = r;
484 delta_r.add(-1.0, r_old);
485 R->push_back(delta_r);
487 VectorType delta_b = delta_d_tilde;
488 delta_b.add(1.0, b_old);
489 delta_b.add(-1.0, b);
490 B.push_back(delta_b);
493 U = std::make_shared<Matrix<Number>>(k);
495 compute_QR_decomposition(Q, *U);
497 std::vector<Number> rhs(k, 0.0);
498 for(
unsigned int i = 0; i < k; ++i)
499 rhs[i] = -Number(Q[i] * r);
502 std::vector<Number> alpha(k, 0.0);
503 backward_substitution(*U, alpha, rhs);
505 for(
unsigned int i = 0; i < k; ++i)
506 d.add(alpha[i], B[i]);
510 d_tilde_old = d_tilde;
514 structure->time_integrator->set_displacement(d);
516 timer_tree->insert({
"IQN-IMVLS"}, timer.wall_time());
527 D_history.push_back(D);
528 R_history.push_back(R);
529 if(D_history.size() > q)
530 D_history.erase(D_history.begin());
531 if(R_history.size() > q)
532 R_history.erase(R_history.begin());
535 std::shared_ptr<std::vector<VectorType>> Z;
536 Z = std::make_shared<std::vector<VectorType>>();
538 backward_substitution_multiple_rhs(*U, *Z, Q);
539 Z_history.push_back(Z);
540 if(Z_history.size() > q)
541 Z_history.erase(Z_history.begin());
543 timer_tree->insert({
"IQN-IMVLS"}, timer.wall_time());
547 AssertThrow(
false, dealii::ExcMessage(
"This AccelerationMethod is not implemented."));
550 partitioned_iterations.first += 1;
551 partitioned_iterations.second += k;
553 print_solver_info_converged(k);
Definition linear_algebra.h:40
Definition structure.h:38
Definition interpolate.h:32
Definition parameters.h:47