Surface-grid operations

Here, we will discuss the various surface-grid operators available in the package. We will start by generating the cache, just as we did in Immersed layer caches

using ImmersedLayers
using Plots
using LinearAlgebra

Set up the grid, shape, and cache

We do this just as we did in Immersed layer caches

Δx = 0.01
Lx = 4.0
xlim = (-Lx/2,Lx/2)
ylim = (-Lx/2,Lx/2)
g = PhysicalGrid(xlim,ylim,Δx)
RadC = 1.0
Δs = 1.4*cellsize(g)
body = Circle(RadC,Δs)
cache = SurfaceScalarCache(body,g,scaling=GridScaling)
Surface cache with scaling of type GridScaling
  448 point data of type ScalarData{448, Float64, Vector{Float64}}
  Grid data of type Nodes{Primal, 405, 405, Float64, Matrix{Float64}}

Basic regularization and interpolation

Let's perform an example in which we regularize the x coordinate of the surface points onto the grid. We can get the surface coordinates by using the points function. This can be applied either to body directly or to the cache. As for any VectorData type, the components of pts are pts.u and pts.v.

pts = points(cache);

Now, we set up some blank grid data on which to regularize onto

gx = zeros_grid(cache);

Now regularize. We will time it to show that it is fast and memory-efficient:

@time regularize!(gx,pts.u,cache);
  0.001504 seconds (35.97 k allocations: 562.078 KiB)

Let's plot this to look at it. This also gives us a chance to highlight the plot recipe for grid data associated with the cache, which is achieved by simply supplying the cache to the plot function in Plots.jl. By default, this plots the immersed points, as well, but this can be suppressed by adding the keyword layers=false.

plot(gx,cache)
Example block output

This shows how the regularization spreads the data over a couple of cells around the surface. In the parlance of potential theory, this is a single layer. If we wish to interpolate, then we do so with the [interpolate!] function. For example, suppose we have a uniform field that we wish to interpolate:

oc = 2.5*ones_grid(cache)
Nodes{Primal, 405, 405, Float64, Matrix{Float64}} data
Printing in grid orientation (lower left is (1,1))
404×404 Matrix{Float64}:
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5  …  2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5  …  2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 ⋮                        ⋮              ⋱                 ⋮              
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5  …  2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5  …  2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5
 2.5  2.5  2.5  2.5  2.5  2.5  2.5  2.5     2.5  2.5  2.5  2.5  2.5  2.5  2.5

Now set up some surface data to receive the interpolated data, and interpolate:

f = zeros_surface(cache);
interpolate!(f,oc,cache)
448 points of scalar-valued Float64 data
448-element Vector{Float64}:
 2.5000000000000004
 2.5
 2.4999999999999982
 2.4999999999999987
 2.5
 2.499999999999998
 2.4999999999999987
 2.4999999999999982
 2.4999999999999987
 2.4999999999999987
 ⋮
 2.499999999999999
 2.4999999999999996
 2.4999999999999987
 2.4999999999999987
 2.4999999999999982
 2.499999999999999
 2.499999999999999
 2.499999999999999
 2.4999999999999996

It is clear that the interpolation preserves the value of the field. This is also true for linearly varying fields, since the DDF is built to ensure this. Let's try this, using the $x$ coordinate of the grid. Here, we use the coordinates function of CartesianGrids.jl, which gets the coordinates of the grid, and set the values of grid data to the $x$ coordinate. We interpolate and plot, comparing to the actual $x$ coordinate of the points on the body:

xg = x_grid(cache)
interpolate!(f,xg,cache)
plot(f,cache,ylim=(-2,2),label="Interpolated from grid",ylabel="x",xlabel="arclength")
plot!(pts.u,cache,label="Actual body coordinate",legend=:true)
Example block output

A double layer

Now we will generate a double layer. Mathematically, this takes the form

\[D_s f = \nabla\cdot \left( \delta(\chi) \mathbf{n} f \right)\]

for some scalar data $f$ on the surface. (See Background for an example.) Notice that it maps scalar data on the surface ($f$) to scalar data in space. So to calculate this using our discrete tools, we set up some grid data to receive the result. Then, we use the function surface_divergence! to compute the double layer. Here, we will demonstrate this on the $y$ coordinate of the surface points:

dl = zeros_grid(cache)
surface_divergence!(dl,pts.v,cache)
plot(dl,cache)
Example block output

If the surface data are vectors, $\mathbf{f}$, then this operation is a little different:

