ExaDG
Loading...
Searching...
No Matches
time_int_interpolate_analytical_solution.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2024 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_AERO_ACOUSTIC_SINGLE_FIELD_SOLVERS_ANALYTICAL_TIME_INT_FLUID_H_
23#define INCLUDE_EXADG_AERO_ACOUSTIC_SINGLE_FIELD_SOLVERS_ANALYTICAL_TIME_INT_FLUID_H_
24
25#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf.h>
26#include <exadg/time_integration/push_back_vectors.h>
27#include <exadg/utilities/print_solver_results.h>
28namespace ExaDG
29{
30namespace IncNS
31{
36template<int dim, typename Number>
37class TimeIntInterpolateAnalyticalSolution : public TimeIntBDF<dim, Number>
38{
39 using Base = TimeIntBDF<dim, Number>;
40 using VectorType = typename Base::VectorType;
41 using Operator = SpatialOperatorBase<dim, Number>;
42
43public:
44 TimeIntInterpolateAnalyticalSolution(
45 std::shared_ptr<Operator> operator_in,
46 std::shared_ptr<HelpersALE<dim, Number> const> helpers_ale_in,
47 std::shared_ptr<PostProcessorInterface<Number>> postprocessor_in,
48 Parameters const & param_in,
49 MPI_Comm const & mpi_comm_in,
50 bool const is_test_in)
51 : Base(operator_in, helpers_ale_in, postprocessor_in, param_in, mpi_comm_in, is_test_in),
52 velocity(this->order),
53 pressure(this->order)
54 {
55 }
56
57 void
58 print_iterations() const final
59 {
60 std::vector<std::string> names = {"Interpolation analytical solution"};
61 std::vector<double> iterations_avg{0.0};
62 print_list_of_iterations(this->pcout, names, iterations_avg);
63 }
64
65 VectorType const &
66 get_velocity() const final
67 {
68 return velocity[0];
69 }
70
71
72 VectorType const &
73 get_pressure() const final
74 {
75 return pressure[0];
76 }
77
78 VectorType const &
79 get_velocity_np() const final
80 {
81 return velocity_np;
82 }
83 VectorType const &
84 get_pressure_np() const final
85 {
86 return pressure_np;
87 }
88
89private:
90 VectorType const &
91 get_velocity(unsigned int i /* t_{n-i} */) const final
92 {
93 return velocity[i];
94 }
95
96
97 VectorType const &
98 get_pressure(unsigned int i /* t_{n-i} */) const final
99 {
100 return pressure[i];
101 }
102
103 void
104 allocate_vectors() final
105 {
106 Base::allocate_vectors();
107
108 // velocity
109 for(unsigned int i = 0; i < velocity.size(); ++i)
110 this->operator_base->initialize_vector_velocity(velocity[i]);
111 this->operator_base->initialize_vector_velocity(velocity_np);
112
113 // pressure
114 for(unsigned int i = 0; i < pressure.size(); ++i)
115 this->operator_base->initialize_vector_pressure(pressure[i]);
116 this->operator_base->initialize_vector_pressure(pressure_np);
117 }
118
119
123 void
124 do_timestep_solve() final
125 {
126 this->operator_base->interpolate_analytical_solution(velocity_np,
127 pressure_np,
128 this->get_next_time());
129 }
130
131 void
132 prepare_vectors_for_next_timestep() final
133 {
134 Base::prepare_vectors_for_next_timestep();
135
136 push_back(velocity);
137 velocity[0].swap(velocity_np);
138
139 push_back(pressure);
140 pressure[0].swap(pressure_np);
141 }
142
143 void
144 solve_steady_problem() final
145 {
146 AssertThrow(false, dealii::ExcMessage("TimeIntAnalytic not implemented for steady problem."));
147 }
148
149
150 void
151 initialize_current_solution() final
152 {
153 if(this->param.ale_formulation)
154 this->helpers_ale->move_grid(this->get_time());
155
156 this->operator_base->interpolate_analytical_solution(velocity[0],
157 pressure[0],
158 this->get_time());
159 }
160
161 void
162 initialize_former_multistep_dof_vectors() final
163 {
164 // note that the loop begins with i=1! (we could also start with i=0 but this is not necessary)
165 for(unsigned int i = 1; i < velocity.size(); ++i)
166 {
167 if(this->param.ale_formulation)
168 this->helpers_ale->move_grid(this->get_previous_time(i));
169
170 this->operator_base->interpolate_analytical_solution(velocity[i],
171 pressure[i],
172 this->get_previous_time(i));
173 }
174 }
175
176 void
177 set_velocity(VectorType const & velocity_in, unsigned int const i /* t_{n-i} */) final
178 {
179 velocity[i] = velocity_in;
180 }
181
182 void
183 set_pressure(VectorType const & pressure_in, unsigned int const i /* t_{n-i} */) final
184 {
185 pressure[i] = pressure_in;
186 }
187
188 std::vector<VectorType> velocity;
189 std::vector<VectorType> pressure;
190 VectorType velocity_np;
191 VectorType pressure_np;
192};
193
194} // namespace IncNS
195} // namespace ExaDG
196
197#endif /*INCLUDE_EXADG_AERO_ACOUSTIC_SINGLE_FIELD_SOLVERS_ANALYTICAL_TIME_INT_FLUID_H_*/
Definition lambda_functions_ale.h:40
Definition parameters.h:46
Definition postprocessor_interface.h:34
Definition spatial_operator_base.h:68
Definition driver.cpp:33