ExaDG
Loading...
Searching...
No Matches
time_int_gen_alpha_base.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 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_TIME_INTEGRATION_TIME_INT_GEN_ALPHA_BASE_H_
23#define INCLUDE_EXADG_TIME_INTEGRATION_TIME_INT_GEN_ALPHA_BASE_H_
24
25// deal.II
26#include <deal.II/lac/la_parallel_vector.h>
27
28// ExaDG
29#include <exadg/time_integration/enum_types.h>
30#include <exadg/time_integration/time_int_base.h>
31
32namespace ExaDG
33{
34template<typename Number>
36{
37public:
38 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
39
40 TimeIntGenAlphaBase(double const & start_time_,
41 double const & end_time_,
42 unsigned int const max_number_of_time_steps_,
43 double const spectral_radius_,
44 GenAlphaType const & gen_alpha_type_,
45 RestartData const & restart_data_,
46 MPI_Comm const & mpi_comm_,
47 bool const is_test_);
48
49 double
50 get_time_step_size() const final;
51
52 void
53 set_current_time_step_size(double const & time_step_size) final;
54
55protected:
56 double
57 get_scaling_factor_acceleration() const;
58
59 double
60 get_scaling_factor_velocity() const;
61
62 double
63 get_mid_time() const;
64
65 /*
66 * Computes the finite element vector corresponding to the *remainder* of the acceleration term,
67 * which depends on past time step data only. Denoting the acceleration by a, the velocity by v
68 * and the displacement by d, this function returns f1(a^n, v^n, d^n) in
69 * D^2/Dt^2(d) = scaling_factor_acceleration * d^(n+1) + f1(a^n, v^n, d^n)
70 */
71 void
72 compute_const_vector_acceleration_remainder(VectorType & const_vector,
73 VectorType const & displacement_n,
74 VectorType const & velocity_n,
75 VectorType const & acceleration_n) const;
76
77 /*
78 * Computes the finite element vector corresponding to the *remainder* of the velocity term, which
79 * depends on past time step data only. Denoting the acceleration by a, the velocity by v and the
80 * displacement by d, this function returns f2(a^n, v^n, d^n) in
81 * D/Dt(d) = scaling_factor_velocity * d^(n+1) + f2(a^n, v^n, d^n).
82 */
83 void
84 compute_const_vector_velocity_remainder(VectorType & const_vector,
85 VectorType const & displacement_n,
86 VectorType const & velocity_n,
87 VectorType const & acceleration_n) const;
88
89 void
90 update_displacement(VectorType & displacement_np, VectorType const & displacement_n) const;
91
92 void
93 update_velocity(VectorType & velocity_np,
94 VectorType const & displacement_np,
95 VectorType const & displacement_n,
96 VectorType const & velocity_n,
97 VectorType const & acceleration_n) const;
98
99 void
100 update_acceleration(VectorType & acceleration_np,
101 VectorType const & displacement_np,
102 VectorType const & displacement_n,
103 VectorType const & velocity_n,
104 VectorType const & acceleration_n) const;
105
106private:
107 void
108 do_timestep_pre_solve(bool const print_header) final;
109
110 void
111 do_timestep_post_solve() final;
112
113 virtual void
114 prepare_vectors_for_next_timestep() = 0;
115
116 /*
117 * returns whether solver info has to be written in the current time step.
118 */
119 virtual bool
120 print_solver_info() const = 0;
121
122 double spectral_radius;
123 double alpha_m;
124 double alpha_f;
125 double beta;
126 double gamma;
127
128 double time_step;
129};
130
131} // namespace ExaDG
132
133#endif /* INCLUDE_EXADG_TIME_INTEGRATION_TIME_INT_GEN_ALPHA_BASE_H_ */
Definition time_int_base.h:46
Definition time_int_gen_alpha_base.h:36
Definition driver.cpp:33
Definition restart_data.h:38