\[D_s \mathbf{f} = \nabla\cdot \left[\delta(\chi) \left( \mathbf{n} \mathbf{f} + \mathbf{f} \mathbf{n} \right) \right]\]

This maps $\mathbf{f}$ to a vector field. We use this in conjunction with a cache generated with SurfaceVectorCache.

The transpose of the double layer, $D_s$, is the operation

\[G_s u = \mathbf{n}\cdot \delta^{T}(\chi) \nabla u\]

for some scalar field data $u$. This operation computes the gradient of the field data, interpolates this gradient onto the surface, and obtains the normal component of that surface data. As such, it represents an important tool for computing the normal derivative of field data on a surface. In the package, we use surface_grad! for this operation.

The vector field version of this is

\[G_s \mathbf{u} = \mathbf{n}\cdot \delta^{T}(\chi) (\nabla \mathbf{u} + \nabla^{T} \mathbf{u})\]

which maps vector field data $\mathbf{u}$ to vector-valued surface data.

A curl layer

We also sometimes need to take the curl of the regularized surface data,

\[C_s f = \nabla\times \left( \delta(\chi) \mathbf{n} f \right)\]

For this, we use the surface_curl! operator. Let's demonstrate this on a uniform field on the surface.

gc = zeros_gridcurl(cache)
f = ones_surface(cache)
surface_curl!(gc,f,cache)
plot(gc,cache)
Example block output

The continuous version of this operation is actually zero. It's not quite zero in discrete form. However, its norm is much smaller than that of the double layer.

norm(gc,cache)/norm(dl,cache)
0.056709335671267

Finally, a pair of operations that are used in support of the previous ones, or occasionally on their own, are

\[R_n f = \delta(\chi)\mathbf{n}\circ f\]

for scalar surface data $f$, which maps to a vector field, effectively a field of doublet strengths; and its transpose

\[R_n^T \mathbf{u} = \mathbf{n}\cdot \delta^{T}(\chi)\mathbf{u}\]

which maps vector field data $\mathbf{u}$ to a scalar surface field, the normal component of the vector field on the surface. These are provided by regularize_normal! and normal_interpolate!, respectively.

Masks

Masks are grid data that take the value 1 in one region (e.g., the interior of a surface) and 0 in the other (e.g., the exterior). The functions mask and complementary_mask achieve this

m = mask(cache)
cm = complementary_mask(cache)
plot(
   surface(m,cache,layers=false),
   surface(cm,cache,layers=false)
   )
Example block output

One can apply a mask to some grid data by multiplying it, using, e.g., the product! function in CartesianGrids.jl. Let's demonstrate that with the grid data of $x$ coordinates:

xmask = zeros_grid(cache)
xcmask = zeros_grid(cache)
product!(xmask,xg,m)
product!(xcmask,xg,cm)
plot(
 plot(xmask,cache),
 plot(xcmask,cache)
 )
Example block output

The mask and complementary mask effectively partition the field into two parts. We can also apply masks in place, using mask! and complementary_mask!:

xmask .= xg
mask!(xmask,cache)
plot(xmask,cache)
Example block output

Surface-grid operator functions

ImmersedLayers.regularize!Function
regularize!(s::Nodes{Primal},f::ScalarData,cache::BasicILMCache)
regularize!(s::Nodes{Primal},f::ScalarData,sys::ILMSystem)

The operation $s = R_c f$, which regularizes scalar surface data f onto the grid in the form of scalar grid data s. This is the adjoint to interpolate!

source
regularize!(s::Nodes{Dual},f::ScalarData,cache::BasicILMCache)
regularize!(s::Nodes{Dual},f::ScalarData,sys::ILMSystem)

The operation $s = R_N f$, which regularizes scalar surface data f onto the grid in the form of scalar grid (curl) data s. Only works if the cache is built for vector data. This is the adjoint to interpolate!

source
regularize!(v::Edges{Primal},vb::VectorData,cache::BasicILMCache)
regularize!(v::Edges{Primal},vb::VectorData,sys::ILMSystem)

The operation $\mathbf{v} = R_f \mathbf{v}_b$, which regularizes vector surface data vb onto the grid in the form of scalar grid data v. This is the adjoint to interpolate!

source
ImmersedLayers.interpolate!Function
interpolate!(f::ScalarData,s::Nodes{Primal},cache::BasicILMCache)
interpolate!(f::ScalarData,s::Nodes{Primal},sys::ILMSystem)

