ExaDG
Loading...
Searching...
No Matches
operator_pressure_correction.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_PRESSURE_CORRECTION_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATOR_PRESSURE_CORRECTION_H_
24
25#include <exadg/incompressible_navier_stokes/spatial_discretization/operator_projection_methods.h>
26
27namespace ExaDG
28{
29namespace IncNS
30{
31template<int dim, typename Number = double>
32class OperatorPressureCorrection : public OperatorProjectionMethods<dim, Number>
33{
34private:
36 typedef OperatorProjectionMethods<dim, Number> ProjectionBase;
37 typedef OperatorPressureCorrection<dim, Number> This;
38
39 typedef typename Base::VectorType VectorType;
40
41 typedef typename Base::scalar scalar;
42 typedef typename Base::vector vector;
43 typedef typename Base::tensor tensor;
44
45 typedef typename Base::Range Range;
46
47 typedef typename Base::FaceIntegratorU FaceIntegratorU;
48 typedef typename Base::FaceIntegratorP FaceIntegratorP;
49
50public:
51 /*
52 * Constructor.
53 */
54 OperatorPressureCorrection(
55 std::shared_ptr<Grid<dim> const> grid,
56 std::shared_ptr<dealii::Mapping<dim> const> mapping,
57 std::shared_ptr<MultigridMappings<dim, Number>> const multigrid_mappings,
58 std::shared_ptr<BoundaryDescriptor<dim> const> boundary_descriptor,
59 std::shared_ptr<FieldFunctions<dim> const> field_functions,
60 Parameters const & parameters,
61 std::string const & field,
62 MPI_Comm const & mpi_comm);
63
64 /*
65 * Destructor.
66 */
67 virtual ~OperatorPressureCorrection();
68
69 /*
70 * Initializes the inverse pressure mass matrix operator needed for the pressure correction
71 * scheme, as well as the pressure mass operator needed in the ALE case only (where the mass
72 * operator may be evaluated at different times depending on the specific ALE formulation chosen).
73 */
74 void
76
77 void
78 update_after_grid_motion(bool const update_matrix_free) final;
79
80 /*
81 * This function evaluates the nonlinear residual of the steady Navier-Stokes equations (momentum
82 * equation and continuity equation).
83 */
84 void
85 evaluate_nonlinear_residual_steady(VectorType & dst_u,
86 VectorType & dst_p,
87 VectorType const & src_u,
88 VectorType const & src_p,
89 double const & time) const;
90
91 /*
92 * This function applies the linearized momentum operator and is used for throughput measurements.
93 */
94 void
95 apply_momentum_operator(VectorType & dst, VectorType const & src);
96
97 /*
98 * Projection step.
99 */
100
101 // rhs pressure gradient
102 void
103 rhs_pressure_gradient_term_dirichlet_bc_from_dof_vector(VectorType & dst,
104 VectorType const & pressure) const;
105
106 void
107 evaluate_pressure_gradient_term_dirichlet_bc_from_dof_vector(VectorType & dst,
108 VectorType const & src,
109 VectorType const & pressure) const;
110
111 /*
112 * Pressure update step.
113 */
114
115 // apply inverse pressure mass operator
116 void
117 apply_inverse_pressure_mass_operator(VectorType & dst, VectorType const & src) const;
118
119 /*
120 * pressure Poisson equation.
121 */
122 unsigned int
123 solve_pressure(VectorType & dst, VectorType const & src, bool const update_preconditioner) const;
124
125 void
126 rhs_ppe_laplace_add(VectorType & dst, double const & time) const;
127
128 void
129 rhs_ppe_laplace_add_dirichlet_bc_from_dof_vector(VectorType & dst, VectorType const & src) const;
130
131 void
132 interpolate_pressure_dirichlet_bc(VectorType & dst, double const & time) const;
133
134private:
135 /*
136 * Setup of inverse mass operator for pressure.
137 */
138 void
139 setup_inverse_mass_operator_pressure();
140
141 void
142 cell_loop_empty(dealii::MatrixFree<dim, Number> const &,
143 VectorType &,
144 VectorType const &,
145 Range const &) const
146 {
147 }
148
149 void
150 face_loop_empty(dealii::MatrixFree<dim, Number> const &,
151 VectorType &,
152 VectorType const &,
153 Range const &) const
154 {
155 }
156
157 void
158 local_interpolate_pressure_dirichlet_bc_boundary_face(
159 dealii::MatrixFree<dim, Number> const & matrix_free,
160 VectorType & dst,
161 VectorType const & src,
162 Range const & face_range) const;
163
164 InverseMassOperator<dim, 1, Number> inverse_mass_pressure;
165};
166
167} // namespace IncNS
168} // namespace ExaDG
169
170#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATOR_PRESSURE_CORRECTION_H_ \
171 */
Definition grid.h:40
void setup_derived() final
Definition operator_pressure_correction.cpp:56
Definition parameters.h:46
Definition spatial_operator_base.h:68
Definition inverse_mass_operator.h:53
Definition grid.h:84
Definition driver.cpp:33
Definition boundary_descriptor.h:240
Definition field_functions.h:31