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>
34class OperatorDualSplitting : public OperatorProjectionMethods<dim, Number>
35{
36private:
38 typedef OperatorProjectionMethods<dim, Number> ProjectionBase;
39 typedef OperatorDualSplitting<dim, Number> This;
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
120 /*
121 * Fill a DoF vector with velocity Dirichlet values on Dirichlet boundaries.
122 *
123 * Note that this function only works as long as one uses a nodal dealii::FE_DGQ element with
124 * Gauss-Lobatto points. Otherwise, the quadrature formula used in this function does not match
125 * the nodes of the element, and the values injected by this function into the DoF vector are not
126 * the degrees of freedom of the underlying finite element space.
127 */
128 void
129 interpolate_velocity_dirichlet_bc(VectorType & dst, double const & time) const;
130
131private:
132 /*
133 * rhs pressure Poisson equation
134 */
135
136 void
137 cell_loop_empty(dealii::MatrixFree<dim, Number> const &,
138 VectorType &,
139 VectorType const &,
140 Range const &) const
141 {
142 }
143
144 void
145 face_loop_empty(dealii::MatrixFree<dim, Number> const &,
146 VectorType &,
147 VectorType const &,
148 Range const &) const
149 {
150 }
151
152 // rhs PPE: velocity divergence term
153
154 // convective term
155 void
156 local_rhs_ppe_div_term_convective_term_boundary_face(
157 dealii::MatrixFree<dim, Number> const & matrix_free,
158 VectorType & dst,
159 VectorType const & src,
160 Range const & face_range) const;
161
162 // body force term
163 void
164 local_rhs_ppe_div_term_body_forces_boundary_face(
165 dealii::MatrixFree<dim, Number> const & matrix_free,
166 VectorType & dst,
167 VectorType const & src,
168 Range const & face_range) const;
169
170 // Neumann boundary condition term
171
172 // dg_u/dt with numerical time derivative
173 void
174 local_rhs_ppe_nbc_numerical_time_derivative_add_boundary_face(
175 dealii::MatrixFree<dim, Number> const & matrix_free,
176 VectorType & dst,
177 VectorType const & src,
178 Range const & face_range) const;
179
180 // body force term
181 void
182 local_rhs_ppe_nbc_body_force_term_add_boundary_face(
183 dealii::MatrixFree<dim, Number> const & matrix_free,
184 VectorType & dst,
185 VectorType const & src,
186 Range const & face_range) const;
187
188 // convective term
189 void
190 local_rhs_ppe_nbc_convective_add_boundary_face(
191 dealii::MatrixFree<dim, Number> const & matrix_free,
192 VectorType & dst,
193 VectorType const & src,
194 Range const & face_range) const;
195
196 // viscous term
197 void
198 local_rhs_ppe_nbc_viscous_add_boundary_face(dealii::MatrixFree<dim, Number> const & matrix_free,
199 VectorType & dst,
200 VectorType const & src,
201 Range const & face_range) const;
202
203 void
204 local_interpolate_velocity_dirichlet_bc_boundary_face(
205 dealii::MatrixFree<dim, Number> const & matrix_free,
206 VectorType & dst,
207 VectorType const & src,
208 Range const & face_range) const;
209};
210
211} // namespace IncNS
212} // namespace ExaDG
213
214#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATOR_DUAL_SPLITTING_H_ \
215 */
Definition grid.h:40
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