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>
41{
43 using VectorType = dealii::LinearAlgebra::distributed::Vector<Number>;
44
45 using CellIntegratorScalar = CellIntegrator<dim, 1, Number>;
46 using CellIntegratorVector = CellIntegrator<dim, dim, Number>;
47
48public:
49 SourceTermCalculator() : matrix_free(nullptr), time(std::numeric_limits<double>::min())
50 {
51 }
52
53 void
54 setup(dealii::MatrixFree<dim, Number> const & matrix_free_in,
55 SourceTermCalculatorData<dim> const & data_in)
56 {
57 matrix_free = &matrix_free_in;
58 data = data_in;
59 }
60
61 void
62 evaluate_integrate(VectorType & dst,
63 VectorType const & velocity_cfd_in,
64 VectorType const & pressure_cfd_in,
65 VectorType const & pressure_cfd_time_derivative_in,
66 double const evaluation_time)
67 {
68 time = evaluation_time;
69
70 dst.zero_out_ghost_values();
71
72 if(data.consider_convection)
73 {
74 velocity_cfd.reset(velocity_cfd_in);
75 velocity_cfd->update_ghost_values();
76
77 pressure_cfd.reset(pressure_cfd_in);
78 pressure_cfd->update_ghost_values();
79 }
80
81 matrix_free->cell_loop(
82 &This::compute_source_term, this, dst, pressure_cfd_time_derivative_in, true);
83 }
84
85private:
86 void
87 compute_source_term(dealii::MatrixFree<dim, Number> const & matrix_free_in,
88 VectorType & dst,
89 VectorType const & dp_cfd_dt,
90 std::pair<unsigned int, unsigned int> const & cell_range) const
91 {
92 CellIntegratorScalar dpdt(matrix_free_in, data.dof_index_pressure, data.quad_index);
93 CellIntegratorScalar p(matrix_free_in, data.dof_index_pressure, data.quad_index);
94 CellIntegratorVector u(matrix_free_in, data.dof_index_velocity, data.quad_index);
95
96 Number rho = static_cast<Number>(data.density);
97
98 // In case we blend in the source term, we check if the scaling is space dependent. Only in that
99 // case we have to evaluate the function in every equadrature point. Otherwise the scaling
100 // is purely temporal and constant during this function.
101 if(data.blend_in)
102 AssertThrow(data.blend_in_function != nullptr,
103 dealii::ExcMessage("No blend-in function provided."));
104
105 bool const space_dependent_scaling =
106 data.blend_in_function != nullptr ? data.blend_in_function->varies_in_space(time) : false;
107 Number const pure_temporal_scaling_factor =
108 (not space_dependent_scaling) ? data.blend_in_function->compute_time_factor(time) : 1.0;
109
110 auto apply_scaling = [&](auto & flux, auto const & q) {
111 if(space_dependent_scaling)
112 flux *= FunctionEvaluator<0, dim, Number>::value(*data.blend_in_function, q, time);
113 else
114 flux *= pure_temporal_scaling_factor;
115 };
116
117 for(unsigned int cell = cell_range.first; cell < cell_range.second; ++cell)
118 {
119 dpdt.reinit(cell);
120 dpdt.gather_evaluate(dp_cfd_dt, dealii::EvaluationFlags::values);
121
122 if(data.consider_convection)
123 {
124 p.reinit(cell);
125 p.gather_evaluate(*pressure_cfd, dealii::EvaluationFlags::gradients);
126 u.reinit(cell);
127 u.gather_evaluate(*velocity_cfd, dealii::EvaluationFlags::values);
128
129 for(unsigned int q = 0; q < dpdt.n_q_points; ++q)
130 {
131 auto flux = -rho * dpdt.get_value(q) + u.get_value(q) * p.get_gradient(q);
132
133 if(data.blend_in)
134 apply_scaling(flux, dpdt.quadrature_point(q));
135
136 dpdt.submit_value(flux, q);
137 }
138 }
139 else
140 {
141 for(unsigned int q = 0; q < dpdt.n_q_points; ++q)
142 {
143 auto flux = -rho * dpdt.get_value(q);
144
145 if(data.blend_in)
146 apply_scaling(flux, dpdt.quadrature_point(q));
147
148 dpdt.submit_value(flux, q);
149 }
150 }
151
152 dpdt.integrate_scatter(dealii::EvaluationFlags::values, dst);
153 }
154 }
155
156 dealii::MatrixFree<dim, Number> const * matrix_free;
157
159
160 lazy_ptr<VectorType> velocity_cfd;
161 lazy_ptr<VectorType> pressure_cfd;
162
163 double time;
164};
165} // namespace AeroAcoustic
166} // namespace ExaDG
167
168#endif /*INCLUDE_EXADG_AERO_ACOUSTIC_SOURCE_TERM_CALCULATOR_H_*/
Definition source_term_calculator.h:41
Definition lazy_ptr.h:29
Definition driver.cpp:33
Definition source_term_calculator.h:13
Definition evaluate_functions.h:40