22#ifndef EXADG_GRID_MAPPING_DEFORMATION_STRUCTURE_H_
23#define EXADG_GRID_MAPPING_DEFORMATION_STRUCTURE_H_
26#include <deal.II/base/timer.h>
29#include <exadg/grid/mapping_deformation_base.h>
30#include <exadg/structure/spatial_discretization/operator.h>
31#include <exadg/utilities/print_solver_results.h>
42template<
int dim,
typename Number>
46 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
57 std::shared_ptr<MaterialDescriptor const> material_descriptor,
59 std::string
const & field,
60 bool const setup_scalar_field,
61 MPI_Comm
const & mpi_comm)
64 pcout(std::cout, dealii::
Utilities::MPI::this_mpi_process(mpi_comm) == 0),
65 iterations({0, {0, 0}})
68 pde_operator = std::make_shared<Operator<dim, Number>>(grid,
70 multigrid_mappings_undeformed,
80 pde_operator->setup();
83 pde_operator->initialize_dof_vector(displacement);
86 std::shared_ptr<Operator<dim, Number>
const>
87 get_pde_operator()
const
92 dealii::MatrixFree<dim, Number>
const &
93 get_matrix_free()
const
95 return *pde_operator->get_matrix_free();
103 bool const print_solver_info,
104 types::time_step time_step_number)
override
109 if(param.large_deformation)
111 VectorType const_vector;
113 bool const update_preconditioner =
114 this->param.update_preconditioner &&
115 time_step_number % this->param.update_preconditioner_every_time_steps == 0;
117 auto const iter = pde_operator->solve_nonlinear(displacement,
122 update_preconditioner);
124 iterations.first += 1;
125 std::get<0>(iterations.second) += std::get<0>(iter);
126 std::get<1>(iterations.second) += std::get<1>(iter);
128 if(print_solver_info)
130 this->pcout << std::endl <<
"Solve moving mesh problem (nonlinear elasticity):";
131 print_solver_info_nonlinear(pcout, std::get<0>(iter), std::get<1>(iter), timer.wall_time());
138 pde_operator->initialize_dof_vector(rhs);
139 pde_operator->rhs(rhs, time);
141 auto const iter = pde_operator->solve_linear(displacement,
148 iterations.first += 1;
149 std::get<1>(iterations.second) += iter;
151 if(print_solver_info)
153 this->pcout << std::endl <<
"Solve moving mesh problem (linear elasticity):";
154 print_solver_info_linear(pcout, iter, timer.wall_time());
160 pde_operator->get_dof_handler());
169 std::vector<std::string> names;
170 std::vector<double> iterations_avg;
172 if(param.large_deformation)
174 names = {
"Nonlinear iterations",
175 "Linear iterations (accumulated)",
176 "Linear iterations (per nonlinear it.)"};
178 iterations_avg.resize(3);
180 (double)std::get<0>(iterations.second) / std::max(1., (
double)iterations.first);
182 (double)std::get<1>(iterations.second) / std::max(1., (
double)iterations.first);
183 if(iterations_avg[0] > std::numeric_limits<double>::min())
184 iterations_avg[2] = iterations_avg[1] / iterations_avg[0];
186 iterations_avg[2] = iterations_avg[1];
190 names = {
"Linear iterations"};
191 iterations_avg.resize(1);
193 (double)std::get<1>(iterations.second) / std::max(1., (
double)iterations.first);
196 print_list_of_iterations(pcout, names, iterations_avg);
201 std::shared_ptr<MatrixFreeData<dim, Number>> matrix_free_data;
202 std::shared_ptr<dealii::MatrixFree<dim, Number>> matrix_free;
205 std::shared_ptr<Operator<dim, Number>> pde_operator;
214 dealii::ConditionalOStream pcout;
218 std::tuple<unsigned long long, unsigned long long> >
void initialize_mapping_from_dof_vector(std::shared_ptr< dealii::Mapping< dim > const > mapping, VectorType const &displacement_vector, dealii::DoFHandler< dim > const &dof_handler)
Definition mapping_dof_vector.h:212
Definition parameters.h:41
Definition interpolate.h:34
Definition boundary_descriptor.h:48
Definition field_functions.h:32