ExaDG
Loading...
Searching...
No Matches
rhs_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 EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_RHS_OPERATOR_H_
23#define EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_RHS_OPERATOR_H_
24
25// ExaDG
26#include <exadg/functions_and_boundary_conditions/evaluate_functions.h>
27#include <exadg/matrix_free/integrators.h>
28#include <exadg/operators/mapping_flags.h>
29
30namespace ExaDG
31{
32namespace IncNS
33{
34namespace Operators
35{
36template<int dim>
37struct RHSKernelData
38{
39 RHSKernelData()
40 : boussinesq_term(false),
41 boussinesq_dynamic_part_only(false),
42 thermal_expansion_coefficient(1.0),
43 reference_temperature(0.0)
44 {
45 }
46
47 std::shared_ptr<dealii::Function<dim>> f;
48
49 // Boussinesq term
50 bool boussinesq_term;
51 bool boussinesq_dynamic_part_only;
52 double thermal_expansion_coefficient;
53 double reference_temperature;
54 std::shared_ptr<dealii::Function<dim>> gravitational_force;
55};
56
57template<int dim, typename Number>
59{
60private:
61 typedef dealii::VectorizedArray<Number> scalar;
62 typedef dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> vector;
63
64 typedef CellIntegrator<dim, dim, Number> Integrator;
65 typedef CellIntegrator<dim, 1, Number> IntegratorScalar;
66
67public:
68 void
69 reinit(RHSKernelData<dim> const & data_in) const
70 {
71 data = data_in;
72 }
73
74 static MappingFlags
75 get_mapping_flags()
76 {
77 MappingFlags flags;
78
79 flags.cells = dealii::update_JxW_values |
80 dealii::update_quadrature_points; // q-points due to rhs function f
81
82 // no face integrals
83
84 return flags;
85 }
86
87 /*
88 * Volume flux, i.e., the term occurring in the volume integral
89 */
90 inline DEAL_II_ALWAYS_INLINE //
91 vector
92 get_volume_flux(Integrator const & integrator,
93 IntegratorScalar const & integrator_temperature,
94 unsigned int const q,
95 Number const & time) const
96 {
97 dealii::Point<dim, scalar> q_points = integrator.quadrature_point(q);
98
99 vector f = FunctionEvaluator<1, dim, Number>::value(*(data.f), q_points, time);
100
101 if(data.boussinesq_term)
102 {
103 vector g =
104 FunctionEvaluator<1, dim, Number>::value(*(data.gravitational_force), q_points, time);
105 scalar T = integrator_temperature.get_value(q);
106 scalar T_ref = data.reference_temperature;
107 // solve only for the dynamic pressure variations
108 if(data.boussinesq_dynamic_part_only)
109 f += g * (-data.thermal_expansion_coefficient * (T - T_ref));
110 else // includes hydrostatic component
111 f += g * (1.0 - data.thermal_expansion_coefficient * (T - T_ref));
112 }
113
114 return f;
115 }
116
117private:
118 mutable RHSKernelData<dim> data;
119};
120
121} // namespace Operators
122
123
124template<int dim>
125struct RHSOperatorData
126{
127 RHSOperatorData() : dof_index(0), quad_index(0), dof_index_scalar(2)
128 {
129 }
130
131 unsigned int dof_index;
132 unsigned int quad_index;
133
134 unsigned int dof_index_scalar;
135
137};
138
139template<int dim, typename Number>
140class RHSOperator
141{
142public:
143 typedef RHSOperator<dim, Number> This;
144
145 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
146
147 typedef std::pair<unsigned int, unsigned int> Range;
148
149 typedef CellIntegrator<dim, dim, Number> Integrator;
150
151 typedef CellIntegrator<dim, 1, Number> IntegratorScalar;
152
153 RHSOperator();
154
155 void
156 initialize(dealii::MatrixFree<dim, Number> const & matrix_free,
157 RHSOperatorData<dim> const & data);
158
159 void
160 evaluate(VectorType & dst, Number const evaluation_time) const;
161
162 void
163 evaluate_add(VectorType & dst, Number const evaluation_time) const;
164
165 void
166 set_temperature(VectorType const & T);
167
168private:
169 void
170 do_cell_integral(Integrator & integrator, IntegratorScalar & integrator_temperature) const;
171
172 void
173 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
174 VectorType & dst,
175 VectorType const & src,
176 Range const & cell_range) const;
177
178 dealii::MatrixFree<dim, Number> const * matrix_free;
179
181
182 mutable double time;
183
185
186 VectorType const * temperature;
187};
188
189} // namespace IncNS
190} // namespace ExaDG
191
192#endif /* EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_RHS_OPERATOR_H_ \
193 */
Definition rhs_operator.h:59
Definition driver.cpp:33
Definition rhs_operator.h:38
Definition rhs_operator.h:126
Definition mapping_flags.h:31