Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/callbacks_step/analysis_surface_integral_2d.jl
2055 views
1
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).
2
# Since these FMAs can increase the performance of many numerical algorithms,
3
# we need to opt-in explicitly.
4
# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.
5
@muladd begin
6
#! format: noindent
7
8
@doc raw"""
9
LiftCoefficientPressure2D(aoa, rho_inf, u_inf, l_inf)
10
11
Compute the lift coefficient
12
```math
13
C_{L,p} \coloneqq \frac{\oint_{\partial \Omega} p \boldsymbol n \cdot \psi_L \, \mathrm{d} S}
14
{0.5 \rho_{\infty} U_{\infty}^2 L_{\infty}}
15
```
16
based on the pressure distribution along a boundary.
17
In 2D, the freestream-normal unit vector ``\psi_L`` is given by
18
```math
19
\psi_L \coloneqq \begin{pmatrix} -\sin(\alpha) \\ \cos(\alpha) \end{pmatrix}
20
```
21
where ``\alpha`` is the angle of attack.
22
Supposed to be used in conjunction with [`AnalysisSurfaceIntegral`](@ref)
23
which stores the the to-be-computed variables (for instance `LiftCoefficientPressure2D`)
24
and boundary information.
25
26
- `aoa::Real`: Angle of attack in radians (for airfoils etc.)
27
- `rho_inf::Real`: Free-stream density
28
- `u_inf::Real`: Free-stream velocity
29
- `l_inf::Real`: Reference length of geometry (e.g. airfoil chord length)
30
"""
31
function LiftCoefficientPressure2D(aoa, rho_inf, u_inf, l_inf)
32
# `psi_lift` is the normal unit vector to the freestream direction.
33
# Note: The choice of the normal vector `psi_lift = (-sin(aoa), cos(aoa))`
34
# leads to positive lift coefficients for positive angles of attack for airfoils.
35
# One could also use `psi_lift = (sin(aoa), -cos(aoa))` which results in the same
36
# value, but with the opposite sign.
37
psi_lift = (-sin(aoa), cos(aoa))
38
return LiftCoefficientPressure(ForceState(psi_lift, rho_inf, u_inf, l_inf))
39
end
40
41
@doc raw"""
42
DragCoefficientPressure2D(aoa, rho_inf, u_inf, l_inf)
43
44
Compute the drag coefficient
45
```math
46
C_{D,p} \coloneqq \frac{\oint_{\partial \Omega} p \boldsymbol n \cdot \psi_D \, \mathrm{d} S}
47
{0.5 \rho_{\infty} U_{\infty}^2 L_{\infty}}
48
```
49
based on the pressure distribution along a boundary.
50
In 2D, the freestream-tangent unit vector ``\psi_D`` is given by
51
```math
52
\psi_D \coloneqq \begin{pmatrix} \cos(\alpha) \\ \sin(\alpha) \end{pmatrix}
53
```
54
where ``\alpha`` is the angle of attack.
55
56
Supposed to be used in conjunction with [`AnalysisSurfaceIntegral`](@ref)
57
which stores the the to-be-computed variables (for instance `DragCoefficientPressure2D`)
58
and boundary information.
59
60
- `aoa::Real`: Angle of attack in radians (for airfoils etc.)
61
- `rho_inf::Real`: Free-stream density
62
- `u_inf::Real`: Free-stream velocity
63
- `l_inf::Real`: Reference length of geometry (e.g. airfoil chord length)
64
"""
65
function DragCoefficientPressure2D(aoa, rho_inf, u_inf, l_inf)
66
# `psi_drag` is the unit vector tangent to the freestream direction
67
psi_drag = (cos(aoa), sin(aoa))
68
return DragCoefficientPressure(ForceState(psi_drag, rho_inf, u_inf, l_inf))
69
end
70
71
@doc raw"""
72
LiftCoefficientShearStress2D(aoa, rho_inf, u_inf, l_inf)
73
74
Compute the lift coefficient
75
```math
76
C_{L,f} \coloneqq \frac{\oint_{\partial \Omega} \boldsymbol \tau_w \cdot \psi_L \, \mathrm{d} S}
77
{0.5 \rho_{\infty} U_{\infty}^2 L_{\infty}}
78
```
79
based on the wall shear stress vector ``\tau_w`` along a boundary.
80
In 2D, the freestream-normal unit vector ``\psi_L`` is given by
81
```math
82
\psi_L \coloneqq \begin{pmatrix} -\sin(\alpha) \\ \cos(\alpha) \end{pmatrix}
83
```
84
where ``\alpha`` is the angle of attack.
85
Supposed to be used in conjunction with [`AnalysisSurfaceIntegral`](@ref)
86
which stores the the to-be-computed variables (for instance `LiftCoefficientShearStress2D`)
87
and boundary information.
88
89
- `aoa::Real`: Angle of attack in radians (for airfoils etc.)
90
- `rho_inf::Real`: Free-stream density
91
- `u_inf::Real`: Free-stream velocity
92
- `l_inf::Real`: Reference length of geometry (e.g. airfoil chord length)
93
"""
94
function LiftCoefficientShearStress2D(aoa, rho_inf, u_inf, l_inf)
95
# `psi_lift` is the normal unit vector to the freestream direction.
96
# Note: The choice of the normal vector `psi_lift = (-sin(aoa), cos(aoa))`
97
# leads to negative lift coefficients for airfoils.
98
# One could also use `psi_lift = (sin(aoa), -cos(aoa))` which results in the same
99
# value, but with the opposite sign.
100
psi_lift = (-sin(aoa), cos(aoa))
101
return LiftCoefficientShearStress(ForceState(psi_lift, rho_inf, u_inf, l_inf))
102
end
103
104
@doc raw"""
105
DragCoefficientShearStress2D(aoa, rho_inf, u_inf, l_inf)
106
107
Compute the drag coefficient
108
```math
109
C_{D,f} \coloneqq \frac{\oint_{\partial \Omega} \boldsymbol \tau_w \cdot \psi_D \, \mathrm{d} S}
110
{0.5 \rho_{\infty} U_{\infty}^2 L_{\infty}}
111
```
112
based on the wall shear stress vector ``\tau_w`` along a boundary.
113
In 2D, the freestream-tangent unit vector ``\psi_D`` is given by
114
```math
115
\psi_D \coloneqq \begin{pmatrix} \cos(\alpha) \\ \sin(\alpha) \end{pmatrix}
116
```
117
where ``\alpha`` is the angle of attack.
118
Supposed to be used in conjunction with [`AnalysisSurfaceIntegral`](@ref)
119
which stores the the to-be-computed variables (for instance `DragCoefficientShearStress2D`)
120
and boundary information.
121
122
- `aoa::Real`: Angle of attack in radians (for airfoils etc.)
123
- `rho_inf::Real`: Free-stream density
124
- `u_inf::Real`: Free-stream velocity
125
- `l_inf::Real`: Reference length of geometry (e.g. airfoil chord length)
126
"""
127
function DragCoefficientShearStress2D(aoa, rho_inf, u_inf, l_inf)
128
# `psi_drag` is the unit vector tangent to the freestream direction
129
psi_drag = (cos(aoa), sin(aoa))
130
return DragCoefficientShearStress(ForceState(psi_drag, rho_inf, u_inf, l_inf))
131
end
132
133
# Compute the three components of the 2D symmetric viscous stress tensor
134
# (tau_11, tau_12, tau_22) based on the gradients of the velocity field.
135
# This is required for drag and lift coefficients based on shear stress,
136
# as well as for the non-integrated quantities such as
137
# skin friction coefficient (to be added).
138
function viscous_stress_tensor(u, normal_direction, equations_parabolic,
139
gradients_1, gradients_2)
140
_, dv1dx, dv2dx, _ = convert_derivative_to_primitive(u, gradients_1,
141
equations_parabolic)
142
_, dv1dy, dv2dy, _ = convert_derivative_to_primitive(u, gradients_2,
143
equations_parabolic)
144
145
# Components of viscous stress tensor
146
# (4/3 * (v1)_x - 2/3 * (v2)_y)
147
tau_11 = 4.0 / 3.0 * dv1dx - 2.0 / 3.0 * dv2dy
148
# ((v1)_y + (v2)_x)
149
# stress tensor is symmetric
150
tau_12 = dv1dy + dv2dx # = tau_21
151
# (4/3 * (v2)_y - 2/3 * (v1)_x)
152
tau_22 = 4.0 / 3.0 * dv2dy - 2.0 / 3.0 * dv1dx
153
154
mu = dynamic_viscosity(u, equations_parabolic)
155
156
return mu .* (tau_11, tau_12, tau_22)
157
end
158
159
# 2D viscous stress vector based on contracting the viscous stress tensor
160
# with the normalized `normal_direction` vector.
161
function viscous_stress_vector(u, normal_direction, equations_parabolic,
162
gradients_1, gradients_2)
163
# Normalize normal direction, should point *into* the fluid => *(-1)
164
n_normal = -normal_direction / norm(normal_direction)
165
166
tau_11, tau_12, tau_22 = viscous_stress_tensor(u, normal_direction,
167
equations_parabolic,
168
gradients_1, gradients_2)
169
170
# Viscous stress vector: Stress tensor * normal vector
171
viscous_stress_vector_1 = tau_11 * n_normal[1] + tau_12 * n_normal[2]
172
viscous_stress_vector_2 = tau_12 * n_normal[1] + tau_22 * n_normal[2]
173
174
return (viscous_stress_vector_1, viscous_stress_vector_2)
175
end
176
177
function (lift_coefficient::LiftCoefficientShearStress{RealT, 2})(u, normal_direction,
178
x, t,
179
equations_parabolic,
180
gradients_1,
181
gradients_2) where {RealT <:
182
Real}
183
visc_stress_vector = viscous_stress_vector(u, normal_direction, equations_parabolic,
184
gradients_1, gradients_2)
185
@unpack psi, rho_inf, u_inf, l_inf = lift_coefficient.force_state
186
return (visc_stress_vector[1] * psi[1] + visc_stress_vector[2] * psi[2]) /
187
(0.5f0 * rho_inf * u_inf^2 * l_inf)
188
end
189
190
function (drag_coefficient::DragCoefficientShearStress{RealT, 2})(u, normal_direction,
191
x, t,
192
equations_parabolic,
193
gradients_1,
194
gradients_2) where {RealT <:
195
Real}
196
visc_stress_vector = viscous_stress_vector(u, normal_direction, equations_parabolic,
197
gradients_1, gradients_2)
198
@unpack psi, rho_inf, u_inf, l_inf = drag_coefficient.force_state
199
return (visc_stress_vector[1] * psi[1] + visc_stress_vector[2] * psi[2]) /
200
(0.5f0 * rho_inf * u_inf^2 * l_inf)
201
end
202
203
# 2D version of the `analyze` function for `AnalysisSurfaceIntegral`, i.e.,
204
# `LiftCoefficientPressure` and `DragCoefficientPressure`.
205
function analyze(surface_variable::AnalysisSurfaceIntegral, du, u, t,
206
mesh::P4estMesh{2},
207
equations, dg::DGSEM, cache, semi)
208
@unpack boundaries = cache
209
@unpack surface_flux_values, node_coordinates, contravariant_vectors = cache.elements
210
@unpack weights = dg.basis
211
212
@unpack variable, boundary_symbols = surface_variable
213
@unpack boundary_symbol_indices = semi.boundary_conditions
214
boundary_indices = get_boundary_indices(boundary_symbols, boundary_symbol_indices)
215
216
surface_integral = zero(eltype(u))
217
index_range = eachnode(dg)
218
for boundary in boundary_indices
219
element = boundaries.neighbor_ids[boundary]
220
node_indices = boundaries.node_indices[boundary]
221
direction = indices2direction(node_indices)
222
223
i_node_start, i_node_step = index_to_start_step_2d(node_indices[1], index_range)
224
j_node_start, j_node_step = index_to_start_step_2d(node_indices[2], index_range)
225
226
i_node = i_node_start
227
j_node = j_node_start
228
for node_index in index_range
229
u_node = Trixi.get_node_vars(cache.boundaries.u, equations, dg,
230
node_index, boundary)
231
# Extract normal direction at nodes which points from the elements outwards,
232
# i.e., *into* the structure.
233
normal_direction = get_normal_direction(direction, contravariant_vectors,
234
i_node, j_node, element)
235
236
# Coordinates at a boundary node
237
x = get_node_coords(node_coordinates, equations, dg,
238
i_node, j_node, element)
239
240
# L2 norm of normal direction (contravariant_vector) is the surface element
241
dS = weights[node_index] * norm(normal_direction)
242
243
# Integral over entire boundary surface. Note, it is assumed that the
244
# `normal_direction` is normalized to be a normal vector within the
245
# function `variable` and the division of the normal scaling factor
246
# `norm(normal_direction)` is then accounted for with the `dS` quantity.
247
surface_integral += variable(u_node, normal_direction, x, t, equations) * dS
248
249
i_node += i_node_step
250
j_node += j_node_step
251
end
252
end
253
return surface_integral
254
end
255
256
# 2D version of the `analyze` function for `AnalysisSurfaceIntegral` of viscous, i.e.,
257
# variables that require gradients of the solution variables.
258
# These are for parabolic equations readily available.
259
# Examples are `LiftCoefficientShearStress` and `DragCoefficientShearStress`.
260
function analyze(surface_variable::AnalysisSurfaceIntegral{Variable}, du, u, t,
261
mesh::P4estMesh{2},
262
equations, equations_parabolic,
263
dg::DGSEM, cache, semi,
264
cache_parabolic) where {Variable <: VariableViscous}
265
@unpack boundaries = cache
266
@unpack surface_flux_values, node_coordinates, contravariant_vectors = cache.elements
267
@unpack weights = dg.basis
268
269
@unpack variable, boundary_symbols = surface_variable
270
@unpack boundary_symbol_indices = semi.boundary_conditions
271
boundary_indices = get_boundary_indices(boundary_symbols, boundary_symbol_indices)
272
273
# Additions for parabolic
274
@unpack viscous_container = cache_parabolic
275
@unpack gradients = viscous_container
276
277
gradients_x, gradients_y = gradients
278
279
surface_integral = zero(eltype(u))
280
index_range = eachnode(dg)
281
for boundary in boundary_indices
282
element = boundaries.neighbor_ids[boundary]
283
node_indices = boundaries.node_indices[boundary]
284
direction = indices2direction(node_indices)
285
286
i_node_start, i_node_step = index_to_start_step_2d(node_indices[1], index_range)
287
j_node_start, j_node_step = index_to_start_step_2d(node_indices[2], index_range)
288
289
i_node = i_node_start
290
j_node = j_node_start
291
for node_index in index_range
292
u_node = Trixi.get_node_vars(cache.boundaries.u, equations, dg,
293
node_index, boundary)
294
# Extract normal direction at nodes which points from the elements outwards,
295
# i.e., *into* the structure.
296
normal_direction = get_normal_direction(direction, contravariant_vectors,
297
i_node, j_node, element)
298
299
# Coordinates at a boundary node
300
x = get_node_coords(node_coordinates, equations, dg,
301
i_node, j_node, element)
302
303
# L2 norm of normal direction (contravariant_vector) is the surface element
304
dS = weights[node_index] * norm(normal_direction)
305
306
gradients_1 = get_node_vars(gradients_x, equations_parabolic, dg,
307
i_node, j_node, element)
308
gradients_2 = get_node_vars(gradients_y, equations_parabolic, dg,
309
i_node, j_node, element)
310
311
# Integral over whole boundary surface. Note, it is assumed that the
312
# `normal_direction` is normalized to be a normal vector within the
313
# function `variable` and the division of the normal scaling factor
314
# `norm(normal_direction)` is then accounted for with the `dS` quantity.
315
surface_integral += variable(u_node, normal_direction, x, t,
316
equations_parabolic,
317
gradients_1, gradients_2) * dS
318
319
i_node += i_node_step
320
j_node += j_node_step
321
end
322
end
323
return surface_integral
324
end
325
end # muladd
326
327