Interpolation
Interpolation types
Tesserae.BSpline
— TypeBSpline(degree)
B-spline kernel. degree
is one of Linear()
, Quadratic()
or Cubic()
.
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.
Tesserae.SteffenBSpline
— TypeSteffenBSpline(degree)
B-spline kernel with boundary correction by Steffen et al.[Steffen] SteffenBSpline
satisfies the partition of unity, $\sum_i w_{ip} = 1$, near boundaries. See also KernelCorrection
.
Tesserae.uGIMP
— TypeuGIMP()
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
Tesserae.CPDI
— TypeCPDI()
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
Tesserae.WLS
— TypeWLS(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
.
Tesserae.KernelCorrection
— TypeKernelCorrection(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
.
Basis function values
Tesserae.MPValue
— TypeMPValue([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
- SteffenSteffen, M., Kirby, R. M., & Berzins, M. (2008). Analysis and reduction of quadrature errors in the material point method (MPM). International journal for numerical methods in engineering, 76(6), 922-948.
- GIMPBardenhagen, S. G., & Kober, E. M. (2004). The generalized interpolation material point method. Computer Modeling in Engineering and Sciences, 5(6), 477-496.
- CPDISadeghirad, A., Brannon, R.M. and Burghardt, J., 2011. A convected particle domain interpolation technique to extend applicability of the material point method for problems involving massive deformations. International Journal for numerical methods in Engineering, 86(12), pp.1435-1456.
- MLSMPMHu, Y., Fang, Y., Ge, Z., Qu, Z., Zhu, Y., Pradhana, A. and Jiang, C., 2018. A moving least squares material point method with displacement discontinuity and two-way rigid body coupling. ACM Transactions on Graphics (TOG), 37(4), pp.1-14.
- KCNakamura, K., Matsumura, S., & Mizutani, T. (2023). Taylor particle-in-cell transfer and kernel correction for material point method. Computer Methods in Applied Mechanics and Engineering, 403, 115720.