Immersed layer caches
This package uses precomputed caches to efficiently implement the immersed layer operators. The starting point for a cache is the specification of the (discretized) body shape and the grid. Let's use an example.
Setting up a cache
using ImmersedLayers
using Plots
using LinearAlgebra
Set up a grid
First, we will set up a grid for performing the operations. We use the PhysicalGrid
constructor of the CartesianGrids.jl
package to create this.
Δx = 0.01
Lx = 4.0
xlim = (-Lx/2,Lx/2)
ylim = (-Lx/2,Lx/2)
g = PhysicalGrid(xlim,ylim,Δx)
PhysicalGrid{2}((405, 405), (202, 202), 0.01, ((-2.0100000000000002, 2.02), (-2.0100000000000002, 2.02)), 1)
Set the shape
Now let's set a shape to immerse into the grid. We will use a circle, but there are a variety of other shapes available. Many of these are in the RigidBodyTools.jl
package. Note that we set the spacing between the points on this shape equal to 1.4 times the grid spacing. This is not critical, but it is generally best to set it to a value between 1 and 2.
RadC = 1.0
Δs = 1.4*cellsize(g)
body = Circle(RadC,Δs)
Circular body with 448 points and radius 1.0
Current position: (0.0,0.0)
Current angle (rad): 0.0
Create the cache
After setting the grid and the surface to immerse, the next step for using the immersed layer tools is to set up a surface cache. This allocates a set of data structures, as well as the critical regularization and interpolation operators that will get used.
There are a few choices to make when setting this up
- What kind of data (scalar or vector) are we dealing with?
We will demonstrate with scalar data, which means we use the SurfaceScalarCache
function. For vector data, use SurfaceVectorCache
.
- What type of scaling (grid or index) do we wish to apply to the operators?
Grid scaling, set with scaling = GridScaling
, means that the various operators are scaled with the physical grid spacing and/or surface point spacing so that they approximate the continuous operators. This means that regularization and interpolation are transposes with respect to inner products that incorporate these physical spacings, rather than the usual linear algebra inner products. Also, differential operations on the grid are true approximations of their continuous counterparts. This choice of scaling is usually the best, and the dot
operator is extended in this package to implement the physically- scaled inner products.
On the other hand, scaling = IndexScaling
does not scale these, but rather, uses pure differencing for the grid differential operators, and regularization is the simple matrix transpose of interpolation.
- What discrete Dirac delta function (DDF) do we wish to use?
This is specified with the ddftype
keyword argument. The default is ddftype=CartesianGrids.Yang3
[1]. However, there are other choices, such as CartesianGrids.Roma
, CartesianGrids.Goza
, CartesianGrids.Witchhat
, CartesianGrids.M3
, CartesianGrids.M4prime
.
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}}
Plotting the immersed points
We can plot the immersed points with the plot
function of the Plots.jl
package, using a special plot recipe. This accepts all of the usual attribute keywords of the basic plot function.
plot(cache,xlims=(-2,2),ylims=(-2,2))
In this plotting demonstration, we have used the points
function to obtain the coordinates of the immersed points. Other useful utilities are normals
, areas
, and arcs
, e.g.
normals(cache)
448 points of vector-valued Float64 data
896-element Vector{Float64}:
0.9999999999988712
0.9999016728287639
0.9996065842578712
0.999115045044428
0.9984267309011199
0.9975423661130551
0.9964613682740765
0.9951848732366676
0.9937120427187925
0.992044421068211
⋮
-0.12589126758565716
-0.11196299239880425
-0.09801862843388753
-0.08404903300348561
-0.07006887224977693
-0.05606894942065292
-0.04206398840993138
-0.02804475491729984
-0.014026009668317694
Some basic utilities
We will see deeper uses of this cache in Surface-grid operations. However, for now we can learn how to get basic copies of the grid and surface data. For example, a copy of the grid data, all initialized to zero:
zeros_grid(cache)
Nodes{Primal, 405, 405, Float64, Matrix{Float64}} data
Printing in grid orientation (lower left is (1,1))
404×404 Matrix{Float64}:
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 … 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 … 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
⋮ ⋮ ⋱ ⋮
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 … 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 … 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
or similarly on the surface
zeros_surface(cache)
448 points of scalar-valued Float64 data
448-element Vector{Float64}:
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0
⋮
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0
0.0
We also might need to initialize data on the grid to accept the curl of a vector field. This is used in conjunction with surface_curl!
, for example.
zeros_gridcurl(cache)
Nodes{Dual, 405, 405, Float64, Matrix{Float64}} data
Printing in grid orientation (lower left is (1,1))
405×405 Matrix{Float64}:
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 … 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 … 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
⋮ ⋮ ⋱ ⋮
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 … 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
If we want them initialized to unity, then use
ones_grid(cache)
Nodes{Primal, 405, 405, Float64, Matrix{Float64}} data
Printing in grid orientation (lower left is (1,1))
404×404 Matrix{Float64}:
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 … 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 … 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
⋮ ⋮ ⋱ ⋮
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 … 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 … 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0
and
ones_surface(cache)
448 points of scalar-valued Float64 data
448-element Vector{Float64}:
1.0
1.0
1.0
1.0
1.0
1.0
1.0
1.0
1.0
1.0
⋮
1.0
1.0
1.0
1.0
1.0
1.0
1.0
1.0
1.0
To evaluate functions on the grid, it is useful to be able to fill grid data with the x and y coordinates. For this, we use
x_grid(cache)
Nodes{Primal, 405, 405, Float64, Matrix{Float64}} data
Printing in grid orientation (lower left is (1,1))
404×404 Matrix{Float64}:
-2.01 -2.0 -1.99 -1.98 -1.97 … 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 … 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
⋮ ⋱ ⋮
-2.01 -2.0 -1.99 -1.98 -1.97 … 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 … 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
-2.01 -2.0 -1.99 -1.98 -1.97 1.97 1.98 1.99 2.0 2.01 2.02
y_grid(cache)
Nodes{Primal, 405, 405, Float64, Matrix{Float64}} data
Printing in grid orientation (lower left is (1,1))
404×404 Matrix{Float64}:
2.02 2.02 2.02 2.02 2.02 … 2.02 2.02 2.02 2.02 2.02
2.01 2.01 2.01 2.01 2.01 2.01 2.01 2.01 2.01 2.01
2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0
1.99 1.99 1.99 1.99 1.99 1.99 1.99 1.99 1.99 1.99
1.98 1.98 1.98 1.98 1.98 1.98 1.98 1.98 1.98 1.98
1.97 1.97 1.97 1.97 1.97 … 1.97 1.97 1.97 1.97 1.97
1.96 1.96 1.96 1.96 1.96 1.96 1.96 1.96 1.96 1.96
1.95 1.95 1.95 1.95 1.95 1.95 1.95 1.95 1.95 1.95
1.94 1.94 1.94 1.94 1.94 1.94 1.94 1.94 1.94 1.94
1.93 1.93 1.93 1.93 1.93 1.93 1.93 1.93 1.93 1.93
⋮ ⋱ ⋮
-1.93 -1.93 -1.93 -1.93 -1.93 … -1.93 -1.93 -1.93 -1.93 -1.93
-1.94 -1.94 -1.94 -1.94 -1.94 -1.94 -1.94 -1.94 -1.94 -1.94
-1.95 -1.95 -1.95 -1.95 -1.95 -1.95 -1.95 -1.95 -1.95 -1.95
-1.96 -1.96 -1.96 -1.96 -1.96 -1.96 -1.96 -1.96 -1.96 -1.96
-1.97 -1.97 -1.97 -1.97 -1.97 -1.97 -1.97 -1.97 -1.97 -1.97
-1.98 -1.98 -1.98 -1.98 -1.98 … -1.98 -1.98 -1.98 -1.98 -1.98
-1.99 -1.99 -1.99 -1.99 -1.99 -1.99 -1.99 -1.99 -1.99 -1.99
-2.0 -2.0 -2.0 -2.0 -2.0 -2.0 -2.0 -2.0 -2.0 -2.0
-2.01 -2.01 -2.01 -2.01 -2.01 -2.01 -2.01 -2.01 -2.01 -2.01
Using norms, inner products, and integrals
It is useful to compute norms and inner products on grid or surface data. These tools are easily accessible and intuitive, e.g., dot(u,v,cache)
and norm(u,cache)
, and they respect the scaling associated with the cache. For example, the following gives an approximation of the circumference of the circle:
os = ones_surface(cache)
dot(os,os,cache)
6.283288300933757
We can also numerically integrate on surfaces or grids. For example, the previous example could have been written
integrate(os,cache)
6.283288300933757
Or, to compute an area integral over the whole domain,
og = ones_grid(cache)
integrate(og,cache)
16.240900000000003
If the underlying data are of vector or tensor form, then the integrate
function operates on every component
ovg = ones_gridgrad(cache)
integrate(ovg,cache)
2-element Vector{Float64}:
16.240900000000003
16.240900000000003
Cache types and constructors
ImmersedLayers.BasicILMCache
— Typestruct BasicILMCache{N, SCA<:AbstractScalingType, GCT, ND, BLT<:BodyList, NT<:VectorData, DST<:ScalarData, REGT<:Regularize, RSNT<:RegularizationMatrix, ESNT<:InterpolationMatrix, RT<:RegularizationMatrix, ET<:InterpolationMatrix, RCT, ECT, RDT, EDT, LT<:Laplacian, GVT, GNT, GDT, SVT, SDT, SST} <: ImmersedLayers.AbstractBasicCache{N, GCT}
A cache of operators and storage data for use in surface operations. Constructed with SurfaceScalarCache
or SurfaceVectorCache
.
ImmersedLayers.SurfaceScalarCache
— FunctionSurfaceScalarCache(g::PhysicalGrid[,scaling=IndexScaling])
Create a cache of type BasicILMCache
with scalar grid data, using the grid specified in g
, with no immersed points. The keyword scaling
can be used to set the scaling in the operations. By default, it is set to IndexScaling
which sets the differential operators to be only differencing operators. By using scaling = GridScaling
, then the grid and spacings are accounted for and differential operators are scaled by this spacing. The keyword phys_params
can be used to supply physical parameters.
SurfaceScalarCache(body::Body/BodyList,g::PhysicalGrid[,ddftype=CartesianGrids.Yang3][,scaling=GridScaling])
Create a cache of type BasicILMCache
, holding operators and storage data for use in immersed layer operations on scalar data. This is sometimes called from withinILMSystem
rather than directly.
The body
can be of type Body
or BodyList
. The keyword scaling
can be used to set the scaling in the operations. By default, it is set to IndexScaling
which sets the regularization and interpolation to be symmetric matrices (i.e., interpolation is the adjoint of regularization with respect to a vector dot product), and the vector calculus operations on the grid are simple differences. By using scaling = GridScaling
, then the grid and point spacings are accounted for. Interpolation and regularization are adjoints with respect to inner products based on discretized surface and volume integrals, and vector calculus operations are scaled by the grid spacing.
SurfaceScalarCache(X::VectorData,g::PhysicalGrid[,ddftype=CartesianGrids.Yang3][,scaling=GridScaling])
Create a cache of type BasicILMCache
, holding operators and storage data for use in immersed layer operations on scalar data. The X
specifies the is assumed to hold the endpoints of the immersed surface segments, and g
the physical grid.
ImmersedLayers.SurfaceVectorCache
— FunctionSurfaceVectorCache(g::PhysicalGrid[,scaling=IndexScaling])
Create a cache of type BasicILMCache
with vector grid data, with no immersed points. See SurfaceScalarCache
for details.
SurfaceVectorCache(body::Body/BodyList,g::PhysicalGrid[,ddftype=CartesianGrids.Yang3][,scaling=Gridcaling])
Create a cache of type BasicILMCache
, holding operators and storage data for use in immersed layer operations on vector data. See SurfaceScalarCache
for details.
SurfaceVectorCache(X::VectorData,g::PhysicalGrid[,ddftype=CartesianGrids.Yang3][,scaling=GridScaling])
Create a cache of type BasicILMCache
, holding operators and storage data for use in immersed layer operations on vector data. See SurfaceScalarCache
for details.
ImmersedLayers.AbstractExtraILMCache
— Typeabstract type AbstractExtraILMCache
When defining problem-specific cache, make it a subtype of this.
Utilities for creating instances of data
ImmersedLayers.similar_grid
— Functionsimilar_grid(::BasicILMCache)
Get a similar
copy of the basic grid data in the cache.
ImmersedLayers.similar_gridgrad
— Functionsimilar_gridgrad(::BasicILMCache)
Get a similar
copy of the gradient of the grid data in the cache.
ImmersedLayers.similar_gridcurl
— Functionsimilar_gridcurl(::BasicILMCache)
Get a similar
copy of the grid curl field data in the cache.
ImmersedLayers.similar_griddiv
— Functionsimilar_griddiv(::BasicILMCache)
Get a similar
copy of the grid div field data in the cache.
ImmersedLayers.similar_gridgradcurl
— Functionsimilar_gridgradcurl(::BasicILMCache)
Get a similar
copy of the grid gradient-of-curl field data in the cache.
ImmersedLayers.similar_surface
— Functionsimilar_surface(::BasicILMCache)
Get a similar
copy of the basic surface point data in the cache.
ImmersedLayers.similar_surfacescalar
— Functionsimilar_surfacescalar(::BasicILMCache)
Get a similar
copy of the surface scalar point data in the cache. This is only used for vector-type caches. Otherwise, similar_surface
should be used.
ImmersedLayers.zeros_grid
— Functionzeros_grid(::BasicILMCache)
Get an instance of the basic grid data in the cache, with values set to zero.
ImmersedLayers.zeros_gridgrad
— Functionzeros_gridgrad(::BasicILMCache)
Get an instance of the gradient of the grid data in the cache, with values set to zero.
ImmersedLayers.zeros_gridcurl
— Functionzeros_gridcurl(::BasicILMCache)
Get an instance of the grid curl field data in the cache, with values set to zero.
ImmersedLayers.zeros_griddiv
— Functionzeros_griddiv(::BasicILMCache)
Get an instance of the grid div field data in the cache, with values set to zero.
ImmersedLayers.zeros_gridgradcurl
— Functionzeros_gridgradcurl(::BasicILMCache)
Get an instance of the grid gradient-of-curl field data in the cache, with values set to zero.
ImmersedLayers.zeros_surface
— Functionzeros_surface(::BasicILMCache)
Get an instance of the basic surface point data in the cache, with values set to zero.
ImmersedLayers.zeros_surfacescalar
— Functionzeros_surfacescalar(::BasicILMCache)
Get an instance of the surface scalar point data in the cache, with values set to zero. This is only used for vector-type caches. Otherwise, zeros_surface
should be used.
ImmersedLayers.ones_grid
— Functionones_grid(::BasicILMCache)
Get an instance of the basic grid data in the cache, with values set to unity.
ones_grid(::BasicILMCache,dim)
For a vector-type cache, get an instance of the basic grid data in the cache, with values set to unity in dimension dim
.
ImmersedLayers.ones_gridgrad
— Functionones_gridgrad(::BasicILMCache)
Get an instance of the gradient of the grid data in the cache, with values set to unity.
ones_gridgrad(::BasicILMCache,dim)
Get an instance of the gradient of the grid data in the cache, in direction dim
, with values set to unity. If the data are of type TensorGridData
, then dim
takes values from 1 to 2^2.
ImmersedLayers.ones_gridcurl
— Functionones_gridcurl(::BasicILMCache)
Get an instance of the grid curl field data in the cache, with values set to unity.
ImmersedLayers.ones_griddiv
— Functionones_griddiv(::BasicILMCache)
Get an instance of the grid div field data in the cache, with values set to unity.
ImmersedLayers.ones_gridgradcurl
— Functionones_gridgradcurl(::BasicILMCache)
Get an instance of the grid gradient-of-curl field data in the cache, with values set to unity.
ImmersedLayers.ones_surface
— Functionones_surface(::BasicILMCache)
Get an instance of the basic surface point data in the cache, with values set to unity.
ones_surface(::BasicILMCache,dim)
Get an instance of the basic surface point data in the cache, with values set to unity in dimension dim
. This only works for a vector-type cache.
ImmersedLayers.ones_surfacescalar
— Functionones_surfacescalar(::BasicILMCache)
Get an instance of the surface scalar point data in the cache, with values set to unity. This is only used for vector-type caches. Otherwise, ones_surface
should be used.
ImmersedLayers.x_grid
— Functionx_grid(::BasicILMCache)
Return basic grid data filled with the grid x
coordinate. If the grid data is scalar, then the output of this function is of the same type. If the grid data is vector, then the output is of type Edges{Primal}
, and the coordinates of each component are in u
, v
fields.
ImmersedLayers.y_grid
— Functiony_grid(::BasicILMCache)
Return basic grid data filled with the grid y
coordinate. If the grid data is scalar, then the output of this function is of the same type. If the grid data is vector, then the output is of type Edges{Primal}
, and the coordinates of each component are in u
, v
fields.
ImmersedLayers.x_gridcurl
— Functionx_gridcurl(::BasicILMCache)
Return basic grid curl field data filled with the grid x
coordinate
ImmersedLayers.y_gridcurl
— Functiony_gridcurl(::BasicILMCache)
Return basic grid curl field data filled with the grid y
coordinate
ImmersedLayers.x_griddiv
— Functionx_griddiv(::BasicILMCache)
Return basic grid div field data filled with the grid x
coordinate
ImmersedLayers.y_griddiv
— Functiony_griddiv(::BasicILMCache)
Return basic grid div field data filled with the grid y
coordinate
ImmersedLayers.x_gridgrad
— Functionx_gridgrad(::BasicILMCache)
Return basic grid gradient field data filled with the grid x
coordinate. If the grid data is scalar, then the output of this function is of Edges{Primal}
type and the coordinate field for each component can be accessed with the u
and v
field, respectively. If the grid data is vector, then the output is of type EdgeGradient{Primal}
, and the coordinates are in dudx
, dudy
, dvdx
, and dvdy
fields.
ImmersedLayers.y_gridgrad
— Functiony_gridgrad(::BasicILMCache)
Return basic grid gradient field data filled with the grid y
coordinate. If the grid data is scalar, then the output of this function is of Edges{Primal}
type and the coordinate field for each component can be accessed with the u
and v
field, respectively. If the grid data is vector, then the output is of type EdgeGradient{Primal}
, and the coordinates are in dudx
, dudy
, dvdx
, and dvdy
fields.
ImmersedLayers.evaluate_field!
— Functionevaluate_field!(f::GridData,s::AbstractSpatialField,cache)
Evaluate a spatial field s
in place as grid data f
.
evaluate_field!(f::GridData,t,s::AbstractSpatialField,cache)
Evaluate a spatial-temporal field s
at time t
in place as grid data f
.
Utilities for accessing surface information
ImmersedLayers.areas
— Methodareas(cache::BasicILMCache)
Return the areas (as ScalarData
) of the surface panels associated with cache
ImmersedLayers.normals
— Methodnormals(cache::BasicILMCache)
Return the normals (as VectorData
) of the surface points associated with cache
ImmersedLayers.points
— Methodpoints(cache::BasicILMCache)
Return the coordinates (as VectorData
) of the surface points associated with cache
ImmersedLayers.arcs
— Methodarcs(cache::BasicILMCache)
Return ScalarData
of arc length coordinates for the body surface(s) in the given cache cache
.
Inner products, norms, and integrals
LinearAlgebra.dot
— Methoddot(u1::GridData,u2::GridData,cache::BasicILMCache)
Calculate the inner product of grid data u1
and u2
, using the scaling associated with cache
.
LinearAlgebra.dot
— Methoddot(u1::PointData,u2::PointData,cache::BasicILMCache)
Calculate the inner product of surface point data u1
and u2
, using the scaling associated with cache
.
LinearAlgebra.norm
— Methodnorm(u::GridData,cache::BasicILMCache)
Calculate the norm of grid data u
, using the scaling associated with cache
.
CartesianGrids.integrate
— Methodintegrate(u::GridData,cache::BasicILMCache)
Calculate the discrete volume integral of grid data u
. If u
is VectorGridData
, then this returns a vector of the integrals in each coordinate direction. This integral produces the same result regardless if GridScaling
or IndexScaling
is used.
LinearAlgebra.norm
— Methodnorm(u::PointData,cache::BasicILMCache)
Calculate the norm of surface point data u
, using the scaling associated with cache
.
CartesianGrids.integrate
— Methodintegrate(u::PointData,cache::BasicILMCache)
Calculate the discrete surface integral of data u
on the immersed points in cache
. This uses trapezoidal rule quadrature. If u
is VectorData
, then this returns a vector of the integrals in each coordinate direction. This operation produces the same effect, regardless if cache
is set up for GridScaling
or IndexScaling
. In both cases, the surface element areas are used.
LinearAlgebra.dot
— Methoddot(u1::PointData,u2::PointData,cache::BasicILMCache,i)
Calculate the inner product of surface point data u1
and u2
for body i
in the cache cache
, scaling as appropriate for this cache.
LinearAlgebra.norm
— Methodnorm(u::PointData,cache::BasicILMCache,i::Int)
Calculate the norm of surface point data u
, using the scaling associated with cache
, for body i
in the body list of cache
.
CartesianGrids.integrate
— Methodintegrate(u::PointData,cache::BasicILMCache,i::Int)
Calculate the discrete surface integral of scalar data u
on the immersed points in cache
, on body i
in the body list in cache
. This uses trapezoidal rule quadrature. If u
is VectorData
, then this returns a vector of the integrals in each coordinate direction. This operation produces the same effect, regardless if cache
is set up for GridScaling
or IndexScaling
. In both cases, the surface element areas are used.
Other cache utilities
Base.view
— Methodview(u::PointData,cache::BasicILMCache,i::Int)
Provide a view
of point data u
corresponding to body i
in the list of bodies in cache
.
Base.copyto!
— Methodcopyto!(u::PointData,v::PointData,cache::BasicILMCache,i::Int)
Copy the data in the elements of v
associated with body i
in the body list in cache
to the corresponding elements in u
. These data must be of the same type (e.g., ScalarData
or VectorData
) and have the same length.
Base.copyto!
— Methodcopyto!(u::ScalarData,v::AbstractVector,cache::BasicILMCache,i::Int)
Copy the data in v
to the elements in u
associated with body i
in the body list in cache
. v
must have the same length as this subarray of u
associated with i
.
CartesianGrids.Laplacian
— MethodLaplacian(cache::BasicILMCache,coeff_factor::Real[,with_inverse=true])
Create an invertible Laplacian operator for the grid in cache
, using the index or grid scaling associated with cache
. The operator is pre-multiplied by the factor coeff_factor
.
CartesianGrids.RegularizationMatrix
— MethodRegularizationMatrix(cache::BasicILMCache,src::PointData,trg::GridData)
Create a regularization matrix for regularizing point data of type src
to grid data of type trg
. (Both src
and trg
must be appropriately sized for the grid and points in cache
.)
CartesianGrids.InterpolationMatrix
— MethodInterpolationMatrix(cache::BasicILMCache,src::GridData,trg::PointData)
Create a interpolation matrix for regularizing grid data of type src
to point data of type trg
. (Both src
and trg
must be appropriately sized for the grid and points in cache
.)
This page was generated using Literate.jl.
- 1Yang, X., et al., (2009) "A smoothing technique for discrete delta functions with application to immersed boundary method in moving boundary simulations," J. Comput. Phys., 228, 7821–7836.