FLOWPanel.FLOWPanel
— ModuleThree-dimensional panel method for high-Reynolds aerodynamics.
AUTHORSHIP
* Created by : Eduardo J. Alvarez
* Email : Edo.AlvarezR@gmail.com
* Date : Jun 2018 originally as MyPanel.jl
* License : MIT License
FLOWPanel.AbstractBody
— TypeAbstract type AbstractBody{N, E<:AbstractElement}
where N
is the number of element types in this body and E
is an Union containing the N
element types.
Implementations of AbstractBody are expected to have the following fields
grid::GeometricTools.GridTriangleSurface
: Paneled geometrynnodes::Int
: Number of nodesncells::Int
: Number of cellsfields::Vector{String}
: Available fields (solutions)Oaxis::Matrix
: Coordinate system of body w.r.t global (3x3 matrix)O::Vector
: Origin of body w.r.t. global (3-dim vector)strength::Matrix
: Strength of each element of each type (ncells x N matrix)CPoffset::Real
: Control point offset in normal directioncharacteristiclength::Function
: Function for computing the characteristic length of each panel used to offset each control pointkerneloffset::Real
: Kernel offset to avoid singularitieskernelcutoff::Real
: Kernel cutoff to avoid singularities
and the following functions
# Imposes boundary conditions to solve for element strengths.
function solve(self::AbstractBody, Uinfs::Array{<:Real, 2}, args...)
.
.
.
end
# Outputs the body as a vtk file
function save(self::AbstractBody, args...; optargs...)
.
.
.
end
# Returns the dimensions of the system of equations solved by `solve(...)`
# as a tuple `(m, n)`, where `m` is the number of equations and `n` is the
# number of unknowns (usually `m==n`).
function _get_Gdims(self::AbstractBody)
.
.
.
end
# Returns the velocity induced by the body on the targets `targets`. It adds
# the velocity at the i-th target to out[:, i].
function _Uind(self::AbstractBody, targets::Array{<:Real, 2},
out::Array{<:Real, 2}, args...; optargs...)
.
.
.
end
# Returns the potential induced by the body on the targets `targets`. It
# adds the potential at the i-th target to out[i].
function _phi(self::AbstractBody, targets::Array{<:Real, 2},
out::Array{<:Real, 1}, args...; optargs...)
.
.
.
end
FLOWPanel.AbstractLiftingBody
— TypeImplementations of AbstractLiftingBody are expected to have the following fields:
U::Array{Int64,1}
: Indices of all panels along the upper side of the trailing edge.L::Array{Int64,1}
: Indices of all panels along the lower side of the trailing edge.ncellsTE::Int64
: Number of cells along trailing edgennodesTE::Int64
: Number of nodes along trailing edge
NOTE: U and L are assumed to have the same number of points.
in addition to the same properties and functions expected from an AbstractBody implementation. The following functions also need to be implemented.
```julia
Impose boundary conditions to solve for element strengths
function solve(self::AbstractLiftingBody, Uinfs::Array{<:Real, 2}, D::Array{<:Real, 2}, args...) . . . end
Outputs a vtk file with the wake
function _savewake(self::AbstractLiftingBody, filename::String; len::Real=1.0, upper::Bool=true, optargs...) . . . end ```
FLOWPanel.NonLiftingBody
— TypeNonLiftingBody{E::AbstractElement, N}(grid::gt.GridTriangleSurface)
Non-lifting body that is solved using a combination of N panel elements. grid
is the grid surface (paneled geometry).
Properties
nnodes::Int
: Number of nodesncells::Int
: Number of cellsfields::Vector{String}
: Available fields (solutions)Oaxis::Matrix
: Coordinate system of body w.r.t. globalO::Vector
: Origin of body w.r.t. global
FLOWPanel.RigidWakeBody
— TypeRigidWakeBody{E::AbstractElement, N}(grid::gt.GridTriangleSurface, shedding::Matrix{Int})
Lifting body that is solver using a combination of N panel elements and a steady rigid wake. grid
is the grid surface (paneled geometry).
shedding[:, i]
contains the information of the i-th edge along which to shed the wake, where shedding[1, i]
is the linear index of the panel shedding the wake, and shedding[2:3, i]
are the indices of the nodes in that panel that make the edge. Since the wake is typically shed at the edge between two panels, shedding[3, i]
is the index of the partner panel (use -1 if none) and shedding[4:5, i]
are the node indices in that panel that make the edge. The user must ensure that both edges are coincident, and the strength of the wake is equal to the difference between the strengths of both panels.
Properties
nnodes::Int
: Number of nodesncells::Int
: Number of cellsfields::Vector{String}
: Available fields (solutions)Oaxis::Matrix
: Coordinate system of body w.r.t. globalO::Vector
: Origin of body w.r.t. globalncellsTE::Int
: Number of cells along trailing edgennodesTE::Int
: Number of nodes along trailing edge
FLOWPanel.U_boundvortex
— MethodComputes the velocity induced by a bound vortex with beginning and end points [pa1, pa2, pa3]
and [pb1, pb2, pb3]
, respectively, and vortex strength strength
on the targets targets
. It adds the velocity at the i-th target to out[i].
FLOWPanel.U_constant_doublet
— FunctionComputes the velocity induced by a panel of vertices nodes[:, panel]
and constant strength doublet strength
on the targets targets
. It adds the velocity at the i-th target to out[i].
FLOWPanel.U_constant_source
— MethodComputes the velocity induced by a panel of vertices nodes[:, panel]
and constant strength source strength
on the targets targets
. It adds the velocity at the i-th target to out[i].
Implementation of equations in Hess & Smith (1967).
FLOWPanel.U_constant_vortexsheet
— MethodComputes the velocity induced by a panel of vertices nodes[:, panel]
and constant vortex sheet (with strength components gammat
and gammao
) on the targets targets
. It adds the velocity at the i-th target to out[i].
Implementation of equations in Pate's 2017 doctoral dissertation, A Surface Vorticity Method for Wake-Body Interactions, Appendix A.
FLOWPanel.U_semiinfinite_doublet
— FunctionComputes the velocity induced by a semi-infinite doublet panel of strength strength
on the targets targets
. The semi-infinite panel starts at the trailing edge between nodes[:, TE[1]]
and nodes[:, TE[2]]
and extends infinitely in the direction [d1, d2, d3]
. It adds the potential at the i-th target to out[i].
FLOWPanel.U_semiinfinite_horseshoe
— MethodComputes the velocity induced by a semi-infinite horseshoe of strength strength
on the targets targets
. The semi-infinite horseshoe comes from ∞ to nodes[:, TE[1]]
and goes out to ∞ from nodes[:, TE[2]]
, with a bound vortex going from nodes[:, TE[1]]
to nodes[:, TE[2]]
. The direction is of the semi-infinite sections is given by [d1, d2, d3]
. It adds the velocity at the i-th target to out[i].
FLOWPanel.U_semiinfinite_horseshoe
— MethodComputes the velocity induced by a semi-infinite horseshoe of strength strength
on the targets targets
. The semi-infinite horseshoe comes from ∞ to pa and goes out to ∞ from pb, with a bound vortex going from pa to pb. The semi-infinite directions of pa and pb are da and db, respectively. It adds the velocity at the i-th target to out[i].
DEFINITIONS
- pa =
nodes[:, TE[1]]
- pb =
nodes[:, TE[2]]
- da =
[da1, da2, da3]
- db =
[db1, db2, db3]
FLOWPanel.U_semiinfinite_vortex
— MethodComputes the velocity induced by a semi-infinite vortex starting at point [p1, p2, p3]
going out in the unitary direction [d1, d2, d3]
with vortex strength strength
on the targets targets
. It adds the velocity at the i-th target to out[i].
A negative strength makes it such that the vortex comes from infinite and ends at the point.
FLOWPanel.U_vortexring
— MethodU_vortexring(nodes::Matrix, panel::Array{Int}, strength::Real, targets::Matrix, out::Matrix; dot_with=nothing, closed_ring=true, cutoff=1e-14, offset=1e-8)
Computes the velocity induced by a vortex ring panel of vertices nodes[:, panel]
and vortex strength strength
on the targets targets
. It adds the velocity at the i-th target to out[i].
FLOWPanel.Uind!
— MethodUind!(self::AbstractBody, targets, out, args...; optargs...)
Returns the velocity induced by the body on the targets targets
, which is a 3xn matrix. It adds the velocity at the i-th target to out[:, i]
.
FLOWPanel._G_U!
— MethodComputes the geometric matrix (left-hand side matrix of the system of equation) and stores it under G
.
ARGUMENTS
G::Array{T,2}
: Pre-allocated output memory.CPs::Array{T,2}
: Control points.normals::Array{T,2}
: Normal associated to every CP.Das::Array{T,2}
: Unitary direction of semi-infinite vortex at pointa
of each trailing edge panel.Dbs::Array{T,2}
: Unitary direction of semi-infinite vortex at point b of each trailing edge panel.
FLOWPanel._G_U!
— MethodComputes the geometric matrix (left-hand side matrix of the system of equation) and stores it under G
.
ARGUMENTS
G::Array{T,2}
: Pre-allocated output memory.CPs::Array{T,2}
: Control points.normals::Array{T,2}
: Normal associated to every CP.
FLOWPanel._G_U!
— MethodComputes the geometric matrix (left-hand side matrix of the system of equation) and stores it under G
.
ARGUMENTS
G::Array{T,2}
: Pre-allocated output memory.CPs::Array{T,2}
: Control points.normals::Array{T,2}
: Normal associated to every CP.
FLOWPanel._G_U!
— MethodComputes the geometric matrix (left-hand side matrix of the system of equation) and stores it under G
.
ARGUMENTS
G::Array{T,2}
: Pre-allocated output memory.CPs::Array{T,2}
: Control points.normals::Array{T,2}
: Normal associated to every CP.
FLOWPanel._G_phi!
— MethodComputes the geometric matrix (left-hand side matrix of the system of equation) and stores it under G
.
ARGUMENTS
G::Array{T,2}
: Pre-allocated output memory.CPs::Array{T,2}
: Control points.normals::Array{T,2}
: Normal associated to every CP.
FLOWPanel._checkTE
— MethodChecks correction definition of trailing edge
FLOWPanel._count
— MethodCount the number of types in an Union type
FLOWPanel._savewake
— MethodOutputs a vtk file with the wake.
FLOWPanel.add_field
— Methodadd_field(self::AbstractBody, field_name::String, field_type::String,
field_data, entry_type::String)
Adds a new field field_name
to the body (overwriting the field if it already existed).
Expected arguments
field_type=="scalar"
, thenfield_data
is a vector of lengthn
.field_type=="vector"
, thenfield_data
is an array of 3-dim vectors of lengthn
.entry_type=="node"
, thenn=length(field_data)
is the number of nodes in the body andfield_data[i]
is the field value at the i-th node.entry_type=="cell"
, thenn=length(field_data)
is the number of cells in the body andfield_data[i]
is the field value at the i-th cell.entry_type=="system"
, thenn=length(field_data)
is any arbritrary number, andfield_data
is a field for the whole body as a system without any data structure.
FLOWPanel.addfields
— Methodaddfields(body::AbstractBody,
sourcefieldname::String, targetfieldname::String)
Adds field sourcefieldname
to field targetfieldname
.
FLOWPanel.calc_Alu!
— Methodcalc_Alu!(A::AbstractMatrix) -> Alu
Returns the LU decomposition of A
. If A
does not carry Dual nor TrackedReal numbers, computation will be done in-place using A
; hence A
should not be reused for multiple solves or for implicit differentiation (use calc_Alu(A)
and calc_Alu!(Apivot, A)
instead).
FLOWPanel.calc_Alu!
— Methodcalc_Alu!(Apivot::AbstractMatrix, A::AbstractMatrix) -> Alu
Returns the LU decomposition of A
using Apivot
as storage memory to pivot leaving A
unchanged.
FLOWPanel.calc_Alu
— Methodcalc_Alu(A::AbstractMatrix) -> Alu
Returns the LU decomposition of A
.
FLOWPanel.calc_Avalue!
— Methodcalc_Avalue!(Avalue::AbstractMatrix, A::AbstractMatrix)
Copy the primal values of A
into Avalue
.
FLOWPanel.calc_Avalue
— Methodcalc_Avalue(A::AbstractMatrix)
Return the primal values of matrix A
, which is simply A
if the elements of A
are not Dual nor TrackedReal numbers.
FLOWPanel.calc_areas
— Functioncalc_areas(self::AbstractBody)
Calculates the area of every cell in grid
, returning an array with all areas.
See calc_areas!
documentation for more details.
FLOWPanel.calc_areas!
— Functioncalc_areas!(body::AbstractBody, areas::Matrix)
Calculates the area of every cell in body
and stores them in the array areas
.
Output: areas[i]
is the area of the i-th cell (linearly indexed).
FLOWPanel.calc_bc_noflowthrough!
— Methodcalc_bc_noflowthrough!(RHS::AbstractVector,
Us::AbstractMatrix, normals::AbstractMatrix)
Given the velocity and normals at/of each control point, it calculates the right-hand side of the linear system of equations that imposes the no-flow-through condition, which is stored under RHS
.
FLOWPanel.calc_bc_noflowthrough
— Methodcalc_bc_noflowthrough(Us::AbstractMatrix, normals::AbstractMatrix) -> RHS
Given the velocity and normals at/of each control point, it calculates and returns the right-hand side of the linear system of equations that imposes the no-flow-through condition.
FLOWPanel.calc_controlpoints
— Functioncalc_controlpoints(body::AbstractBody)
Calculates the control point of every cell in body
returning a 3xN matrix.
See calc_controlpoints!
documentation for more details.
FLOWPanel.calc_controlpoints!
— Functioncalc_controlpoints!(body::AbstractBody, controlpoints::Matrix, normals::Matrix)
Calculates the control point of every cell in body
and stores them in the 3xN matrix controlpoints
. It uses body.CPoffset
, body.charateristiclength
, and normals
to offset the control points off the surface in the normal direction.
Output: controlpoints[:, i]
is the control point of the i-th cell (linearly indexed).
Use normals = calc_normals(body)
to calculate the normals.
FLOWPanel.calc_elprescribe
— Method`calc_elprescribe(multibody::MultiBody; indices=[1], values=[0.0])`
Automatically calculate the elements to prescribe in a MultiBody, recognizing all bodies that are watertight and prescribing the strength of elements of indices indices
in each body the values values
. Returns array elprescribe
that is used by solve!
.
FLOWPanel.calc_minmax_winding
— Method`calc_minmax_winding(body::Union{NonLiftingBody, AbstractLiftingBody}) ->
(minw, maxw)`
Calculates the winding number of each cell (associated to its corresponding control point) in a body with a Meshes.jl mesh object, and returns both minimum and maximum values.
FLOWPanel.calc_normals
— Functioncalc_normals(self::AbstractBody)
Calculates the normal vector of every cell in grid
returning a 3xN matrix.
See calc_normals!
documentation for more details.
FLOWPanel.calc_normals!
— Functioncalc_normals!(body::AbstractBody, normals::Matrix)
Calculates the normal vector of every cell in body
and stores them in the 3xN matrix normals
.
Output: normals[:, i]
is the normal vector of the i-th cell (linearly indexed).
Normals can be accessed through their (i, j) coordinates (or "Cartesian indices") as follows:
coordinates = (i, j)
ndivscells = get_ndivscells(body)
lin = LinearIndices(Tuple(ndivscells))
FLOWPanel.calc_obliques
— Functioncalc_obliques(self::AbstractBody)
Calculates the oblique vector of every cell in grid
returning a 3xN matrix.
See calc_obliques!
documentation for more details.
FLOWPanel.calc_obliques!
— Functioncalc_obliques!(body::AbstractBody, obliques::Matrix)
Calculates the oblique vector of every cell in body
and stores them in the 3xN matrix obliques
.
Output: obliques[:, i]
is the oblique vector of the i-th cell (linearly indexed).
FLOWPanel.calc_shedding
— Method`calc_shedding(grid::GridTriangleSurface{gt.Meshes.SimpleMesh},
trailingedge::Matrix; tolerance=1e2*eps())`
Given an unstructured grid
and a collection of points (line) trailingedge
, it finds the points in grid
that are closer than tolerance
to the line, and automatically builds a shedding
matrix that can be used to shed the wake from this trailing edge.
Note: It is important that the points in trailingedge
have been previously sorted to be contiguous to each other, otherwise the resulting shedding
might have panels that are not contiguous to each other, fail to recognize panels that are at the trailing edge, or unphysically large trailing vortices.
FLOWPanel.calc_tangents
— Functioncalc_tangents(self::AbstractBody)
Calculates the tangent vector of every cell in grid
returning a 3xN matrix.
See calc_tangents!
documentation for more details.
FLOWPanel.calc_tangents!
— Functioncalc_tangents!(body::AbstractBody, tangents::Matrix)
Calculates the tangent vector of every cell in body
and stores them in the 3xN matrix tangents
.
Output: tangents[:, i]
is the tangent vector of the i-th cell (linearly indexed).
FLOWPanel.calcfield_Cp!
— Methodcalcfield_Cp!(out::Vector, body::AbstractBody, Uref;
U_fieldname="U", fieldname="Cp")
Calculate the pressure coefficient $C_p = 1 - \left(\frac{u}{U_\mathrm{ref}}\right)^2}$, where $u$ is the velocity field named U_fieldname
under body
. The $C_p$ is saved as a field named fieldname
.
The field is calculated in-place and added to out
(hence, make sure that out
starts with all zeroes).
FLOWPanel.calcfield_Cp!
— Methodcalcfield_Cp!(out::Vector, body::AbstractBody, Us, Uref;
U_fieldname="U", fieldname="Cp")
Calculate the pressure coefficient $C_p = 1 - \left(\frac{u}{U_\mathrm{ref}}\right)^2}$, where is the velocity Us
of each control point. The $C_p$ is saved as a field named fieldname
.
The field is calculated in-place and added to out
(hence, make sure that out
starts with all zeroes).
FLOWPanel.calcfield_Cp
— Methodcalcfield_Cp(args...; optargs...)
Similar to calcfield_Cp!
but without in-place calculation (out
is not needed).
FLOWPanel.calcfield_F!
— Methodcalcfield_F!(out::Vector, body::AbstractBody,
areas::Vector, normals::Matrix, Cps::Vector,
Uinf::Number, rho::Number;
fieldname="F")
Calculate the force of each element $F = - C_p \frac{\rho U_\infty}{2} A \hat{\mathbf{n}}$, where $C_p$is calculated from the velocity Cps
at each control point, $A$ is the area of each element given in areas
, and $\hat{\mathbf{n}}$ is the normal of each element given in normals
. $F$ is saved as a field named fieldname
.
The field is calculated in-place and added to out
(hence, make sure that out
starts with all zeroes).
FLOWPanel.calcfield_F!
— Methodcalcfield_F!(out::Vector, body::AbstractBody,
Uinf::Number, rho::Number;
Cp_fieldname="Cp", optargs...
)
Calculate the force of each element $F = - C_p \frac{\rho U_\infty}{2} A \hat{\mathbf{n}}$, where $C_p$is fetched from the field of name Cp_fieldname
, $A$ is the area of each element, and $\hat{\mathbf{n}}$ is the normal of each element. $F$ is saved as a field named fieldname
.
The field is calculated in-place and added to out
(hence, make sure that out
starts with all zeroes).
FLOWPanel.calcfield_F
— Methodcalcfield_F(args...; optargs...)
Similar to calcfield_F!
but without in-place calculation (out
is not needed).
FLOWPanel.calcfield_Ftot!
— Methodcalcfield_Ftot!(out::AbstractVector, body::AbstractBody,
Fs::AbstractMatrix; fieldname="Ftot")
Calculate the integrated force of this body, which is a three-dimensional vector. This is calculated from the force of each element given in Fs
and saved as a field named fieldname
.
The field is calculated in-place and added to out
.
FLOWPanel.calcfield_Ftot!
— Methodcalcfield_Ftot!(out::AbstractVector, body::AbstractBody;
F_fieldname="F", optargs...)
Calculate the integrated force of this body, which is a three-dimensional vector. This is calculated from the force field F_fieldname
and saved as a field named fieldname
.
The field is calculated in-place and added to out
.
FLOWPanel.calcfield_Ftot
— Methodcalcfield_Ftot(body, args...; optargs...) = calcfield_Ftot!(zeros(3), body, args...; optargs...)
Similar to calcfield_Ftot!
but without in-place calculation (out
is not needed).
FLOWPanel.calcfield_LDS!
— Methodcalcfield_LDS!(out, body, Lhat, Dhat; optargs...)
Shat
is calculated automatically from Lhat
and Dhat
,
FLOWPanel.calcfield_LDS!
— Methodcalcfield_LDS!(out::Matrix, body::AbstractBody,
Lhat::Vector, Dhat::Vector, Shat::Vector; F_fieldname="F")
Calculate the integrated force decomposed as lift, drag, and sideslip according to the orthonormal basis Lhat
, Dhat
, Shat
. This is calculated from the force field F_fieldname
.
FLOWPanel.calcfield_LDS!
— Methodcalcfield_LDS!(out::Matrix, body::AbstractBody, Fs::Matrix,
Lhat::Vector, Dhat::Vector, Shat::Vector)
Calculate the integrated force decomposed as lift, drag, and sideslip according to the orthonormal basis Lhat
, Dhat
, Shat
. This is calculated from the force of each element given in Fs
. out[:, 1]
is the lift vector and is saved as the field "L". out[:, 2]
is the drag vector and is saved as the field "D". out[:, 3]
is the sideslip vector and is saved as the field "S".
The field is calculated in-place on out
.
FLOWPanel.calcfield_LDS
— Methodcalcfield_LDS(body, args...; optargs...) = calcfield_LDS!(zeros(3, 3), body, args...; optargs...)
Similar to calcfield_LDS!
but without in-place calculation (out
is not needed).
FLOWPanel.calcfield_Mtot!
— Methodcalcfield_Mtot!(out::AbstractVector, body::AbstractBody,
Xac::AbstractVector, controlpoints::AbstractMatrix,
Fs::AbstractMatrix;
fieldname="Ftot", addfield=true)
Calculate the integrated moment of this body (which is a three-dimensional vector) with respect to the aerodynamic center Xac
. This is calculated from the force and position of each element given in Fs
and controlpoints
, respectively, and saved as a field named fieldname
.
The field is calculated in-place and added to out
.
FLOWPanel.calcfield_Mtot!
— Methodcalcfield_Mtot!(out, body, Xac; F_fieldname="F",
offset=nothing, characteristiclength=nothing, optargs...)
Calculate the integrated moment of this body (which is a three-dimensional vector) with respect to the aerodynamic center Xac
. This is calculated from the force field F_fieldname
and saved as a field named fieldname
.
The field is calculated in-place and added to out
.
FLOWPanel.calcfield_Mtot
— Methodcalcfield_Mtot(body, args...; optargs...) = calcfield_Mtot!(zeros(3), body, args...; optargs...)
Similar to calcfield_Mtot!
but without in-place calculation (out
is not needed).
FLOWPanel.calcfield_U!
— Methodcalcfield_U!(out::Matrix,
sourcebody::AbstractBody, targetbody::AbstractBody,
controlpoints::Matrix, Uinfs::Matrix; fieldname="U")
Calculate the velocity induced by sourcebody
on controlpoints
and save it as a field of name fieldname
under targetbody
. The field includes the freestream velocity Uinfs
.
The field is calculated in-place and added to out
(hence, make sure that out
starts with all zeroes).
FLOWPanel.calcfield_U!
— Methodcalcfield_U!(out::Matrix,
sourcebody::AbstractBody, targetbody::AbstractBody;
offset=nothing, characteristiclength=nothing, optargs...)
Calculate the velocity induced by sourcebody
on control points computed using offset
and characteristiclength
, and save it as a field in targetbody
. The field includes the freestream velocity stored as field "Uinf"
in targetbody
.
The field is calculated in-place and added to out
(hence, make sure that out
starts with all zeroes).
FLOWPanel.calcfield_U
— Methodcalcfield_U(args...; optargs...)
Similar to calcfield_U!
but without in-place calculation (out
is not needed).
FLOWPanel.calcfield_Ugradmu_cell!
— Methodcalcfield_Ugradmu_cell!(out::Matrix, body::AbstractBody;
fieldname="Ugradmu")
Calculate the surface velocity on body
due to changes in the constant doublet strength using the Green-Gauss method and save it as a field of name fieldname
.
The field is calculated in-place and added to out
(hence, make sure that out
starts with all zeroes).
TODO: Avoid the large gradient at the trailing edge recognizing the trailing edge and omiting the neighbor that should be the wake.
FLOWPanel.calcfield_Ugradmu_cell
— Methodcalcfield_Ugradmu_cell(body::AbstractBody; fieldname="Ugradmu")
Similar to calcfield_Ugradmu_cell!
but without in-place calculation (out
is not needed).
FLOWPanel.calcfield_Uoff!
— Methodcalcfield_Uoff!(args...; optargs...) = calcfield_U(args...; optargs..., fieldname="Uoff")
See documentation of calcfield_U!(...)
.
FLOWPanel.calcfield_lmn!
— Methodcalcfield_lmn!(out, body, Xac, lhat, mhat; optargs...)
nhat
is calculated automatically from lhat
and mhat
,
FLOWPanel.calcfield_lmn!
— Methodcalcfield_lmn!(out, body, Xac, lhat, mhat, nhat; F_fieldname="F",
offset=nothing, characteristiclength=nothing, optargs...)
Calculate the integrated moment of this body with respect to the aerodynamic center Xac
and decompose it as rolling, pitching, and yawing moments according to the orthonormal basis lhat
, mhat
, nhat
, repsectively. This is calculated from the force field F_fieldname
.
The field is calculated in-place on out
.
FLOWPanel.calcfield_lmn!
— Methodcalcfield_lmn!(out::Matrix, body::AbstractBody,
Xac::AbstractVector, controlpoints::AbstractMatrix,
Fs::Matrix, lhat::Vector, mhat::Vector, nhat::Vector)
Calculate the integrated moment of this body with respect to the aerodynamic center Xac
and decompose it as rolling, pitching, and yawing moments according to the orthonormal basis lhat
, mhat
, nhat
, repsectively. This is calculated from the force and position of each element given in Fs
and controlpoints
, respectively. out[:, 1]
is the rolling moment vector and is saved as the field "Mroll". out[:, 2]
is the pitching moment vector and is saved as the field "Mpitch". out[:, 3]
is the yawing moment vector and is saved as the field "Myaw".
The field is calculated in-place on out
.
FLOWPanel.calcfield_lmn
— Methodcalcfield_lmn(body, args...; optargs...) = calcfield_lmn!(zeros(3, 3), body, args...; optargs...)
Similar to calcfield_lmn!
but without in-place calculation (out
is not needed).
FLOWPanel.calcfield_sectionalforce!
— Methodcalcfield_sectionalforce!(outFs::Matrix, outpos::Vector,
body::Union{NonLiftingBody, AbstractLiftingBody};
F_fieldname="F", optargs...
)
Calculate the sectional force (a vectorial force per unit span) along the span. This is calculated from the force field F_fieldname
and saved as a field named fieldname
.
The field is calculated in-place on outFs
while the spanwise position of each section is stored under outpos
.
FLOWPanel.calcfield_sectionalforce!
— Methodcalcfield_sectionalforce!(outf::Matrix, outpos::Vector,
body::Union{NonLiftingBody, AbstractLiftingBody},
controlpoints::Matrix, Fs::Matrix;
dimspan=2, dimchord=1,
spandirection=[0, 1, 0],
fieldname="sectionalforce"
)
Calculate the sectional force (a vectorial force per unit span) along the span. This is calculated from the force Fs
and the control points controlpoints
and saved as a field named fieldname
.
The field is calculated in-place on outf
while the spanwise position of each section is stored under outpos
.
FLOWPanel.calcfield_sectionalforce
— Methodcalcfield_sectionalforce(args...; optargs...)
Similar to calcfield_sectionalforce!
but without in-place calculation (outFs
nor outpos
are needed).
FLOWPanel.calcfield_winding
— Method`calcfield_winding(body::Union{NonLiftingBody, AbstractLiftingBody})`
Calculate the winding number of a body with a Meshes.jl mesh object. It adds the winding number of each cell (associated to its corresponding control point) as the field "winding"
.
FLOWPanel.characteristiclength_bbox
— Methodcharacteristiclength_bbox(nodes::Matrix, panel::Vector{Int})
Returns the characteristic length of a panel calculated as the diagonal of the minimum bounding box.
FLOWPanel.characteristiclength_maxdist
— Methodcharacteristiclength_maxdist(nodes::Matrix, panel::Vector{Int})
Returns the characteristic length of a panel calculated as the maximum distance between nodes.
FLOWPanel.characteristiclength_sqrtarea
— Methodcharacteristiclength_sqrtarea(nodes::Matrix, panel::Vector{Int})
Returns the characteristic length of a panel calculated as the square-root of its area.
FLOWPanel.characteristiclength_unitary
— Methodcharacteristiclength_unitary(args...) -> 1
FLOWPanel.check_field
— Methodcheck_field(self::AbstractBody, field_name::String)
Returns true
of the body has the field field_name
. Returns false otherwise.
FLOWPanel.check_solved
— Methodcheck_solved(self::AbstractBody)
Returns true
of the body has been solved. Returns false otherwise.
FLOWPanel.decompose!
— Methoddecompose!(out, V, ihat, jhat)
Similar to decompose!(out, V, ihat, jhat, khat)
, but automatically calculates khat
from ihat
and jhat
.
FLOWPanel.decompose!
— Methoddecompose!(out::Matrix, V::Matrix, ihat::Vector, jhat::Vector, khat::Vector)
Project each column of V
onto the orthonormal bases ihat
, jhat
, khat
. The projection is calculated
FLOWPanel.decompose
— Methoddecompose(V, ihat, jhat)
Similar to decompose!(out, V, ihat, jhat)
but without calculating the projection in-place.
FLOWPanel.find_i
— Method`find_i(body, xtarget::Number, gdim::Int, xdim::Int; xdir=nothing)`
Find the row or column of cells in structured grid body
that is the closest to xtarget
in the xdim
spatial dimension. Use gdim=1
to obtain a row, and gdim=2
to obtain a column.
Alternatively, use an arbitrary direction xdir
in place of xdim
, if given.
Returns (itarget, pos, errmin, lin)
where itarget
is the index of the best candidate row/column, pos
is the position of this row/column projected to xdim
or xdir
, and errmin
is the error between pos
and xtarget
. lin
is the LinearIndices
for the user to iterate over the row/column as lin[itarget, j]
if gdim==1
, or lin[j, itarget]
if gdim==2
.
NOTE:
body
cannot be a MultiBody.
FLOWPanel.generate_loft
— Methodgenerate_loft(bodytype::Type{B}, args...; bodyoptargs=(), optargs...) where {B<:NonLiftingBody}
Generates a lofted non-lifting body. See documentation of GeometricTools.generate_loft
for a description of the arguments of this function.
FLOWPanel.generate_loft_liftbody
— Methodgenerate_loft_liftbody(bodytype::Type{<:AbstractLiftingBody}, args...; optargs...)
Generates a lofted lifting body of type bodytype
. See documentation of GeometricTools.generate_loft
for a description of the arguments of this function.
FLOWPanel.generate_revolution
— Methodgenerate_revolution_(bodytype::Type{B}, args...; bodyoptargs=(), optargs...) where {B<:NonLiftingBody}
Generates a non-lifting body of a body of revolution. See documentation of GeometricTools.surface_revolution
for a description of the arguments of this function.
FLOWPanel.generate_revolution_liftbody
— Methodgenerate_revolution_liftbody(bodytype::Type{<:AbstractLiftingBody}, args...; optargs...)
Generates a lifting body type bodytype
of a body of revolution. See documentation of GeometricTools.surface_revolution
for a description of the arguments of this function.
NOTE: For a complete revolution generating a water-tight grid, the system
of equations solved in the body type RigidWakeBody{VortexRing}
becomes singular (this is because in the absence of an open edge, the solution depends only in the difference between adjacent panels rather than the strengths themselves), leading to vortex ring strengths that are astronomically large. To avoid this situation, use RigidWakeBody{VortexRing, 2}
, which indicates the solver to center the solution around 0.
FLOWPanel.get_body
— MethodReturns the requested body
FLOWPanel.get_body
— MethodReturns the requested body
FLOWPanel.get_cart2lin_cells
— Methodget_cart2lin_cells(self::AbstractBody)
Returns a LinearIndices
that converts the coordinates (or "Cartesian index") of a cell to its linear index.
coordinates = (i, j) # (i, j) coordinates of an arbitrary cell
cart2lin = get_cart2lin_cells(body)
index = cart2lin(coordinates...) # Linear index of the cell
FLOWPanel.get_cart2lin_nodes
— Methodget_cart2lin_nodes(self::AbstractBody)
Returns a LinearIndices
that converts the coordinates (or "Cartesian index") of a node to its linear index.
coordinates = (i, j) # (i, j) coordinates of an arbitrary node
cart2lin = get_cart2lin_nodes(body)
index = cart2lin(coordinates...) # Linear index of the node
FLOWPanel.get_field
— Methodget_field(self::AbstractBody, field_name::String)
Returns the requested field.
FLOWPanel.get_fieldval
— Methodget_field(self::AbstractBody, field_name::String, i::Int)
Returns the requested field value of the i-th cell or node (depending of the field type). Give it _check=false
to skip checking logic for faster processing.
FLOWPanel.get_fieldval
— Methodget_fieldval(self::AbstractBody, field_name::String, coor::Array{Int,1})
Returns the requested field value of the cell or node (depending of the field type) of coordinates coor
(1-indexed).
FLOWPanel.get_linearindex
— Methodget_linearindex(body::Union{NonLiftingBody, AbstractLiftingBody})
Return the LinearIndex of the grid of body
.
FLOWPanel.get_ndivscells
— Methodget_ndivscells(body::AbstractBody)
Returns a tuple with the number of cells in each parametric dimension
FLOWPanel.get_ndivsnodes
— Methodget_ndivsnodes(body::AbstractBody)
Returns a tuple with the number of nodes in each parametric dimension
FLOWPanel.get_normal
— Methodget_normal(body::AbstractBody, i::Int64 or coor::Array{Int64,1})
Returns the normal vector the i-th panel.
FLOWPanel.get_unitvectors
— Methodget_unitvectors(body::AbstractBody, i::Int64 or coor::Array{Int64,1})
Returns the unit vectors t
,n
,o
of the i-th panel, with t
the tanget vector, n
normal, and o
oblique.
FLOWPanel.monitor_slice
— Methodmonitor_slice(body::Union{NonLiftingBody, AbstractLiftingBody},
controlpoints, fieldname, sliceposs;
slicenormal=[0, 1, 0],
row=false,
filter=(x, i)->true,
proc_xfun=(points, vals)->getindex.(eachcol(points), 1),
proc_yfun=(points, vals)->vals,
# ------------- PLOT OUTPUTS ---------------------------
disp_plot=true,
_fig=nothing, _axs=nothing,
ttl=nothing, xscaling=1,
stl="-", xlims=nothing, ylims=nothing,
plot_optargs=[],
figext=".png",
figname="",
fignum=nothing,
# ------------- OTHER OUTPUTS --------------------------
save_path=nothing,
filepref="slice",
num=nothing,
out=[],
output_csv=true,
output_vtk=true
)
Plots and outputs a slice of field fieldname
.
ARGUMENTS
sliceposs::Tuple
: Position of each slice.slicenormal::Vector
: Unit vector normal to the slicing plane.row::Bool
: Slice along first grid dimension if true; second if false.filter::Function
: Filter points along the slice according to this logic.proc_xfun::Function
: Process slice and output x-values for plotsproc_yfun::Function
: Process slice and output y-values for plots
NOTE: Current implementation of find_i
does not work on MultiBody.
FLOWPanel.neighbor
— Methodneighbor(body, ni::Int, ci::Int) -> ncoor
neighbor(body, ni::Int, ccoor) -> ncoor
Returns the Cartesian coordinates ncoor
of the ni
-th neighbor of the cell of linear indexing ci
or coordinates ccoor
.
# Calculate all normals
normals = pnl.calc_normals(body)
# Identify the second neighbor of the 10th cell
ncoor = pnl.neighbor(body, 2, 10)
# Convert Cartesian coordinates to linear indexing
ndivscells = Tuple(collect( 1:(d != 0 ? d : 1) for d in body.grid._ndivscells))
lin = LinearIndices(ndivscells)
ni = lin[ncoor...]
# Fetch the normal of such neighbor
normal = normals[:, ni]
FLOWPanel.phi!
— Methodphi!(self::AbstractBody, targets, out, args...; optargs...)
Returns the potential induced by the body on the targets targets
. It adds the potential at the i-th target to out[:, i]
.
FLOWPanel.phi_constant_doublet
— MethodComputes the potential induced by a panel of vertices nodes[:, panel]
and constant strength doublet strength
on the targets targets
. It adds the potential at the i-th target to out[i].
Implementation of equations in Hess & Smith (1967), using ϕµ = −µ*(n⋅∇ϕσ/σ).
FLOWPanel.phi_constant_source
— MethodComputes the potential induced by a panel of vertices nodes[:, panel]
and constant strength source strength
on the targets targets
. It adds the potential at the i-th target to out[i].
Implementation of equations in Hess & Smith (1967), p. 53.
FLOWPanel.phi_semiinfinite_doublet
— MethodComputes the potential induced by a semi-infinite doublet panel of strength strength
on the targets targets
. The semi-infinite panel starts at the trailing edge between nodes[:, TE[1]]
and nodes[:, TE[2]]
and extends infinitely in the direction [d1, d2, d3]
. It adds the potential at the i-th target to out[i].
Implementation of equations in Moran, J. (1984), Appendix F, p. 445.
FLOWPanel.phi_semiinfinite_doublet
— MethodComputes the potential induced by a semi-infinite horseshoe of strength strength
on the targets targets
. The semi-infinite horseshoe comes from ∞ to pa and goes out to ∞ from pb, with a bound vortex going from pa to pb. The semi-infinite directions of pa and pb are da and db, respectively. It adds the velocity at the i-th target to out[i].
DEFINITIONS
- pa =
nodes[:, TE[1]]
- pb =
nodes[:, TE[2]]
- da =
[da1, da2, da3]
- db =
[db1, db2, db3]
Implementation of equations in Moran, J. (1984), Appendix F, p. 445, splitting the non-planar semi-infinite panel into two planar sections.
FLOWPanel.remove_field
— Methodremove_field(self::AbstractBody, field_name)
Removes field from body.
FLOWPanel.rotate!
— Method`rotate!(body::AbstractBody, roll::Number, pitch::Number, yaw::Number;
translation::Vector=zeros(3), reset_fields::Bool=true)`
Rotates the body by the given axial angles. Nomenclature follows the aircraft convention, with
- roll: rotation about local x-axis
- pitch: rotation about local y-axis
- yaw: rotation about local z-axis
Use keyword argument translation
to also translate the body.
FLOWPanel.rotatetranslate!
— Methodrotatetranslate!(body::AbstractBody, M::Matrix, T::Vector; optargs...)
Rotate and translate body
by rotation matrix M
and translation vector T
.
FLOWPanel.save_base
— Methodsave_base(body::AbstractBody, filename::String; optargs...)
Outputs a vtk file of this body. See GeometricTools.save(::GridExtentions, ::String)
for a description of optional arguments optargs...
.
FLOWPanel.save_controlpoints
— Methodsave_controlpoints(body::AbstractBody, filename::String;
suffix::String="_cp", optargs...)
Outputs a vtk file with the control points of the body along with associated normals and surface velocities.
FLOWPanel.set_coordinatesystem
— Methodset_coordinatesystem(body::AbstractBody, O::Vector, Oaxis::Matrix)
Redefines the local coordinate system of the body, where O
is the new origin and Oaxis
is the matrix of new unit vectors.
FLOWPanel.simplewing
— Methodsimplewing(b, ar, tr, twist_root, twist_tip, lambda, gamma; bodytype=RigidWakeBody, span_NDIVS="automatic", rfl_NDIVS="automatic", airfoil_root="naca6412.dat", airfoil_tip="naca6412.dat", airfoil_path=def_rfl_path)
Generates a symmetric single-section wing.
ARGUMENTS
b::Real
: Span.ar::Real
: Aspect ratio defined as b/c_tip.tr::Real
: Taper ratio defined as ctip/croot.twist_root::Real
: (deg) twist of the root.twist_tip::Real
: (deg) twist of the tip.lambda::Real
: (deg) sweep.gamma::Real
: (deg) dihedral.
OPTIONAL ARGUMENTS
bodytype::Type{LBodyTypes}
: Type of lifting body to generate.span_NDIVS::ndivstype
: Spanwise divisions.rfl_NDIVS::ndivstype
: Chordwise divisions.airfoil_root::String
: File to root airfoil contour.airfoil_tip::String
: File to tip airfoil contour.airfoil_path::String
: Path to airfoil files.
NOTE: See gt.multidscretize for a description of arguments of type ndivstype
. NOTE2: In the current implementation, sweep and dihedral are done about the LE.
FLOWPanel.slicefield
— Methodslicefield(body::AbstractBody, fieldname::String,
position::Number, direction::Vector, row::Bool)
Return a slice of the field fieldname
in body
corresponding to the row or column ("row" is the first dimension of the grid, "column" is the second dimension) that is the closest to position
calculated as the projection of the average cell position in the direction direction
.
Example: For a wing with its span aligned along the y-axis, the pressure along a slice of the wing at the spanwise position y=0.5 is obtained as slicefield(wing, "Cp", 0.5, [0, 1, 0], false)
.
FLOWPanel.slicefield
— Methodslicefield(body::AbstractBody, controlpoints::Matrix,
fieldname::String,
position::Number, direction::Vector, row::Bool)
Same thing, but with the option of providing the control points as to save avoid memory allocation.
FLOWPanel.solve
— Method`solve(body::AbstractBody, Uinfs::Array{<:Real, 2})`
Impose boundary conditions to solve for element strengths. Uinds[:, i]
is the velocity at the i-th control point used in the boundary condition.
FLOWPanel.solve
— Method`solve(body::AbstractBody, Uinfs::Array{<:Real, 2})`
Impose boundary conditions to solve for element strengths. Uinds[:, i]
is the velocity at the i-th control point used in the boundary condition.
Das[:, i]
is the unitary vector pointing in the direction that the semi-infinite wake is shed from the first node of the i-th shedding edge, while Dbs[:, i]
is for the second node. NOTE: These directions are expected to point from the node out to infinite.
FLOWPanel.solve_backslash!
— Methodsolve_backslash!(y::AbstractVector, A::AbstractMatrix, b::AbstractVector)
Solves a linear system of equations of the form $Ay = b$ using the \
operator.
The solution is stored under y
.
FLOWPanel.solve_gmres!
— Methodsolve_gmres!(y::AbstractVector,
A::AbstractMatrix, b::AbstractVector; Avalue=nothing,
atol=1e-8, rtol=1e-8, optargs...)
Solves a linear system of equations of the form $Ay = b$ through the generalized minimal residual (GMRES) method, which is an iterative method in the Krylov subspace.
This iterative method is more efficient than a direct method (solve_backslack!
or solve_ludiv!
) when A
is larger than 3000x3000 or so. Also, iterative methods can trade off accuracy for speed by lowering the tolerance (atol
and rtol
). Optional arguments optargs...
will be passed to Krylov.gmres
.
Differentiating through the solver will require extracting the primal values of A
, which can be provided through the argument Avalue
(this is calculated automatically if not already provided).
The solution is stored under y
.
import FLOWPanel as pnl
Avalue = pnl.calc_Avalue(A)
pnl.solve_gmres!(y, A, b; Avalue=Avalue)
FLOWPanel.solve_ludiv!
— Methodsolve_ludiv!(y::AbstractVector,
A::AbstractMatrix, b::AbstractVector; Alu=nothing)
Solves a linear system of equations of the form $Ay = b$ using the LU decomposition of A
provided under Alu
and LinearAlgebra.ldiv!
. If Alu
is not provided, it will be automatically calculated using LinearAlgebra.lu
.
This method is useful when the system needs to be solved multiple times for different b
vectors since Alu
can be precomputed and re-used. We recommend you use FLOWPanel.calc_Alu
to compute Alu
since this function has been overloaded for Dual and TrackedReal numbers. solve_ludiv!
has also been overloaded with ImplicitAD to efficiently differentiate the linear solver as needed.
The solution is stored under y
.
import FLOWPanel as pnl
Alu = pnl.calc_Alu(A)
pnl.solve_ludiv!(y, A, b; Alu=Alu)