ExaDG
Loading...
Searching...
No Matches
time_int_abm.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_ACOUSTIC_CONSERVATION_EQUATIONS_TIME_INTEGRATION_TIME_INT_ABM_H_
23#define EXADG_ACOUSTIC_CONSERVATION_EQUATIONS_TIME_INTEGRATION_TIME_INT_ABM_H_
24
25// ExaDG
26#include <exadg/acoustic_conservation_equations/postprocessor/postprocessor_interface.h>
27#include <exadg/acoustic_conservation_equations/spatial_discretization/interface.h>
28#include <exadg/acoustic_conservation_equations/user_interface/parameters.h>
29#include <exadg/operators/grid_related_time_step_restrictions.h>
30#include <exadg/time_integration/time_int_abm_base.h>
31#include <exadg/time_integration/time_step_calculation.h>
32
33namespace ExaDG
34{
35namespace Acoustics
36{
37template<typename Number>
38class TimeIntAdamsBashforthMoulton
39 : public TimeIntAdamsBashforthMoultonBase<Interface::SpatialOperator<Number>,
40 dealii::LinearAlgebra::distributed::BlockVector<Number>>
41{
42public:
43 TimeIntAdamsBashforthMoulton(std::shared_ptr<Interface::SpatialOperator<Number>> pde_operator_in,
44 Parameters const & param_in,
45 std::shared_ptr<PostProcessorInterface<Number>> postprocessor_in,
46 MPI_Comm const & mpi_comm_in,
47 bool const is_test_in)
48 : TimeIntAdamsBashforthMoultonBase<Interface::SpatialOperator<Number>,
49 dealii::LinearAlgebra::distributed::BlockVector<Number>>(
50 pde_operator_in,
51 param_in.start_time,
52 param_in.end_time,
53 param_in.max_number_of_time_steps,
54 param_in.order_time_integrator,
55 param_in.start_with_low_order,
56 param_in.adaptive_time_stepping,
57 param_in.restart_data,
58 mpi_comm_in,
59 is_test_in),
60 param(param_in),
61 postprocessor(postprocessor_in),
62 pcout(std::cout, dealii::Utilities::MPI::this_mpi_process(mpi_comm_in) == 0),
63 initial_time_step_size(std::numeric_limits<double>::max())
64 {
65 }
66
67 bool
68 print_solver_info() const final
69 {
70 return param.solver_info_data.write(this->global_timer.wall_time(),
71 this->time - this->start_time,
72 this->time_step_number);
73 }
74
75private:
76 double
77 calculate_time_step_size() final
78 {
79 pcout << std::endl << "Calculation of time step size:" << std::endl << std::endl;
80
81 if(param.calculation_of_time_step_size == TimeStepCalculation::UserSpecified)
82 {
83 initial_time_step_size = calculate_const_time_step(param.time_step_size, param.n_refine_time);
84
85 print_parameter(pcout, "time step size", initial_time_step_size);
86 }
87 else if(param.calculation_of_time_step_size == TimeStepCalculation::CFL)
88 {
89 double const cfl = param.cfl / std::pow(2.0, param.n_refine_time);
90
91 initial_time_step_size = cfl * this->get_underlying_operator().calculate_time_step_cfl();
92
93 this->pcout << std::endl
94 << "Calculation of time step size according to CFL condition:" << std::endl
95 << std::endl;
96 print_parameter(this->pcout, "CFL", cfl);
97 print_parameter(this->pcout, "time step size", initial_time_step_size);
98 }
99
100 return initial_time_step_size;
101 }
102
103 double
104 recalculate_time_step_size() const final
105 {
106 // Currently the time step sice can not vary since the
107 // it depends only on the speed of sound that is
108 // constant over time. This changes once ALE is used
109 // or additional convective terms are considered.
110 return initial_time_step_size;
111 }
112
113 void
114 postprocessing() const final
115 {
116 dealii::Timer timer;
117 timer.restart();
118
119 postprocessor->do_postprocessing(this->get_solution(),
120 this->get_time(),
121 this->time_step_number);
122
123 this->timer_tree->insert({"Timeloop", "Postprocessing"}, timer.wall_time());
124 }
125
126 Parameters const param;
127
128 std::shared_ptr<PostProcessorInterface<Number>> postprocessor;
129
130 dealii::ConditionalOStream pcout;
131
132 // Store initially compute time step size. We do this because the value
133 // will not change, but we need adaptive time-stepping for efficient
134 // sub-time stepping methods.
135 double initial_time_step_size;
136};
137
138} // namespace Acoustics
139} // namespace ExaDG
140
141#endif /* EXADG_ACOUSTIC_CONSERVATION_EQUATIONS_TIME_INTEGRATION_TIME_INT_ABM_H_ */
Definition parameters.h:41
Definition postprocessor_interface.h:34
Definition driver.cpp:33