Path: blob/main/src/callbacks_step/analysis_surface_integral_2d.jl
2055 views
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).1# Since these FMAs can increase the performance of many numerical algorithms,2# we need to opt-in explicitly.3# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.4@muladd begin5#! format: noindent67@doc raw"""8LiftCoefficientPressure2D(aoa, rho_inf, u_inf, l_inf)910Compute the lift coefficient11```math12C_{L,p} \coloneqq \frac{\oint_{\partial \Omega} p \boldsymbol n \cdot \psi_L \, \mathrm{d} S}13{0.5 \rho_{\infty} U_{\infty}^2 L_{\infty}}14```15based on the pressure distribution along a boundary.16In 2D, the freestream-normal unit vector ``\psi_L`` is given by17```math18\psi_L \coloneqq \begin{pmatrix} -\sin(\alpha) \\ \cos(\alpha) \end{pmatrix}19```20where ``\alpha`` is the angle of attack.21Supposed to be used in conjunction with [`AnalysisSurfaceIntegral`](@ref)22which stores the the to-be-computed variables (for instance `LiftCoefficientPressure2D`)23and boundary information.2425- `aoa::Real`: Angle of attack in radians (for airfoils etc.)26- `rho_inf::Real`: Free-stream density27- `u_inf::Real`: Free-stream velocity28- `l_inf::Real`: Reference length of geometry (e.g. airfoil chord length)29"""30function LiftCoefficientPressure2D(aoa, rho_inf, u_inf, l_inf)31# `psi_lift` is the normal unit vector to the freestream direction.32# Note: The choice of the normal vector `psi_lift = (-sin(aoa), cos(aoa))`33# leads to positive lift coefficients for positive angles of attack for airfoils.34# One could also use `psi_lift = (sin(aoa), -cos(aoa))` which results in the same35# value, but with the opposite sign.36psi_lift = (-sin(aoa), cos(aoa))37return LiftCoefficientPressure(ForceState(psi_lift, rho_inf, u_inf, l_inf))38end3940@doc raw"""41DragCoefficientPressure2D(aoa, rho_inf, u_inf, l_inf)4243Compute the drag coefficient44```math45C_{D,p} \coloneqq \frac{\oint_{\partial \Omega} p \boldsymbol n \cdot \psi_D \, \mathrm{d} S}46{0.5 \rho_{\infty} U_{\infty}^2 L_{\infty}}47```48based on the pressure distribution along a boundary.49In 2D, the freestream-tangent unit vector ``\psi_D`` is given by50```math51\psi_D \coloneqq \begin{pmatrix} \cos(\alpha) \\ \sin(\alpha) \end{pmatrix}52```53where ``\alpha`` is the angle of attack.5455Supposed to be used in conjunction with [`AnalysisSurfaceIntegral`](@ref)56which stores the the to-be-computed variables (for instance `DragCoefficientPressure2D`)57and boundary information.5859- `aoa::Real`: Angle of attack in radians (for airfoils etc.)60- `rho_inf::Real`: Free-stream density61- `u_inf::Real`: Free-stream velocity62- `l_inf::Real`: Reference length of geometry (e.g. airfoil chord length)63"""64function DragCoefficientPressure2D(aoa, rho_inf, u_inf, l_inf)65# `psi_drag` is the unit vector tangent to the freestream direction66psi_drag = (cos(aoa), sin(aoa))67return DragCoefficientPressure(ForceState(psi_drag, rho_inf, u_inf, l_inf))68end6970@doc raw"""71LiftCoefficientShearStress2D(aoa, rho_inf, u_inf, l_inf)7273Compute the lift coefficient74```math75C_{L,f} \coloneqq \frac{\oint_{\partial \Omega} \boldsymbol \tau_w \cdot \psi_L \, \mathrm{d} S}76{0.5 \rho_{\infty} U_{\infty}^2 L_{\infty}}77```78based on the wall shear stress vector ``\tau_w`` along a boundary.79In 2D, the freestream-normal unit vector ``\psi_L`` is given by80```math81\psi_L \coloneqq \begin{pmatrix} -\sin(\alpha) \\ \cos(\alpha) \end{pmatrix}82```83where ``\alpha`` is the angle of attack.84Supposed to be used in conjunction with [`AnalysisSurfaceIntegral`](@ref)85which stores the the to-be-computed variables (for instance `LiftCoefficientShearStress2D`)86and boundary information.8788- `aoa::Real`: Angle of attack in radians (for airfoils etc.)89- `rho_inf::Real`: Free-stream density90- `u_inf::Real`: Free-stream velocity91- `l_inf::Real`: Reference length of geometry (e.g. airfoil chord length)92"""93function LiftCoefficientShearStress2D(aoa, rho_inf, u_inf, l_inf)94# `psi_lift` is the normal unit vector to the freestream direction.95# Note: The choice of the normal vector `psi_lift = (-sin(aoa), cos(aoa))`96# leads to negative lift coefficients for airfoils.97# One could also use `psi_lift = (sin(aoa), -cos(aoa))` which results in the same98# value, but with the opposite sign.99psi_lift = (-sin(aoa), cos(aoa))100return LiftCoefficientShearStress(ForceState(psi_lift, rho_inf, u_inf, l_inf))101end102103@doc raw"""104DragCoefficientShearStress2D(aoa, rho_inf, u_inf, l_inf)105106Compute the drag coefficient107```math108C_{D,f} \coloneqq \frac{\oint_{\partial \Omega} \boldsymbol \tau_w \cdot \psi_D \, \mathrm{d} S}109{0.5 \rho_{\infty} U_{\infty}^2 L_{\infty}}110```111based on the wall shear stress vector ``\tau_w`` along a boundary.112In 2D, the freestream-tangent unit vector ``\psi_D`` is given by113```math114\psi_D \coloneqq \begin{pmatrix} \cos(\alpha) \\ \sin(\alpha) \end{pmatrix}115```116where ``\alpha`` is the angle of attack.117Supposed to be used in conjunction with [`AnalysisSurfaceIntegral`](@ref)118which stores the the to-be-computed variables (for instance `DragCoefficientShearStress2D`)119and boundary information.120121- `aoa::Real`: Angle of attack in radians (for airfoils etc.)122- `rho_inf::Real`: Free-stream density123- `u_inf::Real`: Free-stream velocity124- `l_inf::Real`: Reference length of geometry (e.g. airfoil chord length)125"""126function DragCoefficientShearStress2D(aoa, rho_inf, u_inf, l_inf)127# `psi_drag` is the unit vector tangent to the freestream direction128psi_drag = (cos(aoa), sin(aoa))129return DragCoefficientShearStress(ForceState(psi_drag, rho_inf, u_inf, l_inf))130end131132# Compute the three components of the 2D symmetric viscous stress tensor133# (tau_11, tau_12, tau_22) based on the gradients of the velocity field.134# This is required for drag and lift coefficients based on shear stress,135# as well as for the non-integrated quantities such as136# skin friction coefficient (to be added).137function viscous_stress_tensor(u, normal_direction, equations_parabolic,138gradients_1, gradients_2)139_, dv1dx, dv2dx, _ = convert_derivative_to_primitive(u, gradients_1,140equations_parabolic)141_, dv1dy, dv2dy, _ = convert_derivative_to_primitive(u, gradients_2,142equations_parabolic)143144# Components of viscous stress tensor145# (4/3 * (v1)_x - 2/3 * (v2)_y)146tau_11 = 4.0 / 3.0 * dv1dx - 2.0 / 3.0 * dv2dy147# ((v1)_y + (v2)_x)148# stress tensor is symmetric149tau_12 = dv1dy + dv2dx # = tau_21150# (4/3 * (v2)_y - 2/3 * (v1)_x)151tau_22 = 4.0 / 3.0 * dv2dy - 2.0 / 3.0 * dv1dx152153mu = dynamic_viscosity(u, equations_parabolic)154155return mu .* (tau_11, tau_12, tau_22)156end157158# 2D viscous stress vector based on contracting the viscous stress tensor159# with the normalized `normal_direction` vector.160function viscous_stress_vector(u, normal_direction, equations_parabolic,161gradients_1, gradients_2)162# Normalize normal direction, should point *into* the fluid => *(-1)163n_normal = -normal_direction / norm(normal_direction)164165tau_11, tau_12, tau_22 = viscous_stress_tensor(u, normal_direction,166equations_parabolic,167gradients_1, gradients_2)168169# Viscous stress vector: Stress tensor * normal vector170viscous_stress_vector_1 = tau_11 * n_normal[1] + tau_12 * n_normal[2]171viscous_stress_vector_2 = tau_12 * n_normal[1] + tau_22 * n_normal[2]172173return (viscous_stress_vector_1, viscous_stress_vector_2)174end175176function (lift_coefficient::LiftCoefficientShearStress{RealT, 2})(u, normal_direction,177x, t,178equations_parabolic,179gradients_1,180gradients_2) where {RealT <:181Real}182visc_stress_vector = viscous_stress_vector(u, normal_direction, equations_parabolic,183gradients_1, gradients_2)184@unpack psi, rho_inf, u_inf, l_inf = lift_coefficient.force_state185return (visc_stress_vector[1] * psi[1] + visc_stress_vector[2] * psi[2]) /186(0.5f0 * rho_inf * u_inf^2 * l_inf)187end188189function (drag_coefficient::DragCoefficientShearStress{RealT, 2})(u, normal_direction,190x, t,191equations_parabolic,192gradients_1,193gradients_2) where {RealT <:194Real}195visc_stress_vector = viscous_stress_vector(u, normal_direction, equations_parabolic,196gradients_1, gradients_2)197@unpack psi, rho_inf, u_inf, l_inf = drag_coefficient.force_state198return (visc_stress_vector[1] * psi[1] + visc_stress_vector[2] * psi[2]) /199(0.5f0 * rho_inf * u_inf^2 * l_inf)200end201202# 2D version of the `analyze` function for `AnalysisSurfaceIntegral`, i.e.,203# `LiftCoefficientPressure` and `DragCoefficientPressure`.204function analyze(surface_variable::AnalysisSurfaceIntegral, du, u, t,205mesh::P4estMesh{2},206equations, dg::DGSEM, cache, semi)207@unpack boundaries = cache208@unpack surface_flux_values, node_coordinates, contravariant_vectors = cache.elements209@unpack weights = dg.basis210211@unpack variable, boundary_symbols = surface_variable212@unpack boundary_symbol_indices = semi.boundary_conditions213boundary_indices = get_boundary_indices(boundary_symbols, boundary_symbol_indices)214215surface_integral = zero(eltype(u))216index_range = eachnode(dg)217for boundary in boundary_indices218element = boundaries.neighbor_ids[boundary]219node_indices = boundaries.node_indices[boundary]220direction = indices2direction(node_indices)221222i_node_start, i_node_step = index_to_start_step_2d(node_indices[1], index_range)223j_node_start, j_node_step = index_to_start_step_2d(node_indices[2], index_range)224225i_node = i_node_start226j_node = j_node_start227for node_index in index_range228u_node = Trixi.get_node_vars(cache.boundaries.u, equations, dg,229node_index, boundary)230# Extract normal direction at nodes which points from the elements outwards,231# i.e., *into* the structure.232normal_direction = get_normal_direction(direction, contravariant_vectors,233i_node, j_node, element)234235# Coordinates at a boundary node236x = get_node_coords(node_coordinates, equations, dg,237i_node, j_node, element)238239# L2 norm of normal direction (contravariant_vector) is the surface element240dS = weights[node_index] * norm(normal_direction)241242# Integral over entire boundary surface. Note, it is assumed that the243# `normal_direction` is normalized to be a normal vector within the244# function `variable` and the division of the normal scaling factor245# `norm(normal_direction)` is then accounted for with the `dS` quantity.246surface_integral += variable(u_node, normal_direction, x, t, equations) * dS247248i_node += i_node_step249j_node += j_node_step250end251end252return surface_integral253end254255# 2D version of the `analyze` function for `AnalysisSurfaceIntegral` of viscous, i.e.,256# variables that require gradients of the solution variables.257# These are for parabolic equations readily available.258# Examples are `LiftCoefficientShearStress` and `DragCoefficientShearStress`.259function analyze(surface_variable::AnalysisSurfaceIntegral{Variable}, du, u, t,260mesh::P4estMesh{2},261equations, equations_parabolic,262dg::DGSEM, cache, semi,263cache_parabolic) where {Variable <: VariableViscous}264@unpack boundaries = cache265@unpack surface_flux_values, node_coordinates, contravariant_vectors = cache.elements266@unpack weights = dg.basis267268@unpack variable, boundary_symbols = surface_variable269@unpack boundary_symbol_indices = semi.boundary_conditions270boundary_indices = get_boundary_indices(boundary_symbols, boundary_symbol_indices)271272# Additions for parabolic273@unpack viscous_container = cache_parabolic274@unpack gradients = viscous_container275276gradients_x, gradients_y = gradients277278surface_integral = zero(eltype(u))279index_range = eachnode(dg)280for boundary in boundary_indices281element = boundaries.neighbor_ids[boundary]282node_indices = boundaries.node_indices[boundary]283direction = indices2direction(node_indices)284285i_node_start, i_node_step = index_to_start_step_2d(node_indices[1], index_range)286j_node_start, j_node_step = index_to_start_step_2d(node_indices[2], index_range)287288i_node = i_node_start289j_node = j_node_start290for node_index in index_range291u_node = Trixi.get_node_vars(cache.boundaries.u, equations, dg,292node_index, boundary)293# Extract normal direction at nodes which points from the elements outwards,294# i.e., *into* the structure.295normal_direction = get_normal_direction(direction, contravariant_vectors,296i_node, j_node, element)297298# Coordinates at a boundary node299x = get_node_coords(node_coordinates, equations, dg,300i_node, j_node, element)301302# L2 norm of normal direction (contravariant_vector) is the surface element303dS = weights[node_index] * norm(normal_direction)304305gradients_1 = get_node_vars(gradients_x, equations_parabolic, dg,306i_node, j_node, element)307gradients_2 = get_node_vars(gradients_y, equations_parabolic, dg,308i_node, j_node, element)309310# Integral over whole boundary surface. Note, it is assumed that the311# `normal_direction` is normalized to be a normal vector within the312# function `variable` and the division of the normal scaling factor313# `norm(normal_direction)` is then accounted for with the `dS` quantity.314surface_integral += variable(u_node, normal_direction, x, t,315equations_parabolic,316gradients_1, gradients_2) * dS317318i_node += i_node_step319j_node += j_node_step320end321end322return surface_integral323end324end # muladd325326327