ExaDG
Loading...
Searching...
No Matches
mapping_deformation_structure.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef INCLUDE_EXADG_GRID_MAPPING_DEFORMATION_STRUCTURE_H_
23#define INCLUDE_EXADG_GRID_MAPPING_DEFORMATION_STRUCTURE_H_
24
25// deal.II
26#include <deal.II/base/timer.h>
27
28// ExaDG
29#include <exadg/grid/mapping_deformation_base.h>
30#include <exadg/structure/spatial_discretization/operator.h>
31#include <exadg/utilities/print_solver_results.h>
32
33namespace ExaDG
34{
35namespace Structure
36{
42template<int dim, typename Number>
43class DeformedMapping : public DeformedMappingBase<dim, Number>
44{
45public:
46 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
47
52 std::shared_ptr<Grid<dim> const> grid,
53 std::shared_ptr<dealii::Mapping<dim> const> mapping_undeformed,
54 std::shared_ptr<MultigridMappings<dim, Number>> const multigrid_mappings_undeformed,
55 std::shared_ptr<BoundaryDescriptor<dim> const> boundary_descriptor,
56 std::shared_ptr<FieldFunctions<dim> const> field_functions,
57 std::shared_ptr<MaterialDescriptor const> material_descriptor,
58 Parameters const & param,
59 std::string const & field,
60 MPI_Comm const & mpi_comm)
61 : DeformedMappingBase<dim, Number>(mapping_undeformed, param.degree, *grid->triangulation),
62 param(param),
63 pcout(std::cout, dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0),
64 iterations({0, {0, 0}})
65 {
66 // initialize PDE operator
67 pde_operator = std::make_shared<Operator<dim, Number>>(grid,
69 multigrid_mappings_undeformed,
70 boundary_descriptor,
71 field_functions,
72 material_descriptor,
73 param,
74 field,
75 mpi_comm);
76
77 // setup PDE operator and solver
78 pde_operator->setup();
79
80 // finally, initialize dof vector
81 pde_operator->initialize_dof_vector(displacement);
82 }
83
84 std::shared_ptr<Operator<dim, Number> const>
85 get_pde_operator() const
86 {
87 return pde_operator;
88 }
89
90 dealii::MatrixFree<dim, Number> const &
91 get_matrix_free() const
92 {
93 return *pde_operator->get_matrix_free();
94 }
95
99 void
100 update(double const time,
101 bool const print_solver_info,
102 types::time_step time_step_number) override
103 {
104 dealii::Timer timer;
105 timer.restart();
106
107 if(param.large_deformation) // nonlinear problem
108 {
109 VectorType const_vector;
110
111 bool const update_preconditioner =
112 this->param.update_preconditioner &&
113 time_step_number % this->param.update_preconditioner_every_time_steps == 0;
114
115 auto const iter = pde_operator->solve_nonlinear(displacement,
116 const_vector,
117 0.0 /* no acceleration term */,
118 0.0 /* no damping term */,
119 time,
120 update_preconditioner);
121
122 iterations.first += 1;
123 std::get<0>(iterations.second) += std::get<0>(iter);
124 std::get<1>(iterations.second) += std::get<1>(iter);
125
126 if(print_solver_info)
127 {
128 this->pcout << std::endl << "Solve moving mesh problem (nonlinear elasticity):";
129 print_solver_info_nonlinear(pcout, std::get<0>(iter), std::get<1>(iter), timer.wall_time());
130 }
131 }
132 else // linear problem
133 {
134 // calculate right-hand side vector
135 VectorType rhs;
136 pde_operator->initialize_dof_vector(rhs);
137 pde_operator->rhs(rhs, time);
138
139 auto const iter = pde_operator->solve_linear(displacement,
140 rhs,
141 0.0 /* no acceleration term */,
142 0.0 /* no damping term */,
143 time,
144 false /* do not update preconditioner */);
145
146 iterations.first += 1;
147 std::get<1>(iterations.second) += iter;
148
149 if(print_solver_info)
150 {
151 this->pcout << std::endl << "Solve moving mesh problem (linear elasticity):";
152 print_solver_info_linear(pcout, iter, timer.wall_time());
153 }
154 }
155
156 this->initialize_mapping_from_dof_vector(this->mapping_undeformed,
157 displacement,
158 pde_operator->get_dof_handler());
159 }
160
164 void
165 print_iterations() const override
166 {
167 std::vector<std::string> names;
168 std::vector<double> iterations_avg;
169
170 if(param.large_deformation)
171 {
172 names = {"Nonlinear iterations",
173 "Linear iterations (accumulated)",
174 "Linear iterations (per nonlinear it.)"};
175
176 iterations_avg.resize(3);
177 iterations_avg[0] =
178 (double)std::get<0>(iterations.second) / std::max(1., (double)iterations.first);
179 iterations_avg[1] =
180 (double)std::get<1>(iterations.second) / std::max(1., (double)iterations.first);
181 if(iterations_avg[0] > std::numeric_limits<double>::min())
182 iterations_avg[2] = iterations_avg[1] / iterations_avg[0];
183 else
184 iterations_avg[2] = iterations_avg[1];
185 }
186 else // linear
187 {
188 names = {"Linear iterations"};
189 iterations_avg.resize(1);
190 iterations_avg[0] =
191 (double)std::get<1>(iterations.second) / std::max(1., (double)iterations.first);
192 }
193
194 print_list_of_iterations(pcout, names, iterations_avg);
195 }
196
197private:
198 // matrix-free
199 std::shared_ptr<MatrixFreeData<dim, Number>> matrix_free_data;
200 std::shared_ptr<dealii::MatrixFree<dim, Number>> matrix_free;
201
202 // PDE operator
203 std::shared_ptr<Operator<dim, Number>> pde_operator;
204
205 Parameters const & param;
206
207 // store solution of previous time step / iteration so that a good initial
208 // guess is available in the next step, easing convergence or reducing computational
209 // costs by allowing larger tolerances
210 VectorType displacement;
211
212 dealii::ConditionalOStream pcout;
213
214 std::pair<
215 unsigned int /* calls */,
216 std::tuple<unsigned long long, unsigned long long> /* iteration counts {Newton, linear}*/>
217 iterations;
218};
219
220} // namespace Structure
221} // namespace ExaDG
222
223#endif /* INCLUDE_EXADG_GRID_MAPPING_DEFORMATION_STRUCTURE_H_ */
Definition mapping_deformation_base.h:49
std::shared_ptr< dealii::Mapping< dim > const > mapping_undeformed
Definition mapping_deformation_base.h:106
Definition grid.h:40
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 grid.h:84
Definition mapping_deformation_structure.h:44
DeformedMapping(std::shared_ptr< Grid< dim > const > grid, std::shared_ptr< dealii::Mapping< dim > const > mapping_undeformed, std::shared_ptr< MultigridMappings< dim, Number > > const multigrid_mappings_undeformed, std::shared_ptr< BoundaryDescriptor< dim > const > boundary_descriptor, std::shared_ptr< FieldFunctions< dim > const > field_functions, std::shared_ptr< MaterialDescriptor const > material_descriptor, Parameters const &param, std::string const &field, MPI_Comm const &mpi_comm)
Definition mapping_deformation_structure.h:51
void print_iterations() const override
Definition mapping_deformation_structure.h:165
void update(double const time, bool const print_solver_info, types::time_step time_step_number) override
Definition mapping_deformation_structure.h:100
Definition parameters.h:40
Definition driver.cpp:33
Definition boundary_descriptor.h:48
Definition field_functions.h:32