ExaDG
Loading...
Searching...
No Matches
fluid.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2023 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_FLUID_H_
23#define INCLUDE_EXADG_AERO_ACOUSTIC_SINGLE_FIELD_SOLVERS_FLUID_H_
24
25// IncNS
26#include <exadg/incompressible_navier_stokes/postprocessor/postprocessor.h>
27#include <exadg/incompressible_navier_stokes/spatial_discretization/create_operator.h>
28#include <exadg/incompressible_navier_stokes/spatial_discretization/operator_coupled.h>
29#include <exadg/incompressible_navier_stokes/spatial_discretization/operator_dual_splitting.h>
30#include <exadg/incompressible_navier_stokes/spatial_discretization/operator_pressure_correction.h>
31#include <exadg/incompressible_navier_stokes/time_integration/create_time_integrator.h>
32#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf_coupled_solver.h>
33#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf_dual_splitting.h>
34#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf_pressure_correction.h>
35
36// utilities
37#include <exadg/utilities/timer_tree.h>
38
39// application
40#include <exadg/aero_acoustic/user_interface/application_base.h>
41
42namespace ExaDG
43{
44namespace AeroAcoustic
45{
46template<int dim, typename Number>
48{
49public:
50 using VectorType = dealii::LinearAlgebra::distributed::Vector<Number>;
51
53 : timer_tree(std::make_shared<TimerTree>()),
54 adaptive_time_stepping_limiting_factor(std::numeric_limits<double>::min())
55 {
56 }
57
58 void
59 setup(std::shared_ptr<FluidAeroAcoustic::ApplicationBase<dim, Number>> application,
60 MPI_Comm const mpi_comm,
61 bool const is_test)
62 {
63 // setup application
64 application->setup(grid, mapping, multigrid_mappings);
65
66 // initialize pde_operator
67 pde_operator = IncNS::create_operator<dim, Number>(grid,
68 mapping,
69 multigrid_mappings,
70 application->get_boundary_descriptor(),
71 application->get_field_functions(),
72 application->get_parameters(),
73 "fluid",
74 mpi_comm);
75
76 // setup Navier-Stokes operator
77 pde_operator->setup();
78
79 // setup postprocessor
80 postprocessor = application->create_postprocessor();
81 postprocessor->setup(*pde_operator);
82
83 // setup time integrator before calling setup_solvers (this is necessary since the setup
84 // of the solvers depends on quantities such as the time_step_size or gamma0!)
85 AssertThrow(application->get_parameters().solver_type == IncNS::SolverType::Unsteady,
86 dealii::ExcMessage("Invalid parameter in context of fluid-structure interaction."));
87
88 // initialize time_integrator
89 time_integrator = IncNS::create_time_integrator<dim, Number>(pde_operator,
90 nullptr /*no ALE*/,
91 postprocessor,
92 application->get_parameters(),
93 mpi_comm,
94 is_test);
95
96 time_integrator->setup(application->get_parameters().restarted_simulation);
97
98 // initialize vector that stores the pressure time derivative
99 pde_operator->initialize_vector_pressure(pressure_time_derivative);
100
101 // store the adaptive time-stepping limiting factor
102 adaptive_time_stepping_limiting_factor =
103 application->get_parameters().adaptive_time_stepping_limiting_factor;
104 }
105
106 double
107 max_next_time_step_size() const
108 {
109 const auto dt = time_integrator->get_time_step_size();
110 return dt * adaptive_time_stepping_limiting_factor;
111 }
112
113 void
114 advance_one_timestep_and_compute_pressure_time_derivative(bool const update_dpdt)
115 {
116 time_integrator->advance_one_timestep_pre_solve(true);
117 time_integrator->advance_one_timestep_solve();
118
119 // The pressure time derivative has to be computed before the push back
120 // of pressure vectors that is triggered in advance_one_timestep_post_solve()
121 if(update_dpdt)
122 {
123 std::vector<VectorType const *> pressures;
124 std::vector<double> times;
125 time_integrator->get_pressures_and_times_np(pressures, times);
126
127 compute_bdf_time_derivative(pressure_time_derivative, pressures, times);
128 }
129
130 time_integrator->advance_one_timestep_post_solve();
131 }
132
133 VectorType const &
134 get_pressure_time_derivative() const
135 {
136 return pressure_time_derivative;
137 }
138
139 // grid and mapping
140 std::shared_ptr<Grid<dim>> grid;
141 std::shared_ptr<dealii::Mapping<dim>> mapping;
142
143 std::shared_ptr<MultigridMappings<dim, Number>> multigrid_mappings;
144
145 // spatial discretization
146 std::shared_ptr<IncNS::SpatialOperatorBase<dim, Number>> pde_operator;
147
148 // temporal discretization
149 std::shared_ptr<IncNS::TimeIntBDF<dim, Number>> time_integrator;
150
151 // Postprocessor
152 std::shared_ptr<IncNS::PostProcessorBase<dim, Number>> postprocessor;
153
154 // Computation time (wall clock time).
155 std::shared_ptr<TimerTree> timer_tree;
156
157private:
158 // The aeroacoustic source term needs the pressure time derivative.
159 // The update of the vector is performed using
160 // advance_one_timestep_and_compute_pressure_time_derivative(true).
161 VectorType pressure_time_derivative;
162
163 // To be able to estimate if dp/dt has to be evaluated we have to
164 // estimate the maximum time-step size that can be performed after
165 // the current time step. For that we need the limiting factor.
166 double adaptive_time_stepping_limiting_factor;
167};
168
169} // namespace AeroAcoustic
170} // namespace ExaDG
171
172
173
174#endif /* INCLUDE_EXADG_AERO_ACOUSTIC_SINGLE_FIELD_SOLVERS_FLUID_H_ */
Definition application_base.h:186
Definition driver.cpp:33