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