ExaDG
Loading...
Searching...
No Matches
acoustics.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_ACOUSTICS_H_
23#define INCLUDE_EXADG_AERO_ACOUSTIC_SINGLE_FIELD_SOLVERS_ACOUSTICS_H_
24
25// Acoustics
26#include <exadg/acoustic_conservation_equations/time_integration/time_int_abm.h>
27#include <exadg/acoustic_conservation_equations/user_interface/application_base.h>
28
29// application
30#include <exadg/aero_acoustic/user_interface/application_base.h>
31
32namespace ExaDG
33{
34namespace AeroAcoustic
35{
36template<int dim, typename Number>
38{
39public:
40 void
41 setup(std::shared_ptr<AcousticsAeroAcoustic::ApplicationBase<dim, Number>> application,
42 MPI_Comm const mpi_comm,
43 bool const is_test)
44 {
45 // setup application
46 application->setup(grid, mapping);
47
48 // setup spatial operator
49 pde_operator = std::make_shared<Acoustics::SpatialOperator<dim, Number>>(
50 grid,
51 mapping,
52 application->get_boundary_descriptor(),
53 application->get_field_functions(),
54 application->get_parameters(),
55 "acoustic",
56 mpi_comm);
57
58 pde_operator->setup();
59
60 // initialize postprocessor
61 postprocessor = application->create_postprocessor();
62 postprocessor->setup(*pde_operator);
63
64 // initialize time integrator
65 time_integrator = std::make_shared<Acoustics::TimeIntAdamsBashforthMoulton<Number>>(
66 pde_operator, application->get_parameters(), postprocessor, mpi_comm, is_test);
67
68 time_integrator->setup(application->get_parameters().restarted_simulation);
69 }
70
71 void
72 advance_multiple_timesteps(double const macro_dt)
73 {
74 // in case the macro time step is smaller than the acoustic timestep, simply run
75 // with macro time step, otherwise, ensure that n_sub_time_steps * sub_dt = macro_dt
76 double const dt = time_integrator->get_time_step_size();
77 double const sub_dt =
78 (macro_dt < dt) ? macro_dt : adjust_time_step_to_hit_end_time(0.0, macro_dt, dt);
79
80
81 // define epsilon dependent on time-step size to avoid
82 // wrong terminations of the sub-timestepping loop due to
83 // round-off errors.
84 double const eps = 0.01 * sub_dt;
85
86 // perform time-steps until we advanced the global time step
87 unsigned int n_sub_time_steps = 0;
88 while(n_sub_time_steps * sub_dt + eps < macro_dt)
89 {
90 time_integrator->set_current_time_step_size(sub_dt);
91 time_integrator->advance_one_timestep();
92 ++n_sub_time_steps;
93 }
94
95 if(time_integrator->started())
96 {
97 ++sub_time_steps.first;
98 sub_time_steps.second += n_sub_time_steps;
99 }
100 }
101
102 double
103 get_average_number_of_sub_time_steps() const
104 {
105 return get_number_of_sub_time_steps() / get_number_of_macro_time_steps();
106 }
107
108 unsigned int
109 get_number_of_macro_time_steps() const
110 {
111 return sub_time_steps.first;
112 }
113
114 unsigned int
115 get_number_of_sub_time_steps() const
116 {
117 return sub_time_steps.second;
118 }
119
120 // grid and mapping
121 std::shared_ptr<Grid<dim>> grid;
122 std::shared_ptr<dealii::Mapping<dim>> mapping;
123
124 // spatial discretization
125 std::shared_ptr<Acoustics::SpatialOperator<dim, Number>> pde_operator;
126
127 // temporal discretization
128 std::shared_ptr<Acoustics::TimeIntAdamsBashforthMoulton<Number>> time_integrator;
129
130 // postprocessor
131 std::shared_ptr<Acoustics::PostProcessorBase<dim, Number>> postprocessor;
132
133private:
134 std::pair<unsigned int /* n_macro_dt */, unsigned long long /* n_sub_dt */> sub_time_steps;
135};
136
137
138} // namespace AeroAcoustic
139} // namespace ExaDG
140
141#endif /* INCLUDE_EXADG_AERO_ACOUSTIC_SINGLE_FIELD_SOLVERS_ACOUSTICS_H_ */
Definition application_base.h:60
Definition acoustics.h:38
Definition driver.cpp:33