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