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