The operation $f = R_c^T s$, which interpolates scalar grid data s onto the surface points in the form of scalar point data f. This is the adjoint to regularize!

source
interpolate!(f::ScalarData,s::Nodes{Dual},cache::BasicILMCache)
interpolate!(f::ScalarData,s::Nodes{Dual},sys::ILMSystem)

The operation $f = R_N^T s$, which interpolates scalar grid (curl) data s onto the surface points in the form of scalar point data f. Only works if the cache is built for vector data. This is the adjoint to regularize!

source
interpolate!(vb::VectorData,v::Edges{Primal},cache::BasicILMCache)
interpolate!(vb::VectorData,v::Edges{Primal},sys::ILMSystem)

The operation $\mathbf{v}_b = R_c^T \mathbf{v}$, which interpolates vector grid data v onto the surface points in the form of scalar point data vb. This is the adjoint to regularize!

source
ImmersedLayers.regularize_normal!Function
regularize_normal!(v::Edges{Primal},f::ScalarData,cache::BasicILMCache)
regularize_normal!(v::Edges{Primal},f::ScalarData,sys::ILMSystem)

The operation $\mathbf{v} = R_f \mathbf{n}\circ f$, which maps scalar surface data f (like a jump in scalar potential) to grid data v (like velocity). This is the adjoint to normal_interpolate!.

source
regularize_normal!(qt::EdgeGradient{Primal},v::VectorData,cache::BasicILMCache)
regularize_normal!(qt::EdgeGradient{Primal},v::VectorData,sys::ILMSystem)

The operation $\mathbf{q}_t = R_t \mathbf{n}\circ \mathbf{v}$, which maps scalar vector data v (like a jump in velocity) to grid data qt (like velocity-normal tensor). This is the adjoint to normal_interpolate!.

source
ImmersedLayers.regularize_normal_symm!Function
regularize_normal_symm!(qt::EdgeGradient{Primal},v::VectorData,cache::BasicILMCache)
regularize_normal_symm!(qt::EdgeGradient{Primal},v::VectorData,sys::ILMSystem)

The operation $\mathbf{q}_t = R_t (\mathbf{n}\circ \mathbf{v}+\mathbf{v}\circ \mathbf{n})$, which maps scalar vector data v (like a jump in velocity) to grid data qt (like velocity-normal tensor). This is the adjoint to normal_interpolate_symm!.

source
ImmersedLayers.normal_interpolate!Function
normal_interpolate!(vn::ScalarData,v::Edges{Primal},cache::BasicILMCache)
normal_interpolate!(vn::ScalarData,v::Edges{Primal},sys::ILMSystem)

The operation $v_n = \mathbf{n} \cdot R_f^T \mathbf{v}$, which maps grid data v (like velocity) to scalar surface data vn (like normal component of surface velocity). This is the adjoint to regularize_normal!.

source
normal_interpolate!(τ::VectorData,A::EdgeGradient{Primal},cache::BasicILMCache)
normal_interpolate!(τ::VectorData,A::EdgeGradient{Primal},sys::ILMSystem)

The operation $\mathbf{\tau} = \mathbf{n} \cdot R_t^T \mathbf{A}$, which maps grid tensor data A (like velocity gradient tensor) to vector surface data τ (like traction). This is the adjoint to regularize_normal!.

source
ImmersedLayers.normal_interpolate_symm!Function
normal_interpolate_symm!(τ::VectorData,A::EdgeGradient{Primal},cache::BasicILMCache)
normal_interpolate_symm!(τ::VectorData,A::EdgeGradient{Primal},sys::ILMSystem)

The operation $\mathbf{\tau} = \mathbf{n} \cdot R_t^T (\mathbf{A} +\mathbf{A}^T - (\mathrm{tr}\mathbf{A})\mathbf{I})$, which maps grid tensor data A (like velocity gradient tensor) to vector surface data τ (like traction). This is the adjoint to regularize_normal_symm!.

source
ImmersedLayers.regularize_normal_cross!Function
regularize_normal_cross!(v::Edges{Primal},f::ScalarData,cache::BasicILMCache)
regularize_normal_cross!(v::Edges{Primal},f::ScalarData,sys::ILMSystem)

The operation $\mathbf{v} = R_f \mathbf{n}\times f \mathbf{e}_z$, which maps scalar surface data f (like a jump in streamfunction, endowed with the out-of-plane unit vector $\mathbf{e}_z$) to grid data v (like velocity). This is the negative adjoint to normal_cross_interpolate!.

