Elements

The library currently has these built-in potential flow elements:

Most functions in the library that act on elements can take either a single element, or a collection of elements. These collections can be represented as an array or a tuple. Arrays should be used when the elements are the same type, for example:

julia> points = Vortex.Point.(rand(ComplexF64, 5), rand(5))
5-element Vector{PotentialFlow.Points.Point{Float64, Float64}}:
 Vortex.Point(0.23603334566204692 + 0.34651701419196046im, 0.5557510873245723)
 Vortex.Point(0.3127069683360675 + 0.00790928339056074im, 0.43710797460962514)
 Vortex.Point(0.4886128300795012 + 0.21096820215853596im, 0.42471785049513144)
 Vortex.Point(0.951916339835734 + 0.9999046588986136im, 0.773223048457377)
 Vortex.Point(0.25166218303197185 + 0.9866663668987996im, 0.2811902322857298)

julia> Elements.impulse(points)
1.3362266530178137 - 1.2821936908564113im

julia> blobs = [Vortex.Blob(rand(ComplexF64), rand(), 0.1) for i in 1:5]
5-element Vector{PotentialFlow.Blobs.Blob{Float64, Float64}}:
 Vortex.Blob(0.20947237319807077 + 0.25137920979222494im, 0.02037486871266725, 0.1)
 Vortex.Blob(0.2877015122756894 + 0.859512136087661im, 0.07695088688120899, 0.1)
 Vortex.Blob(0.6403962459899388 + 0.8735441302706854im, 0.27858242002877853, 0.1)
 Vortex.Blob(0.7513126327861701 + 0.6448833539420931im, 0.07782644396003469, 0.1)
 Vortex.Blob(0.8481854810000327 + 0.0856351682044918im, 0.5532055454580578, 0.1)

julia> Elements.impulse(blobs)
0.41217890550975256 - 0.7325028967929701im

Knowing that every element has the same type allows the compiler to perform more aggressive optimizations. Tuples are used when we want to mix and match different element types. For example:

julia> sys = (points, blobs);

julia> Elements.impulse(sys)
1.7484055585275664 - 2.0146965876493814im

This rest of this page documents the data types that represent these elements and some key functions that act on them. For more detailed examples, please refer to the Jupyter notebooks.

Built-in Types

PotentialFlow.Vortex.PointType
Vortex.Point(z::ComplexF64, Γ::Float64)

A point vortex located at z with circulation Γ.

A new point vortex can be created from an existing one by treating the existing point vortex as a function and passing in the parameter you want to change as keyword arguments. For example,

julia> p = Vortex.Point(1.0, 1.0)
Vortex.Point(1.0 + 0.0im, 1.0)

julia> p()
Vortex.Point(1.0 + 0.0im, 1.0)

julia> p(Γ = 2.0)
Vortex.Point(1.0 + 0.0im, 2.0)
source
PotentialFlow.Vortex.BlobType
Vortex.Blob(z::ComplexF64, Γ::Float64, δ::Float64)

A regularized point vortex located at z with circulation Γ and blob radius δ.

A new vortex blob can be created from an existing one by treating the existing blob as a function and passing in the parameter you want to change as keyword arguments. For example,

julia> b = Vortex.Blob(1.0, 1.0, 0.1)
Vortex.Blob(1.0 + 0.0im, 1.0, 0.1)

julia> b()
Vortex.Blob(1.0 + 0.0im, 1.0, 0.1)

julia> b(Γ = 2.0, δ = 0.01)
Vortex.Blob(1.0 + 0.0im, 2.0, 0.01)
source
PotentialFlow.Vortex.SheetType
Vortex.Sheet <: Elements.Element

A vortex sheet represented by vortex blob control points

Fields

  • blobs: the underlying array of vortex blobs
  • Ss: the cumulated sum of circulation starting from the first control point
  • δ: the blob radius of all the vortex blobs
  • zs: a mapped array that accesses the position of each control point

Constructors:

  • Vortex.Sheet(blobs, Γs, δ)
  • Vortex.Sheet(zs, Γs, δ) where zs is an array of positions for the control points
