Interpolation

Interpolation types

Tesserae.BSplineType
BSpline(degree)

B-spline kernel. degree is one of Linear(), Quadratic() or Cubic().

Warning

BSpline(Quadratic()) and BSpline(Cubic()) cannot handle boundaries correctly because the kernel values are merely truncated, which leads to unstable behavior. Therefore, it is recommended to use either SteffenBSpline or KernelCorrection in cases where proper handling of boundaries is necessary.

source
Tesserae.uGIMPType
uGIMP()

A kernel for the unchanged generalized interpolation material point (uGIMP) [GIMP]. uGIMP requires the initial particle length l in the particle property as follows:

ParticleProp = @NamedTuple begin
    < variables... >
    l :: Float64
end
source
Tesserae.CPDIType
CPDI()

A kernel for convected particle domain interpolation (CPDI) [CPDI]. CPDI requires the initial particle length l and the deformation gradient F in the particle property. For example, in two dimensions, the property is likely to be as follows:

ParticleProp = @NamedTuple begin
    < variables... >
    l :: Float64
    F :: Mat{2, 2, Float64, 4}
end
source
Tesserae.WLSType
WLS(kernel)

WLS performs a local weighted least squares fit for the kernel. This results in the same kernel used in moving least squares MPM[MLSMPM]. kernel is one of BSpline and uGIMP.

source
Tesserae.KernelCorrectionType
KernelCorrection(kernel)

KernelCorrection[KC] modifies kernel to achieve stable simulations near boundaries. The corrected kernel satisfies not only the partition of unity, $\sum_i w_{ip} = 1$, but also the linear field reproduction, $\sum_i w_{ip} \bm{x}_i = \bm{x}_p$, near boundaries. In the implementation, this simply applies WLS near boundaries. kernel is one of BSpline and uGIMP. See also SteffenBSpline.

source

Basis function values

Tesserae.MPValueType
MPValue([T,] interpolation, mesh)

MPValue stores properties for interpolation, such as the basis function values and its derivatives.

julia> mesh = CartesianMesh(1.0, (0,5), (0,5));

julia> xₚ = Vec(2.2, 3.4); # particle position

julia> mp = MPValue(BSpline(Quadratic()), mesh);

julia> update!(mp, xₚ, mesh) # update `mp` at position `xₚ` in `mesh`
MPValue:
  Interpolation: BSpline(Quadratic())
  Property names: w::Matrix{Float64}, ∇w::Matrix{Vec{2, Float64}}
  Neighboring nodes: CartesianIndices((2:4, 3:5))

julia> sum(mp.w) ≈ 1 # partition of unity
true

julia> nodeindices = neighboringnodes(mp) # grid indices within a particles' local domain
CartesianIndices((2:4, 3:5))

julia> sum(eachindex(nodeindices)) do ip # linear field reproduction
           i = nodeindices[ip]
           mp.w[ip] * mesh[i]
       end ≈ xₚ
true
source