source
regularize_normal_cross!(w::Nodes{Dual},vs::VectorData,cache::BasicILMCache)
regularize_normal_cross!(w::Nodes{Dual},vs::VectorData,sys::ILMSystem)

The operation $\omega = R_N \mathbf{n}\times \mathbf{v}_s$, which maps surface vector data vs (like a jump in velocity) to grid dual nodal data w (like vorticity). This is for use only with vector-data caches. This is the negative adjoint to normal_cross_interpolate!.

source
ImmersedLayers.normal_cross_interpolate!Function
normal_cross_interpolate!(wn::ScalarData,v::Edges{Primal},cache::BasicILMCache)
normal_cross_interpolate!(wn::ScalarData,v::Edges{Primal},sys::ILMSystem)

The operation $w_n = e_z\cdot (\mathbf{n} \times R_f^T \mathbf{v})$, which maps grid data v (like velocity) to scalar surface data wn (like vorticity in the surface). This is the negative adjoint to regularize_normal_cross!.

source
normal_cross_interpolate!(vs::VectorData,s::Nodes{Dual},cache::BasicILMCache)
normal_cross_interpolate!(vs::VectorData,s::Nodes{Dual},sys::ILMSystem)

The operation $\mathbf{v}_s = \mathbf{n}\times R_N^T \psi\mathbf{e}_z)$, which maps grid data s (like streamfunction) to vector surface data vs (like velocity in the surface). It only works for a vector-type cache. This is the negative adjoint to regularize_normal_cross!.

source
ImmersedLayers.regularize_normal_dot!Function
regularize_normal_dot!(f::Nodes{Primal},vs::VectorData,cache::BasicILMCache)
regularize_normal_dot!(f::Nodes{Primal},vs::VectorData,sys::ILMSystem)

The operation $\phi = R_c \mathbf{n}\cdot \mathbf{v}_s$, which maps surface vector data vs (like a jump in velocity) to grid nodal data f (like scalar potential). This is the negative adjoint to normal_dot_interpolate!.

source
regularize_normal_dot!(v::Edges{Primal},taus::TensorData,cache::BasicILMCache)
regularize_normal_dot!(v::Edges{Primal},taus::TensorData,sys::ILMSystem)

The operation $\mathbf{v} = R_f \mathbf{n}\cdot \mathbf{\tau}_s$, which maps surface tensor data taus (like a stress) to grid edge data v (like velocity). This is the negative adjoint to normal_dot_interpolate!.

source
ImmersedLayers.normal_dot_interpolate!Function
normal_dot_interpolate!(vs::VectorData,f::Nodes{Primal},cache::BasicILMCache)
normal_dot_interpolate!(vs::VectorData,f::Nodes{Primal},sys::ILMSystem)

The operation $\mathbf{v}_s = \mathbf{n} R_c^T \phi$, which maps grid nodal data f (like scalar potential) to surface vector data vs (like velocity). This is the negative adjoint to regularize_normal_dot!.

source
normal_dot_interpolate!(taus::TensorData,v::Edges{Primal},cache::BasicILMCache)
normal_dot_interpolate!(taus::TensorData,v::Edges{Primal},sys::ILMSystem)

The operation $\mathbf{\tau}_s = \mathbf{n}\circ R_f^T\mathbf{v}$, which maps grid edge data v (like velocity) to surface tensor data taus (like stress). This is the negative adjoint to regularize_normal_dot!.

source
ImmersedLayers.surface_divergence!Function
surface_divergence!(Θ::Nodes{Primal},f::ScalarData,cache::BasicILMCache)
surface_divergence!(Θ::Nodes{Primal},f::ScalarData,sys::ILMSystem)

The operation $\theta = D_s f = D R_f \mathbf{n} \circ f$, which maps surface scalar data f (like jump in scalar potential) to grid data Θ (like dilatation, i.e. divergence of velocity). This is the negative adjoint of surface_grad!. Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
surface_divergence!(v::Edges{Primal},dv::VectorData,cache::BasicILMCache)
surface_divergence!(v::Edges{Primal},dv::VectorData,sys::ILMSystem)

