ExaDG
Loading...
Searching...
No Matches
operator_dual_splitting.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_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATOR_DUAL_SPLITTING_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATOR_DUAL_SPLITTING_H_
24
25// base class
26#include <exadg/incompressible_navier_stokes/spatial_discretization/curl_compute.h>
27#include <exadg/incompressible_navier_stokes/spatial_discretization/operator_projection_methods.h>
28
29namespace ExaDG
30{
31namespace IncNS
32{
33template<int dim, typename Number = double>
35{
36private:
40
41 typedef typename Base::VectorType VectorType;
42
43 typedef typename Base::scalar scalar;
44 typedef typename Base::vector vector;
45 typedef typename Base::tensor tensor;
46
47 typedef typename Base::Range Range;
48
49 typedef typename Base::FaceIntegratorU FaceIntegratorU;
50 typedef typename Base::FaceIntegratorP FaceIntegratorP;
51
52public:
53 /*
54 * Constructor.
55 */
56 OperatorDualSplitting(std::shared_ptr<Grid<dim> const> grid,
57 std::shared_ptr<dealii::Mapping<dim> const> mapping,
58 std::shared_ptr<MultigridMappings<dim, Number>> const multigrid_mappings,
59 std::shared_ptr<BoundaryDescriptor<dim> const> boundary_descriptor,
60 std::shared_ptr<FieldFunctions<dim> const> field_functions,
61 Parameters const & parameters,
62 std::string const & field,
63 MPI_Comm const & mpi_comm);
64
65 /*
66 * Destructor.
67 */
68 virtual ~OperatorDualSplitting();
69
70 /*
71 * Pressure Poisson equation.
72 */
73
74 // rhs pressure: velocity divergence
75 void
76 apply_velocity_divergence_term(VectorType & dst, VectorType const & src) const;
77
78 void
79 rhs_velocity_divergence_term_dirichlet_bc_from_dof_vector(VectorType & dst,
80 VectorType const & velocity) const;
81
82 // rhs pressure Poisson equation: velocity divergence term: body force term
83 void
84 rhs_ppe_div_term_body_forces_add(VectorType & dst, double const & time) const;
85
86 // rhs pressure Poisson equation: velocity divergence term: convective term
87 void
88 rhs_ppe_div_term_convective_term_add(VectorType & dst, VectorType const & src) const;
89
90 // rhs pressure Poisson equation: Neumann BC body force term
91 void
92 rhs_ppe_nbc_body_force_term_add(VectorType & dst, double const & time) const;
93
94 // rhs pressure Poisson equation: Neumann BC numerical time derivative term
95 void
96 rhs_ppe_nbc_numerical_time_derivative_add(VectorType & dst, VectorType const & src) const;
97
98 // rhs pressure Poisson equation: Neumann BC convective term
99 void
100 rhs_ppe_nbc_convective_add(VectorType & dst, VectorType const & src) const;
101
102 // rhs pressure Poisson equation: Neumann BC viscous term
103 void
104 rhs_ppe_nbc_viscous_add(VectorType & dst, VectorType const & src) const;
105
106 void
107 rhs_ppe_laplace_add(VectorType & dst, double const & time) const;
108
109 unsigned int
110 solve_pressure(VectorType & dst, VectorType const & src, bool const update_preconditioner) const;
111
112 /*
113 * Viscous step.
114 */
115
116 void
117 apply_helmholtz_operator(VectorType & dst, VectorType const & src) const;
118
119 void
120 rhs_add_viscous_term(VectorType & dst, double const time) const;
121
122 unsigned int
123 solve_viscous(VectorType & dst,
124 VectorType const & src,
125 bool const & update_preconditioner,
126 double const & scaling_factor_mass);
127
128 /*
129 * Fill a DoF vector with velocity Dirichlet values on Dirichlet boundaries.
130 *
131 * Note that this function only works as long as one uses a nodal dealii::FE_DGQ element with
132 * Gauss-Lobatto points. Otherwise, the quadrature formula used in this function does not match
133 * the nodes of the element, and the values injected by this function into the DoF vector are not
134 * the degrees of freedom of the underlying finite element space.
135 */
136 void
137 interpolate_velocity_dirichlet_bc(VectorType & dst, double const & time) const;
138
139private:
140 void
141 setup_preconditioners_and_solvers() final;
142
143 /*
144 * Setup of Helmholtz solver.
145 */
146 void
147 setup_helmholtz_preconditioner();
148
149 void
150 setup_helmholtz_solver();
151
152 /*
153 * rhs pressure Poisson equation
154 */
155
156 void
157 cell_loop_empty(dealii::MatrixFree<dim, Number> const &,
158 VectorType &,
159 VectorType const &,
160 Range const &) const
161 {
162 }
163
164 void
165 face_loop_empty(dealii::MatrixFree<dim, Number> const &,
166 VectorType &,
167 VectorType const &,
168 Range const &) const
169 {
170 }
171
172 // rhs PPE: velocity divergence term
173
174 // convective term
175 void
176 local_rhs_ppe_div_term_convective_term_boundary_face(
177 dealii::MatrixFree<dim, Number> const & matrix_free,
178 VectorType & dst,
179 VectorType const & src,
180 Range const & face_range) const;
181
182 // body force term
183 void
184 local_rhs_ppe_div_term_body_forces_boundary_face(
185 dealii::MatrixFree<dim, Number> const & matrix_free,
186 VectorType & dst,
187 VectorType const & src,
188 Range const & face_range) const;
189
190 // Neumann boundary condition term
191
192 // dg_u/dt with numerical time derivative
193 void
194 local_rhs_ppe_nbc_numerical_time_derivative_add_boundary_face(
195 dealii::MatrixFree<dim, Number> const & matrix_free,
196 VectorType & dst,
197 VectorType const & src,
198 Range const & face_range) const;
199
200 // body force term
201 void
202 local_rhs_ppe_nbc_body_force_term_add_boundary_face(
203 dealii::MatrixFree<dim, Number> const & matrix_free,
204 VectorType & dst,
205 VectorType const & src,
206 Range const & face_range) const;
207
208 // convective term
209 void
210 local_rhs_ppe_nbc_convective_add_boundary_face(
211 dealii::MatrixFree<dim, Number> const & matrix_free,
212 VectorType & dst,
213 VectorType const & src,
214 Range const & face_range) const;
215
216 // viscous term
217 void
218 local_rhs_ppe_nbc_viscous_add_boundary_face(dealii::MatrixFree<dim, Number> const & matrix_free,
219 VectorType & dst,
220 VectorType const & src,
221 Range const & face_range) const;
222
223 void
224 local_interpolate_velocity_dirichlet_bc_boundary_face(
225 dealii::MatrixFree<dim, Number> const & matrix_free,
226 VectorType & dst,
227 VectorType const & src,
228 Range const & face_range) const;
229
230
231 /*
232 * Viscous step (Helmholtz-like equation).
233 */
234 std::shared_ptr<PreconditionerBase<Number>> helmholtz_preconditioner;
235
236 std::shared_ptr<Krylov::SolverBase<VectorType>> helmholtz_solver;
237};
238
239} // namespace IncNS
240} // namespace ExaDG
241
242#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATOR_DUAL_SPLITTING_H_ \
243 */
Definition grid.h:40
Definition operator_dual_splitting.h:35
Definition operator_projection_methods.h:37
Definition parameters.h:46
Definition spatial_operator_base.h:68
Definition grid.h:84
Definition driver.cpp:33
Definition boundary_descriptor.h:240
Definition field_functions.h:31