source
PotentialFlow.Source.PointType
Source.Point(z::ComplexF64, S::Float64)

A point source located at z with strength S.

A new point source can be created from an existing one by treating the existing source as a function and passing in the parameter you want to change as keyword arguments. For example,

julia> p = Source.Point(1.0, 1.0)
Source.Point(1.0 + 0.0im, 1.0)

julia> p()
Source.Point(1.0 + 0.0im, 1.0)

julia> p(S = 2.0)
Source.Point(1.0 + 0.0im, 2.0)
source
PotentialFlow.Source.BlobType
Source.Blob(z::ComplexF64, S::Float64, δ::Float64)

A regularized point source located at z with strength S and blob radius δ.

A new blob source can be created from an existing one by treating the existing blob as a function and passing in the parameter you want to change as keyword arguments. For example,

julia> b = Source.Blob(1.0, 1.0, 0.1)
Source.Blob(1.0 + 0.0im, 1.0, 0.1)

julia> b()
Source.Blob(1.0 + 0.0im, 1.0, 0.1)

julia> b(S = 2.0, δ = 0.01)
Source.Blob(1.0 + 0.0im, 2.0, 0.01)
source
PotentialFlow.Plates.PlateType
Plate <: Elements.Element

An infinitely thin, flat plate, represented as a bound vortex sheet

Constructors

  • Plate(N, L, c, α)
source
PotentialFlow.Bodies.ConformalBodyType
ConformalBody <: Elements.Element

Generates a body from a conformal map. This might be a Schwarz-Christoffel map, in which case the constructor is supplied a polygon, or it might be a power- series map, in which case the constructor is given a set of complex coefficients.

Example

julia> p = Bodies.Polygon([-1.0,0.2,1.0,-1.0],[-1.0,-1.0,0.5,1.0])
Polygon with 4 vertices at
             (-1.0,-1.0) (0.2,-1.0) (1.0,0.5) (-1.0,1.0)
             interior angles/π = [0.5, 0.656, 0.422, 0.422]

julia> Bodies.ConformalBody(p)
Body generated by: Schwarz-Christoffel map of unit circle to exterior of polygon with 4 vertices

  centroid at 0.0 + 0.0im
  angle 0.0

julia> a1 = 1; b1 = 0.1; ccoeff = ComplexF64[0.5(a1+b1),0,0.5(a1-b1)];

julia> Bodies.ConformalBody(ccoeff,ComplexF64(1.0),π/4)
Body generated by: Power series map

  centroid at 1.0 + 0.0im
  angle 0.7854
source

Element Properties

PotentialFlow.Elements.positionFunction
Elements.position(src::Element)

Returns the complex position of a potential flow element. This is a required method for all Element types.

Example

julia> point = Vortex.Point(1.0 + 0.0im, 1.0);

julia> Elements.position(point)
1.0 + 0.0im

julia> points = Vortex.Point.([1.0im, 2.0im], 1.0);

julia> Elements.position.(points)
2-element Array{Complex{Float64},1}:
 0.0 + 1.0im
 0.0 + 2.0im
source
PotentialFlow.Elements.circulationFunction
Elements.circulation(src)

Returns the total circulation contained in src.

Example

julia> points = Vortex.Point.([1.0im, 2.0im], [1.0, 2.0]);

julia> blobs = Vortex.Blob.([1.0im, 2.0im], [1.0, 2.0], 0.1);

julia> Elements.circulation(points[1])
1.0

julia> Elements.circulation(points)
3.0

julia> Elements.circulation((points, blobs))
6.0

julia> Elements.circulation.(points)
2-element Array{Float64,1}:
 1.0
 2.0

julia> Elements.circulation.((points, blobs))
(3.0, 3.0)

julia> Elements.circulation(Source.Point(rand(), rand()))
0.0

julia> Elements.circulation(Source.Blob(rand(), rand(), rand()))
0.0
source
PotentialFlow.Elements.fluxFunction
Elements.flux(src)

Returns the flux through a unit circle induced by src.

Example

julia> points = Source.Point.([1.0im, 2.0im], [1.0, 2.0]);

julia> blobs = Source.Blob.([1.0im, 2.0im], [1.0, 2.0], 0.1);