The operation $\mathbf{v} = D_s d\mathbf{v} = D R_t (\mathbf{n} \circ d\mathbf{v})$, which maps surface vector data dv (like jump in velocity) to grid data v (like velocity). This is the negative adjoint of surface_grad!. Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
ImmersedLayers.surface_grad!Function
surface_grad!(vn::ScalarData,ϕ::Nodes{Primal},cache::BasicILMCache)
surface_grad!(vn::ScalarData,ϕ::Nodes{Primal},sys::ILMSystem)

The operation $v_n = G_s\phi = \mathbf{n} \cdot R_f^T G\phi$, which maps grid data ϕ (like scalar potential) to scalar surface data vn (like normal component of velocity). This is the negative adjoint of surface_divergence!. Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
surface_grad!(τ::VectorData,v::Edges{Primal},cache::BasicILMCache)
surface_grad!(τ::VectorData,v::Edges{Primal},sys::ILMSystem)

The operation $\mathbf{\tau} = G_s v = \mathbf{n} \cdot R_t^T G \mathbf{v}$, which maps grid vector data v (like velocity) to vector surface data τ (like traction). This is the negative adjoint of surface_divergence!. Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
ImmersedLayers.surface_curl!Function
surface_curl!(w::Nodes{Dual},f::ScalarData,cache::BasicILMCache)
surface_curl!(w::Nodes{Dual},f::ScalarData,sys::ILMSystem)

The operation $w = C_s^T f = C^T R_f \mathbf{n}\circ f$, which maps scalar surface data f (like a jump in scalar potential) to grid data w (like vorticity). This is the adjoint to $C_s$, also given by surface_curl! (but with arguments switched). Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
surface_curl!(w::Nodes{Dual},v::VectorData,cache::BasicILMCache)
surface_curl!(w::Nodes{Dual},v::VectorData,sys::ILMSystem)

The operation $w = C_s^T \mathbf{v} = C^T R_f \mathbf{v}$, which maps vector surface data v (like velocity) to grid data w (like vorticity). This is the adjoint to $C_s$, also given by surface_curl! (but with arguments switched). Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
surface_curl!(vn::ScalarData,s::Nodes{Dual},cache::BasicILMCache)
surface_curl!(vn::ScalarData,s::Nodes{Dual},sys::ILMSystem)

The operation $v_n = C_s s = \mathbf{n} \cdot R_f^T C s$, which maps grid data s (like streamfunction) to scalar surface data vn (like normal component of velocity). This is the adjoint to $C_s^T$, also given by surface_curl!, but with arguments switched. Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
surface_curl!(v::VectorData,s::Nodes{Dual},cache::BasicILMCache)
surface_curl!(v::VectorData,s::Nodes{Dual},sys::ILMSystem)

The operation $\mathbf{v} = C_s s = R_f^T C s$, which maps grid data s (like streamfunction) to vector surface data v (like velocity). This is the adjoint to $C_s^T$, also given by surface_curl!, but with arguments switched. Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
ImmersedLayers.surface_divergence_symm!Function
surface_divergence_symm!(v::Edges{Primal},dv::VectorData,cache::BasicILMCache)
surface_divergence_symm!(v::Edges{Primal},dv::VectorData,sys::ILMSystem)

The operation $\mathbf{v} = D_s d\mathbf{v} = D R_t (\mathbf{n} \circ d\mathbf{v} + d\mathbf{v} \circ \mathbf{n})$, which maps surface vector data dv (like jump in velocity) to grid data v (like velocity). This is the negative adjoint of surface_grad_symm!. Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
ImmersedLayers.surface_divergence_cross!Function
surface_divergence_cross!(Θ::Nodes{Primal},f::ScalarData,cache::BasicILMCache)
surface_divergence_cross!(Θ::Nodes{Primal},f::ScalarData,sys::ILMSystem)

The operation $\theta = \hat{D}_s f = D R_f \mathbf{n} \times f \mathbf{e}_z$, which maps surface scalar data f (like jump in streamfunction) to grid data Θ (like dilatation, i.e. divergence of velocity). This is the adjoint of surface_grad_cross!. Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
ImmersedLayers.surface_grad_symm!Function
surface_grad_symm!(τ::VectorData,v::Edges{Primal},cache::BasicILMCache)
surface_grad_symm!(τ::VectorData,v::Edges{Primal},sys::ILMSystem)

The operation $\mathbf{\tau} = G_s v = \mathbf{n} \cdot R_t^T (G \mathbf{v} + (G \mathbf{v})^T)$, which maps grid vector data v (like velocity) to vector surface data τ (like traction). This is the negative adjoint of surface_divergence_symm!. Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
ImmersedLayers.surface_grad_cross!Function
surface_grad_cross!(γ::ScalarData,ϕ::Nodes{Primal},cache::BasicILMCache)
surface_grad_cross!(γ::ScalarData,ϕ::Nodes{Primal},sys::ILMSystem)

