ExaDG
Loading...
Searching...
No Matches
source_term_calculator.h
1#ifndef INCLUDE_EXADG_AERO_ACOUSTIC_SOURCE_TERM_CALCULATOR_H_
2#define INCLUDE_EXADG_AERO_ACOUSTIC_SOURCE_TERM_CALCULATOR_H_
3
4#include <exadg/matrix_free/integrators.h>
5#include <exadg/utilities/lazy_ptr.h>
6
7namespace ExaDG
8{
9namespace AeroAcoustic
10{
11template<int dim>
13{
14 unsigned int dof_index_pressure;
15 unsigned int dof_index_velocity;
16 unsigned int quad_index;
17
18 // density of the underlying fluid.
19 double density;
20
21 // use material or partial temporal derivative of pressure as source term.
22 bool consider_convection;
23
24 // function if blend in is required.
25 bool blend_in;
26 std::shared_ptr<Utilities::SpatialAwareFunction<dim>> blend_in_function;
27};
28
39template<int dim, typename Number>
40class SourceTermCalculator
41{
42 using This = SourceTermCalculator<dim, Number>;
43 using VectorType = dealii::LinearAlgebra::distributed::Vector<Number>;
44
45 using CellIntegratorScalar = CellIntegrator<dim, 1, Number>;
46 using CellIntegratorVector = CellIntegrator<dim, dim, Number>;
47
48 using scalar = dealii::VectorizedArray<Number>;
49 using qpoint = dealii::Point<dim, dealii::VectorizedArray<Number>>;
50
51public:
52 SourceTermCalculator() : matrix_free(nullptr), time(std::numeric_limits<double>::min())
53 {
54 }
55
56 void
57 setup(dealii::MatrixFree<dim, Number> const & matrix_free_in,
58 SourceTermCalculatorData<dim> const & data_in)
59 {
60 matrix_free = &matrix_free_in;
61 data = data_in;
62 }
63
64 void
65 evaluate_integrate(VectorType & dst,
66 dealii::Function<dim> & analytical_source_term,
67 double const evaluation_time)
68 {
69 time = evaluation_time;
70
71 dst.zero_out_ghost_values();
72
73 analytical_source_term.set_time(time);
74
75 matrix_free->cell_loop(&This::compute_source_term, this, dst, analytical_source_term, true);
76 }
77
78
79 void
80 evaluate_integrate(VectorType & dst,
81 VectorType const & velocity_cfd_in,
82 VectorType const & pressure_cfd_in,
83 VectorType const & pressure_cfd_time_derivative_in,
84 double const evaluation_time)
85 {
86 time = evaluation_time;
87
88 dst.zero_out_ghost_values();
89
90 if(data.consider_convection)
91 {
92 velocity_cfd.reset(velocity_cfd_in);
93 velocity_cfd->update_ghost_values();
94
95 pressure_cfd.reset(pressure_cfd_in);
96 pressure_cfd->update_ghost_values();
97 }
98
99 matrix_free->cell_loop(
100 &This::compute_source_term, this, dst, pressure_cfd_time_derivative_in, true);
101 }
102
103private:
104 void
105 compute_source_term(dealii::MatrixFree<dim, Number> const & matrix_free_in,
106 VectorType & dst,
107 VectorType const & dp_cfd_dt,
108 std::pair<unsigned int, unsigned int> const & cell_range) const
109 {
110 CellIntegratorScalar dpdt(matrix_free_in, data.dof_index_pressure, data.quad_index);
111 CellIntegratorScalar p(matrix_free_in, data.dof_index_pressure, data.quad_index);
112 CellIntegratorVector u(matrix_free_in, data.dof_index_velocity, data.quad_index);
113
114 Number rho = static_cast<Number>(data.density);
115
116 auto get_scaling_factor = get_scaling_function();
117
118 for(unsigned int cell = cell_range.first; cell < cell_range.second; ++cell)
119 {
120 dpdt.reinit(cell);
121 dpdt.gather_evaluate(dp_cfd_dt, dealii::EvaluationFlags::values);
122
123 if(data.consider_convection)
124 {
125 p.reinit(cell);
126 p.gather_evaluate(*pressure_cfd, dealii::EvaluationFlags::gradients);
127 u.reinit(cell);
128 u.gather_evaluate(*velocity_cfd, dealii::EvaluationFlags::values);
129
130 for(unsigned int q = 0; q < dpdt.n_q_points; ++q)
131 {
132 scalar flux = -rho * dpdt.get_value(q) + u.get_value(q) * p.get_gradient(q);
133
134 if(data.blend_in)
135 flux *= get_scaling_factor(dpdt.quadrature_point(q));
136
137 dpdt.submit_value(flux, q);
138 }
139 }
140 else
141 {
142 for(unsigned int q = 0; q < dpdt.n_q_points; ++q)
143 {
144 scalar flux = -rho * dpdt.get_value(q);
145
146 if(data.blend_in)
147 flux *= get_scaling_factor(dpdt.quadrature_point(q));
148
149 dpdt.submit_value(flux, q);
150 }
151 }
152
153 dpdt.integrate_scatter(dealii::EvaluationFlags::values, dst);
154 }
155 }
156
157 void
158 compute_source_term(dealii::MatrixFree<dim, Number> const & matrix_free_in,
159 VectorType & dst,
160 dealii::Function<dim> const & analytical_source_term,
161 std::pair<unsigned int, unsigned int> const & cell_range) const
162 {
163 CellIntegratorScalar dpdt(matrix_free_in, data.dof_index_pressure, data.quad_index);
164
165 Number rho = static_cast<Number>(data.density);
166
167 auto get_scaling_factor = get_scaling_function();
168
169 for(unsigned int cell = cell_range.first; cell < cell_range.second; ++cell)
170 {
171 dpdt.reinit(cell);
172
173 for(unsigned int q = 0; q < dpdt.n_q_points; ++q)
174 {
175 scalar flux = -rho * FunctionEvaluator<0, dim, Number>::value(analytical_source_term,
176 dpdt.quadrature_point(q));
177
178 if(data.blend_in)
179 flux *= get_scaling_factor(dpdt.quadrature_point(q));
180
181 dpdt.submit_value(flux, q);
182 }
183
184 dpdt.integrate_scatter(dealii::EvaluationFlags::values, dst);
185 }
186 }
187
188 std::function<scalar(qpoint const &)>
189 get_scaling_function() const
190 {
191 // In case we blend in the source term, we check if the scaling is space dependent. Only in that
192 // case we have to evaluate the function in every equadrature point. Otherwise the scaling
193 // is purely temporal and constant during this function.
194 if(data.blend_in)
195 {
196 AssertThrow(data.blend_in_function != nullptr,
197 dealii::ExcMessage("No blend-in function provided."));
198 }
199
200 bool const space_dependent_scaling =
201 data.blend_in_function != nullptr ? data.blend_in_function->varies_in_space(time) : false;
202 Number const pure_temporal_scaling_factor =
203 (not space_dependent_scaling) ? data.blend_in_function->compute_time_factor(time) : 1.0;
204
205 if(space_dependent_scaling)
206 {
207 return [&](qpoint const & q) {
208 return FunctionEvaluator<0, dim, Number>::value(*data.blend_in_function, q, time);
209 };
210 }
211 else
212 {
213 // capture scaling factor by copy
214 return
215 [pure_temporal_scaling_factor](qpoint const &) { return pure_temporal_scaling_factor; };
216 }
217 }
218
219
220 dealii::MatrixFree<dim, Number> const * matrix_free;
221
223
224 lazy_ptr<VectorType> velocity_cfd;
225 lazy_ptr<VectorType> pressure_cfd;
226
227 double time;
228};
229} // namespace AeroAcoustic
230} // namespace ExaDG
231
232#endif /*INCLUDE_EXADG_AERO_ACOUSTIC_SOURCE_TERM_CALCULATOR_H_*/
Definition lazy_ptr.h:29
Definition driver.cpp:33
Definition source_term_calculator.h:13