Grid operations
Basic differential operations
There are a variety of (purely) grid-based operators that are useful for carrying out calculations in immersed layer problems. We will demonstrate a few of them here. We will start by generating the cache, just as we did in Immersed layer caches
using ImmersedLayers
using Plots
Set up a grid and cache
Δ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)
We still generate a cache for these operations, but now, we only supply the grid. There are no immersed surfaces for this demonstration.
cache = SurfaceScalarCache(g,scaling=GridScaling)
Surface cache with scaling of type GridScaling
0 point data of type ScalarData{0, Float64, Vector{Float64}}
Grid data of type Nodes{Primal, 405, 405, Float64, Matrix{Float64}}
To demonstrate, let's generate a Gaussian
p = zeros_grid(cache)
xg, yg = x_grid(cache), y_grid(cache)
p .= exp.(-(xg∘xg)-(yg∘yg))
Nodes{Primal, 405, 405, Float64, Matrix{Float64}} data
Printing in grid orientation (lower left is (1,1))
404×404 Matrix{Float64}:
0.00029738 0.000309547 0.000322148 … 0.00029738 0.000285634
0.000309609 0.000322277 0.000335396 0.000309609 0.00029738
0.000322277 0.000335463 0.000349118 0.000322277 0.000309547
0.000335396 0.000349118 0.00036333 0.000335396 0.000322148
0.000348979 0.000363257 0.000378044 0.000348979 0.000335194
0.000363039 0.000377893 0.000393276 … 0.000363039 0.0003487
0.000377591 0.00039304 0.000409039 0.000377591 0.000362676
0.000392647 0.000408712 0.000425349 0.000392647 0.000377138
0.000408222 0.000424924 0.000442221 0.000408222 0.000392098
0.00042433 0.000441691 0.000459671 0.00042433 0.000407569
⋮ ⋱
0.00042433 0.000441691 0.000459671 … 0.00042433 0.000407569
0.000408222 0.000424924 0.000442221 0.000408222 0.000392098
0.000392647 0.000408712 0.000425349 0.000392647 0.000377138
0.000377591 0.00039304 0.000409039 0.000377591 0.000362676
0.000363039 0.000377893 0.000393276 0.000363039 0.0003487
0.000348979 0.000363257 0.000378044 … 0.000348979 0.000335194
0.000335396 0.000349118 0.00036333 0.000335396 0.000322148
0.000322277 0.000335463 0.000349118 0.000322277 0.000309547
0.000309609 0.000322277 0.000335396 0.000309609 0.00029738
Now, let's generate the gradient of these data
v = zeros_gridgrad(cache)
grad!(v,p,cache)
plot(v,cache)
We can then compute the derivative of this data
divv = zeros_grid(cache)
divergence!(divv,v,cache)
plot(divv,cache)
Convective derivatives
And finally, let's compute convective derivatives. First, we will compute
\[\mathbf{v}\cdot\nabla p\]
For this operation, we create a special additional cache using ConvectiveDerivativeCache
. This extra cache holds additional memory for making the calculation of the convective derivative faster if we compute it often.
cdcache = ConvectiveDerivativeCache(cache)
vdp = zeros_grid(cache)
@time convective_derivative!(vdp,v,p,cache,cdcache)
0.001086 seconds (18 allocations: 3.078 KiB)
Plot it
plot(vdp,cache)
Now, let's compute
\[\mathbf{v}\cdot\nabla\mathbf{v}\]
For this, we create a cache for VectorGridData
, and a new instance of ConvectiveDerivativeCache
to go along with it.
vcache = SurfaceVectorCache(g,scaling=GridScaling)
vdv = zeros_grid(vcache)
cdvcache = ConvectiveDerivativeCache(vcache)
@time convective_derivative!(vdv,v,vcache,cdvcache)
0.002948 seconds (55 allocations: 8.875 KiB)
Plot it
plot(vdv,vcache)
Finally, let's compute
\[(\curl\mathbf{v})\times\mathbf{v}\]
For this, we create a cache called RotConvectiveDerivativeCache
to go along with it.
w = zeros_gridcurl(vcache)
curl!(w,v,vcache)
wv = zeros_grid(vcache)
cdrcache = RotConvectiveDerivativeCache(vcache)
@time w_cross_v!(wv,w,v,vcache,cdrcache)
Edges{Primal, 405, 405, Float64, Vector{Float64}} data
u (in grid orientation)
404×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
v (in grid orientation)
405×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
Plot it
plot(wv,vcache)
Grid operator functions
CartesianGrids.divergence!
— Functiondivergence!(p::Nodes{Primal/Dual},v::Edges{Primal/Dual},cache::BasicILMCache)
divergence!(p::Nodes{Primal/Dual},v::Edges{Primal/Dual},sys::ILMSystem)
Compute the discrete divergence of v
and return it in p
, scaling it by the grid spacing if cache
(or sys
) is of GridScaling
type, or leaving it as a simple differencing if cache
(or sys
) is of IndexScaling
type.
CartesianGrids.grad!
— Functiongrad!(v::Edges{Primal/Dual},p::Nodes{Primal/Dual},cache::BasicILMCache)
grad!(v::Edges{Primal/Dual},p::Nodes{Primal/Dual},sys::ILMSystem)
Compute the discrete gradient of p
and return it in v
, scaling it by the grid spacing if cache
(or sys
) is of GridScaling
type, or leaving it as a simple differencing if cache
(or sys
) is of IndexScaling
type.
grad!(v::EdgeGradient{Primal/Dual,Dual/Primal},p::Edges{Primal/Dual},cache::BasicILMCache)
grad!(v::EdgeGradient{Primal/Dual,Dual/Primal},p::Edges{Primal/Dual},sys::ILMSystem)
Compute the discrete gradient of edge data p
and return it in v
, scaling it by the grid spacing if cache
(or sys
) is of GridScaling
type, or leaving it as a simple differencing if cache
(or sys
) is of IndexScaling
type.
CartesianGrids.curl!
— Functioncurl!(v::Edges{Primal},s::Nodes{Dual},cache::BasicILMCache)
curl!(v::Edges{Primal},s::Nodes{Dual},sys::ILMSystem)
Compute the discrete curl of s
and return it in v
, scaling it by the grid spacing if cache
(or sys
) is of GridScaling
type, or leaving it as a simple differencing if cache
(or sys
) is of IndexScaling
type.
curl!(w::Nodes{Dual},v::Edges{Primal},cache::BasicILMCache)
curl!(w::Nodes{Dual},v::Edges{Primal},sys::ILMSystem)
Compute the discrete curl of v
and return it in w
, scaling it by the grid spacing if cache
(or sys
) is of GridScaling
type, or leaving it as a simple differencing if cache
(or sys
) is of IndexScaling
type.
CartesianGrids.convective_derivative!
— Functionconvective_derivative!(vdp::Nodes{Primal},v::Edges,p::Nodes{Primal},base_cache::BasicILMCache,extra_cache::ConvectiveDerivativeCache)
Compute the convective derivative of p
with velocity v
, i.e., $v\cdot \nabla p$, and return the result in vdp
. The result is either divided by unity or the grid cell size depending on whether base_cache
is of type IndexScaling
or GridScaling
. This version of the method uses extra_cache
of type ConvectiveDerivativeCache
.
convective_derivative!(vdw::Nodes{Dual},v::Edges,w::Nodes{Dual},base_cache::BasicILMCache,extra_cache::ConvectiveDerivativeCache)
Compute the convective derivative of w
with velocity v
, i.e., $v\cdot \nabla w$, and return the result in vdw
. The result is either divided by unity or the grid cell size depending on whether base_cache
is of type IndexScaling
or GridScaling
. This version of the method uses extra_cache
of type ConvectiveDerivativeCache
.
convective_derivative!(vdu::Edges,v::Edges,u::Edges,base_cache::BasicILMCache,extra_cache::ConvectiveDerivativeCache)
Compute the convective derivative of u
, i.e., $v\cdot \nabla u$, and return the result in vdu
. The result is either divided by unity or the grid cell size depending on whether base_cache
is of type IndexScaling
or GridScaling
. This version of the method uses extra_cache
of type ConvectiveDerivativeCache
.
convective_derivative!(vdv::Edges,v::Edges,base_cache::BasicILMCache,extra_cache::ConvectiveDerivativeCache)
Compute the convective derivative of v
, i.e., $v\cdot \nabla v$, and return the result in vdv
. The result is either divided by unity or the grid cell size depending on whether base_cache
is of type IndexScaling
or GridScaling
. This version of the method uses extra_cache
of type ConvectiveDerivativeCache
.
ImmersedLayers.convective_derivative
— Functionconvective_derivative(v::Edges{Primal},p::Nodes{Primal},base_cache::BasicILMCache)
Compute the convective derivative of p
with velocity v
, i.e., $v\cdot \nabla p$. The result is either divided by unity or the grid cell size depending on whether base_cache
is of type IndexScaling
or GridScaling
.
convective_derivative(v::Edges{Dual},w::Nodes{Dual},base_cache::BasicILMCache)
Compute the convective derivative of w
with velocity v
, i.e., $v\cdot \nabla w$. The result is either divided by unity or the grid cell size depending on whether base_cache
is of type IndexScaling
or GridScaling
.
convective_derivative(v::Edges,base_cache::BasicILMCache)
Compute the convective derivative of v
, i.e., $v\cdot \nabla v$. The result is either divided by unity or the grid cell size depending on whether base_cache
is of type IndexScaling
or GridScaling
.
ImmersedLayers.ConvectiveDerivativeCache
— TypeConvectiveDerivativeCache(dv::EdgeGradient)
Create a cache (a subtype of AbstractExtraILMCache
) for computing the convective derivative, using dv
to define the cache data.
ImmersedLayers.w_cross_v!
— Functionw_cross_v!(vw::Edges{Primal},w::Nodes{Dual},v::Edges{Primal},base_cache::BasicILMCache,extra_cache::RotConvectiveDerivativeCache)
Compute the term w \times v
, with vorticity w
and velocity v
, and return the result in vw
. This version of the method uses extra_cache
of type RotConvectiveDerivativeCache
.
ImmersedLayers.w_cross_v
— Functionw_cross_v(w::Nodes{Dual},v::Edges{Primal},base_cache::BasicILMCache)
Compute the term w \times v
, with vorticity w
and velocity v
.
ImmersedLayers.RotConvectiveDerivativeCache
— TypeRotConvectiveDerivativeCache(w::Nodes{Dual})
Create a cache (a subtype of AbstractExtraILMCache
) for computing the convective derivative in rotational form, using w
to define the cache data.
ImmersedLayers.laplacian!
— Functionlaplacian!(w::GridData,s::GridData,sys::ILMSystem)
Compute the Laplacian of grid data s
, and divide the result by unity or by the grid cell size, depending on whether sys
has IndexScaling
or GridScaling
, respectively, and return the result as w
.
laplacian!(w::GridData,s::GridData,cache::BasicILMCache)
Compute the Laplacian of grid data s
, and divide the result by unity or by the grid cell size, depending on whether cache
has IndexScaling
or GridScaling
, respectively, and return the result as w
.
ImmersedLayers.inverse_laplacian!
— Functioninverse_laplacian!(w::GridData,sys::ILMSystem)
Compute the in-place inverse Laplacian of grid data w
, and multiply the result by unity or by the grid cell size, depending on whether sys
has IndexScaling
or GridScaling
, respectively.
inverse_laplacian!(w::GridData,cache::BasicILMCache)
Compute the in-place inverse Laplacian of grid data w
, and multiply the result by unity or by the grid cell size, depending on whether cache
has IndexScaling
or GridScaling
, respectively.
inverse_laplacian!(sol::GridData,rhs::GridData,cache::BasicILMCache)
Compute the iinverse Laplacian of grid data rhs
, multiply the result by unity or by the grid cell size, depending on whether cache
has IndexScaling
or GridScaling
, respectively, and return the result as sol
.
This page was generated using Literate.jl.