ExaDG
Loading...
Searching...
No Matches
projection_operator.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_OPERATORS_PROJECTION_OPERATOR_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_PROJECTION_OPERATOR_H_
24
25#include <exadg/incompressible_navier_stokes/spatial_discretization/operators/continuity_penalty_operator.h>
26#include <exadg/incompressible_navier_stokes/spatial_discretization/operators/divergence_penalty_operator.h>
27#include <exadg/operators/operator_base.h>
28
29namespace ExaDG
30{
31namespace IncNS
32{
33/*
34 * Combined divergence and continuity penalty operator: applies the operation
35 *
36 * mass operator + dt * divergence penalty operator + dt * continuity penalty operator .
37 *
38 * In detail
39 *
40 * Mass operator: ( v_h , u_h )_Omega^e where
41 *
42 * Divergence penalty operator: ( div(v_h) , tau_div * div(u_h) )_Omega^e
43 *
44 * Continuity penalty operator: ( v_h , tau_conti * jump(u_h) )_dOmega^e, where
45 *
46 * jump(u_h) = u_h^{-} - u_h^{+} or ( (u_h^{-} - u_h^{+})*normal ) * normal
47 *
48 * and
49 *
50 * v_h : test function
51 * u_h : solution
52 */
53
54/*
55 * Operator data.
56 */
57template<int dim>
58struct ProjectionOperatorData : public OperatorBaseData
59{
60 ProjectionOperatorData()
61 : OperatorBaseData(),
62 use_divergence_penalty(true),
63 use_continuity_penalty(true),
64 use_boundary_data(false)
65 {
66 }
67
68 // specify which penalty terms are used
69 bool use_divergence_penalty, use_continuity_penalty;
70
71 bool use_boundary_data;
72
73 std::shared_ptr<BoundaryDescriptorU<dim> const> bc;
74};
75
76template<int dim, typename Number>
77class ProjectionOperator : public OperatorBase<dim, Number, dim>
78{
79private:
80 typedef OperatorBase<dim, Number, dim> Base;
81
82 typedef dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> vector;
83
84 typedef typename Base::VectorType VectorType;
85 typedef typename Base::IntegratorCell IntegratorCell;
86 typedef typename Base::IntegratorFace IntegratorFace;
87
90
91public:
92 typedef Number value_type;
93
94 ProjectionOperator() : velocity(nullptr), time_step_size(1.0)
95 {
96 }
97
98 void
99 initialize(dealii::MatrixFree<dim, Number> const & matrix_free,
100 dealii::AffineConstraints<Number> const & affine_constraints,
101 ProjectionOperatorData<dim> const & data,
102 Operators::DivergencePenaltyKernelData const & div_kernel_data,
103 Operators::ContinuityPenaltyKernelData const & conti_kernel_data);
104
105 void
106 initialize(dealii::MatrixFree<dim, Number> const & matrix_free,
107 dealii::AffineConstraints<Number> const & affine_constraints,
108 ProjectionOperatorData<dim> const & data,
109 std::shared_ptr<DivKernel> div_penalty_kernel,
110 std::shared_ptr<ContiKernel> conti_penalty_kernel);
111
113 get_data() const;
114
116 get_divergence_kernel_data() const;
117
119 get_continuity_kernel_data() const;
120
121 double
122 get_time_step_size() const;
123
124 dealii::LinearAlgebra::distributed::Vector<Number> const &
125 get_velocity() const;
126
127 void
128 update(VectorType const & velocity, double const & dt);
129
130private:
131 void
132 reinit_cell_derived(IntegratorCell & integrator, unsigned int const cell) const final;
133
134 void
135 reinit_face_derived(IntegratorFace & integrator_m,
136 IntegratorFace & integrator_p,
137 unsigned int const face) const final;
138
139 void
140 reinit_boundary_face_derived(IntegratorFace & integrator_m, unsigned int const face) const final;
141
142 void
143 reinit_face_cell_based_derived(IntegratorFace & integrator_m,
144 IntegratorFace & integrator_p,
145 unsigned int const cell,
146 unsigned int const face,
147 dealii::types::boundary_id const boundary_id) const final;
148
149 void
150 do_cell_integral(IntegratorCell & integrator) const final;
151
152 void
153 do_face_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const final;
154
155 void
156 do_face_int_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const final;
157
158 void
159 do_face_ext_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const final;
160
161 void
162 do_boundary_integral(IntegratorFace & integrator_m,
163 OperatorType const & operator_type,
164 dealii::types::boundary_id const & boundary_id) const final;
165
166 ProjectionOperatorData<dim> operator_data;
167
168 VectorType const * velocity;
169 double time_step_size;
170
171 std::shared_ptr<Operators::DivergencePenaltyKernel<dim, Number>> div_kernel;
172 std::shared_ptr<Operators::ContinuityPenaltyKernel<dim, Number>> conti_kernel;
173};
174
175} // namespace IncNS
176} // namespace ExaDG
177
178#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_PROJECTION_OPERATOR_H_ \
179 */
Definition continuity_penalty_operator.h:95
Definition divergence_penalty_operator.h:84
Definition driver.cpp:33
Definition continuity_penalty_operator.h:66
Definition divergence_penalty_operator.h:60
Definition projection_operator.h:59