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