julia> Elements.flux(points[1])
1.0

julia> Elements.flux((points, blobs))
6.0

julia> Elements.flux.(points)
2-element Array{Float64,1}:
 1.0
 2.0

julia> Elements.flux.((points, blobs))
(3.0, 3.0)

julia> Elements.flux(Vortex.Point(rand(), rand()))
0.0

julia> Elements.flux(Vortex.Blob(rand(), rand(), rand()))
0.0
source
PotentialFlow.Elements.impulseFunction
Elements.impulse(src)

Return the aerodynamic impulse of src about (0,0):

\[P := \int \boldsymbol{x} \times \boldsymbol{\omega}\,\mathrm{d}A.\]

This is a required method for all vortex types.

Example

julia> sys = (Vortex.Point(1.0im, π), Vortex.Blob(1.0im, -π, 0.1));

julia> Elements.impulse(sys[1])
3.141592653589793 + 0.0im

julia> Elements.impulse(sys)
0.0 + 0.0im
source

Methods on Vortex Sheets

PotentialFlow.Sheets.append_segment!Function
Sheets.append_segment!(sheet::Sheet, z, Γ)

Append a new segment with circulation Γ extending from the end of the sheet to z.

Example

julia> sheet = Vortex.Sheet(0:0.1:1, 0.0:10, 0.2)
Vortex Sheet: L ≈ 1.0, Γ = 10.0, δ = 0.2

julia> sheet.blobs[end]
Vortex.Blob(1.0 + 0.0im, 0.5, 0.2)

julia> Sheets.append_segment!(sheet, 1.1, 2.0)

julia> sheet
Vortex Sheet: L ≈ 1.1, Γ = 12.0, δ = 0.2

julia> sheet.blobs[end]
Vortex.Blob(1.1 + 0.0im, 1.0, 0.2)
source
PotentialFlow.Sheets.truncate!Function
Sheets.truncate!(sheet, n::Int)

Remove segments 0:n from sheet, and return the circulation in those segments.

Example

julia> sheet = Vortex.Sheet(0:0.1:1, 0.0:10, 0.2)
Vortex Sheet: L ≈ 1.0, Γ = 10.0, δ = 0.2

julia> Sheets.truncate!(sheet, 5)
4.0
source
PotentialFlow.Sheets.redistribute_points!Function
Sheets.redistribute_points!(sheet, zs, Γs)

Returns the modified sheet with replacement control points at positions zs and strength Γs.

julia> sheet = Vortex.Sheet(0:0.1:1, 0.0:10, 0.2)
Vortex Sheet: L ≈ 1.0, Γ = 10.0, δ = 0.2

julia> sys = (sheet,)
(Vortex Sheet: L ≈ 1.0, Γ = 10.0, δ = 0.2,)

julia> Sheets.redistribute_points!(sheet, 0:0.2:2, 0.0:0.5:5)
Vortex Sheet: L ≈ 2.0, Γ = 5.0, δ = 0.2

julia> sys[1]
Vortex Sheet: L ≈ 2.0, Γ = 5.0, δ = 0.2
source
PotentialFlow.Sheets.remeshFunction
Sheets.remesh(sheet, Δs::Float64 , params::Tuple = ())

Uniformly redistribute the control points of the sheet to have a nominal spacing of Δs. Material quantities that should be redistributed along with the control points can be passed in as elements of params.

Returns the tuple (z₌, Γ₌, L [, p₌]) where

  • z₌ is an array with the positions of the uniformly distributed points
  • Γ₌ is circulation interpolated onto z₌
  • L is total length of the sheet
  • p₌ is a tuple containing the material quantities from params interpolated onto z₌

Example

julia> sheet = Vortex.Sheet(0:0.1:1, 0.0:10, 0.2)
Vortex Sheet: L ≈ 1.0, Γ = 10.0, δ = 0.2

julia> age = collect(10.0:-1:0);

