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