The operation $\gamma = \hat{G}_s\phi = \mathbf{e}_z \cdot (\mathbf{n} \times R_f^T G\phi)$, which maps grid data ϕ (like scalar potential) to scalar surface data γ (like surface vorticity). This is the adjoint of surface_divergence_cross!. Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
ImmersedLayers.surface_curl_cross!Function
surface_curl_cross!(w::Nodes{Dual},f::ScalarData,cache::BasicILMCache)
surface_curl_cross!(w::Nodes{Dual},f::ScalarData,sys::ILMSystem)

The operation $w = \hat{C}_s^T f = C^T R_f \mathbf{n}\times f \mathbf{e}_z$, which maps scalar surface data f (like a jump in streamfunction, multiplied by out-of-plane unit vector $\mathbf{e}_z$) to grid data w (like vorticity). This is the adjoint to $\hat{C}_s$, also given by surface_curl_cross! (but with arguments switched). Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
surface_curl_cross!(γ::ScalarData,s::Nodes{Dual},cache::BasicILMCache)
surface_curl_cross!(γ::ScalarData,s::Nodes{Dual},sys::ILMSystem)

The operation $\gamma = \hat{C}_s s = \mathbf{e}_z\cdot (\mathbf{n} \times R_f^T C s)$, which maps grid data s (like streamfunction) to scalar surface data γ (like vorticity in the surface). This is the adjoint to $\hat{C}_s^T$, also given by surface_curl_cross!, but with arguments switched. Note that the differential operations are divided either by 1 or by the grid cell size, depending on whether sys has been designated with IndexScaling or GridScaling, respectively.

source
ImmersedLayers.mask!Function
mask!(w::GridData,cache::BasicILMCache)
mask!(w::GridData,sys::ILMSystem)

Mask the data w in place by multiplying it by 1s inside of a surface (i.e., on a side opposite the normal vectors) and 0s outside. The grid data w must be of the same type as the output data type of cache. Only allows cache to have GridScaling.

source
mask!(w::GridData,surface::Body,grid::PhysicalGrid)

Mask the data w in place by multiplying it by 1s inside of surface surface (i.e., on a side opposite the normal vectors) and 0s outside.

source
ImmersedLayers.maskFunction
mask(cache::BasicILMCache) -> GridData
mask(sys::ILMSystem) -> GridData

Create grid data that consist of 1s inside of a surface (i.e., on a side opposite the normal vectors) and 0s outside. The grid data are the same type as the output data type of sys. Only allows sys to have GridScaling. If the cache has no immersed surface, then it reverts to 1s everywhere.

source
mask(w::GridData,surface::Body,grid::PhysicalGrid)

Create grid data that consist of 1s inside of surface surface (i.e., on a side opposite the normal vectors) and 0s outside. The grid data are the same type as w.

source
mask(ar::AreaRegionCache)

Return the mask for the given area region ar.

source
ImmersedLayers.complementary_mask!Function
complementary_mask!(w::GridData,cache::BasicILMCache)
complementary_mask!(w::GridData,sys::ILMSystem)

Mask the data w in place by multiplying it by 0s inside of a surface (i.e., on a side opposite the normal vectors) and 1s outside. The grid data w must be of the same type as the output data type of cache. Only allows cache to have GridScaling.

source
complementary_mask!(w::GridData,surface::Body,grid::PhysicalGrid)

Mask the data w in place by multiplying it by 0s inside of surface surface (i.e., on a side opposite the normal vectors) and 1s outside.

source
ImmersedLayers.complementary_maskFunction
complementary_mask(cache::BasicILMCache) -> GridData
complementary_mask(sys::ILMSystem) -> GridData

Create grid data that consist of 0s inside of a surface (i.e., on a side opposite the normal vectors) and 1s outside. The grid data are the same type as the output data type of cache. Only allows cache to have GridScaling. If the cache has no immersed surface, then it reverts to 0s everywhere.

source
complementary_mask(w::GridData,surface::Body,grid::PhysicalGrid)

Create grid data that consist of 0s inside of surface surface (i.e., on a side opposite the normal vectors) and 1s outside. The grid data are the same type as w.

source

This page was generated using Literate.jl.