julia> Sheets.remesh(sheet, 0.2, (age, ))
(Complex{Float64}[0.0+0.0im, 0.25+0.0im, 0.5+0.0im, 0.75+0.0im, 1.0+0.0im], [0.0, 2.5, 5.0, 7.5, 10.0], 1.0, ([10.0, 7.5, 5.0, 2.5, 0.0],))
source
PotentialFlow.Sheets.remesh!Function
Sheets.remesh!(sheet::Sheet, Δs::Float64, params::Tuple = ())

Same as Sheets.remesh, except sheet is replaced internally by a uniformly interpolated control points. Returns the tuple (sheet, L, p₌) where

  • sheet is the modified sheet
  • L is total length of the sheet
  • p₌ is a tuple containing the material quantities from params interpolated onto the new control points of sheet
julia> sheet = Vortex.Sheet(0:0.1:1, 0.0:10, 0.2)
Vortex Sheet: L ≈ 1.0, Γ = 10.0, δ = 0.2

julia> age = collect(10.0:-1:0);

julia> Sheets.remesh!(sheet, 0.2, (age,));

julia> Elements.position.(sheet.blobs)
5-element Array{Complex{Float64},1}:
  0.0 + 0.0im
 0.25 + 0.0im
  0.5 + 0.0im
 0.75 + 0.0im
  1.0 + 0.0im

julia> age
5-element Array{Float64,1}:
 10.0
  7.5
  5.0
  2.5
  0.0
source
PotentialFlow.Sheets.split!Function
Sheets.split!(sheet, n::Int)

Remove segments 0:n from sheet, and return those segments as a new sheet.

Example

julia> sheet = Vortex.Sheet(0:0.1:1, 0.0:10, 0.2)
Vortex Sheet: L ≈ 1.0, Γ = 10.0, δ = 0.2

julia> sheet₋ = Sheets.split!(sheet, 5)
Vortex Sheet: L ≈ 0.4, Γ = 4.0, δ = 0.2

julia> sheet
Vortex Sheet: L ≈ 0.6, Γ = 6.0, δ = 0.2
source
PotentialFlow.Sheets.filter!Function
Sheets.filter!(sheet, Δs, Δf[, params])

Redistribute and filter the control points of a vortex sheet

Arguments

  • sheet: the vortex sheet to be modified
  • Δs: the nominal spacing between the uniform points
  • Δf: the minimum length scale that the filter should allow to pass through
  • params: an optional tuple of vectors containing material properties

Returns

If params is passed in, then its vectors will be overwritten by their interpolated values on the new control points, and the function returns the tuple (sheet, params). Otherwise, it returns (sheet, ())

source
PotentialFlow.Sheets.filter_position!Function
Sheets.filter_position!(s, Δf, L = arclength(z₌))

Filter out any length scales in s that is smaller than Δf, storing the result back in s. s can be either a vector of complex positions, or a Vortex.Sheet.

source
PotentialFlow.Sheets.arclengthFunction
Sheets.arclength(s)

Compute the polygonal arc length of s, where s can be either an vector of complex numbers or a Vortex.Sheet.

Example

```jldoctest julia> sheet = Vortex.Sheet(0:0.1:1, 0.0:10, 0.2) Vortex Sheet: L ≈ 1.0, Γ = 10.0, δ = 0.2

julia> Sheets.arclength(sheet) 1.0

source
PotentialFlow.Sheets.arclengthsFunction
Sheets.arclengths(s)

Cumulative sum of the polygonal arc length of s, where s can be either an vector of complex numbers or a Vortex.Sheet.

Example

julia> sheet = Vortex.Sheet(0:0.1:1, 0.0:10, 0.2)
Vortex Sheet: L ≈ 1.0, Γ = 10.0, δ = 0.2

julia> Sheets.arclengths(sheet)
11-element Array{Float64,1}:
 0.0
 0.1
 0.2
 0.3
 0.4
 0.5
 0.6
 0.7
 0.8
 0.9
 1.0
source

Methods on Plates

PotentialFlow.Plates.edgesFunction
edges(plate)

Return the coordinates of the leading and trailing edges

Example

julia> p = Plate(128, 1.0, 0, π/4)
Plate: N = 128, L = 1.0, c = 0.0 + 0.0im, α = 45.0ᵒ
       LESP = 0.0, TESP = 0.0

julia> Plates.edges(p)
(0.3535533905932738 + 0.35355339059327373im, -0.3535533905932738 - 0.35355339059327373im)
source
PotentialFlow.Plates.enforce_no_flow_through!Function
enforce_no_flow_through!(p::Plate, motion, elements, t)

Update the plate, p, to enforce the no-flow-through condition given ambient vortex elements, elements, and while moving with kinematics specified by motion.

Example

julia> plate = Plate(128, 2.0, 0.0, π/3)
Plate: N = 128, L = 2.0, c = 0.0 + 0.0im, α = 60.0ᵒ
       LESP = 0.0, TESP = 0.0

julia> motion = Plates.RigidBodyMotion(1.0, 0.0);

julia> point = Vortex.Point(0.0 + 2im, 1.0);

julia> Plates.enforce_no_flow_through!(plate, motion, point, 0.0)

julia> plate
Plate: N = 128, L = 2.0, c = 0.0 + 0.0im, α = 60.0ᵒ
       LESP = 1.27, TESP = -1.93
source
PotentialFlow.Plates.vorticity_fluxFunction
vorticity_flux(p::Plate, v₁, v₂,
               lesp = 0.0, tesp = 0.0,
               ∂C₁ = Vector{ComplexF64}(undef, plate.N),
               ∂C₂ = Vector{ComplexF64}(undef, plate.N)[,clamp_constraints=false])

Return strengths of new vortex elements that satisfies edge suction parameters. For a given edge, if the current suction parameter is less than the criticial suction parameter, then no vorticity is released. If it is higher, however, vorticity will be released so that the suction parameter equals the critical value. If clamp_constraints=true, and if one of the constraints is currently satisfied, then it fixes the constraint at that critical value while it solves for new strengths. The default is false, in which it does not assign any new strength to a vortex element if the associated constraint is satisfied.

Arguments

  • p: the plate
  • v₁, v₂: the vortex elements (with unit circulation) that the vorticity flux is going into
  • lesp, tesp: the critical leading and trailing edge suction parameters we want to enforce. By default, both parameters are set to 0.0 to enforce the Kutta condition on both edges. We can disable vortex shedding from an edge by setting the its critical suction parameter to Inf

Returns

  • Γ₁, Γ₂: the strengths that the vortex element should have in order to satisfy the edge suction parameters
  • ∂C₁, ∂C₂: Chebyshev coefficients of the normal velocity induced by the vortex elements Instead of running enforce_bc! with the new vortex elements, we can use this matrix to directly update the Chebyshev coefficients associated with the bound vortex sheet without recomputing all the velocities.

Example

Enforcing the trailing edge Kutta condition with an point vortex at negative infinity:

julia> plate = Plate(128, 2.0, 0.0, π/6)
Plate: N = 128, L = 2.0, c = 0.0 + 0.0im, α = 30.0ᵒ
       LESP = 0.0, TESP = 0.0

julia> motion = Plates.RigidBodyMotion(1.0, 0.0);

julia> Plates.enforce_no_flow_through!(plate, motion, (), 0.0)

julia> point = Vortex.Point(-Inf, 1.0);

julia> _, Γ, _, _ = Plates.vorticity_flux(plate, (), point, 0.0, Inf);

julia> Γ # should equal -πULsin(α) = -π
-3.1415926535897927
source
PotentialFlow.Plates.vorticity_flux!Function
vorticity_flux!(p::Plate, v₁, v₂,
                lesp = 0.0, tesp = 0.0,
                ∂C₁ = Vector{ComplexF64}(undef,plate.N),
                ∂C₂ = Vector{ComplexF64}(undef,plate.N))

In-place version of vorticity_flux, except instead of just returning the possible changes in plate Chebyshev coefficients, we modify plate.C with those changes so that no-flow-through is enforced in the presence of v₁ and v₂ with strengths that satisfy the suction parameters.

source
PotentialFlow.Plates.bound_circulationFunction
bound_circulation(plate[, s])

Compute the bound circulation between the trailing edge of the plate to s.

s can be either a single normalized arc length coordinate (between -1 and 1), or a whole array of coordinates.

source
PotentialFlow.Plates.bound_circulation!Function
bound_circulation!(Γs, plate[, ss])

Compute the bound circulation between the trailing edge of the plate to ss, then store it in Γs.

If an array, ss, with normalized arc length coordinates is omitted, then the circulation will be computed at the plate's Chebyshev nodes.

source
PotentialFlow.Plates.rate_of_impulseFunction
rate_of_impulse(plate, motion, elements::Source, velocities::Source)

Compute the rate of change of impulse of a vortex element and its image relative to a plate.

Note that this is not just the rate of impulse of the vortex element itself, but also includes the rate of impulse of the bound vortex sheet generated in response to the vortex element.

source
PotentialFlow.Plates.forceFunction
force(plate, motion, elements, velocities, newelements = ())

Compute the force on plate, given its motion and the state of the ambient vorticity.

Arguments

  • plate: the plate
  • motion: a structure that contains the velocity, acceleration, and angular velocity of the plate.
  • elements: vortex elements representing the ambient vorticity
  • velocities: the velocities of the vortex elements
  • newelements: an optional argument listing vortex elements that are just added to the flow field (it can be an element that is contained in elements)
  • Δt: this is only required if newelements is not empty, we assume that the new vortex elements are created over the span of Δt

Returns

  • F: the force exerted on the plate in complex coordinates
source
PotentialFlow.Plates.surface_pressureFunction
surface_pressure(plate, motion, te_sys, Γs₋, Δt)

Compute the pressure difference across the plate along Chebyshev nodes.

Note

The pressure difference across the bound vortex sheet is given by:

\[ [p]_-^+ = -\rho \left[ \frac{1}{2}(\boldsymbol{v}^+ + \boldsymbol{v}^-) - \boldsymbol{v}_b \right] \cdot ( \boldsymbol{\gamma} \times \boldsymbol{\hat{n}}) +\rho \frac{\mathrm{d}\Gamma}{\mathrm{d}t}\]

where $\rho$ is the fluid density, $\boldsymbol{v}^\pm$ is the velocity on either side of the plate, $\boldsymbol{v}_b$ is the local velocity of the plate, $\boldsymbol{\gamma}$ is the bound vortex sheet strength, and $\Gamma$ is the integrated circulation. We will compute $\frac{\mathrm{d}\Gamma}{\mathrm{d}t}$ using finite differences. So we will need the circulation along the plate from a previous time-step in order to compute the current pressure distribution. We assume that value of circulation at the trailing edge of the plate is equal the the net circulation of all the vorticity that has been shed from the trailing edge.

Arguments

  • plate: we assume that the Plate structure that is passed in already enforces the no-flow-through condition
  • motion: the motion of the plate used to compute $\boldsymbol{v}_b$
  • te_sys: the system of vortex elements representing the vorticity shed from the trailing edge of the plate
  • Γs₋: the circulation along the plate's Chebyshev nodes, this should be equivalent to calling Vortex.circulation(te_sys) .+ Vortex.bound_circulation(plate) from a previous time-step.
  • Δt: time-step used to compute $\frac{\mathrm{d}\Gamma}{\mathrm{d}t}$ using finite differences

Returns

  • Δp: the pressure difference across the plate along Chebyshev nodes
  • Γs₊: the circulation along the plate at the current time-step (this value is used in computing the current Δp and can be used as the Γs₋ for computing pressure differences at the next time-step)
source

Methods on Conformally-Mapped Bodies

PotentialFlow.Bodies.enforce_no_flow_through!Function
enforce_no_flow_through!(b::ConformalBody, motion, elements, t)

Update the body, b, to enforce the no-flow-through condition given ambient vortex elements, elements, and while moving with kinematics specified by motion.

Example

julia> p = Bodies.Polygon([-1.0,0.2,1.0,-1.0],[-1.0,-1.0,0.5,1.0])
Polygon with 4 vertices at
             (-1.0,-1.0) (0.2,-1.0) (1.0,0.5) (-1.0,1.0)
             interior angles/π = [0.5, 0.656, 0.422, 0.422]

julia> b = Bodies.ConformalBody(p)
Body generated by: Schwarz-Christoffel map of unit circle to exterior of polygon with 4 vertices

  centroid at 0.0 + 0.0im
  angle 0.0

julia> motion = RigidBodyMotion(1.0, 0.0);

julia> point = Vortex.Point(0.0 + 2im, 1.0);

julia> Bodies.enforce_no_flow_through!(b, motion, point, 0.0)

julia> b.img
1-element Array{Element,1}:
 Vortex.Point(0.0 + 0.5im, -1.0)
source
PotentialFlow.Bodies.normalFunction
normal(ζ,v,b::ConformalBody)

Returns the normal component of the complex vector(s) v in the physical plane at a point(s) on the surface of body b. Each surface point is specified by its pre-image ζ on the unit circle. v and ζ can be arrays of points.

Example

julia> p = Bodies.Polygon([-1.0,1.0,1.0,-1.0],[-1.0,-1.0,1.0,1.0]);

julia> b = Bodies.ConformalBody(p);

julia> Bodies.normal(exp(im*0),exp(im*π/4),b)
0.7071067811865472
source
PotentialFlow.Bodies.tangentFunction
tangent(ζ,v,b::ConformalBody)

Returns the (counter-clockwise) tangent component of the complex vector(s) v in the physical plane at a point(s) on the surface of body b. Each surface point is specified by its pre-image ζ on the unit circle. v and ζ can be arrays of points.

Example

julia> p = Bodies.Polygon([-1.0,1.0,1.0,-1.0],[-1.0,-1.0,1.0,1.0]);

julia> b = Bodies.ConformalBody(p);

julia> Bodies.tangent(exp(im*0),exp(im*π/4),b)
0.7071067811865478
source
PotentialFlow.Bodies.transform_velocity!Function
transform_velocity!(wout, win, targets, body::ConformalBody)

Transforms the velocity win in the circle plane of a conformal mapping to a velocity wout that can actually be used to transport the pre-images of elements in targets in this circle plane. This transformation applies the Routh correction and subtracts the relative motion of the body.

Example

julia> a1 = 1; b1 = 0.1; ccoeff = ComplexF64[0.5(a1+b1),0,0.5(a1-b1)];

julia> body = Bodies.ConformalBody(ccoeff);

julia> motion = RigidBodyMotion(0,0);

julia> points = Vortex.Point.([-2, 2], 1.0);

julia> Bodies.enforce_no_flow_through!(body, motion, points, 0);

julia> sys = (body,points);

julia> ẋ = (motion, allocate_velocity(points));

julia> self_induce_velocity!(ẋ, sys, 0)
(Rigid Body Motion:
  ċ = 0.0 + 0.0im
  c̈ = 0.0 + 0.0im
  α̇ = 0.0
  α̈ = 0.0
  Constant (ċ = 0 + 0im, α̇ = 0), Complex{Float64}[0.0+0.129977im, 0.0-0.129977im])

julia> Bodies.transform_velocity!(ẋ, ẋ, sys, body)
(Rigid Body Motion:
  ċ = 0.0 + 0.0im
  c̈ = 0.0 + 0.0im
  α̇ = 0.0
  α̈ = 0.0
  Constant (ċ = 0 + 0im, α̇ = 0), Complex{Float64}[0.0+0.785969im, 0.0-0.785969im])
transform_velocity(win, target::ComplexF64, body::ConformalBody)

Returns the velocity in the physical plane from the velocity win in the circle plane.

julia> a1 = 1; b1 = 0.1; ccoeff = ComplexF64[0.5(a1+b1),0,0.5(a1-b1)];

julia> body = Bodies.ConformalBody(ccoeff,0.0+0.0im,π/4);

julia> motion = RigidBodyMotion(1,0);

julia> points = Vortex.Point.([-2, 2], 1.0);

julia> Bodies.enforce_no_flow_through!(body, motion, points, 0);

julia> sys = (body,points);

julia> ζ = exp(-im*π/4);

julia> w̃ = induce_velocity(ζ,sys,0);

julia> w = Bodies.transform_velocity(w̃,ζ,body)
0.7497272298496697 - 0.3058889412948484im
source

Index