Matrix operators

Many solutions of partial differential equations with immersed layers, particularly Poisson problems, lead to saddle-point problems, in which the Schur complement operator is a matrix composed from some of the surface-grid operators discussed in Surface-grid operations. The package provides some convenience tools for constructing these matrices. By their nature, the construction of these matrices is slow, since each column involves the application of the same set of operations. However, the point of this construction is to do it once and store it for repeated application.

One common saddle point system is

\[A = \begin{bmatrix} L & R \\ R^T & 0 \end{bmatrix}\]

where $L$ is the discrete Laplacian and $R$ and $R^T$ are the regularization and interpolation operators (regularize! and interpolate!), respectively. This system arises in the solution of the Poisson equation with Dirichlet boundary conditions on the immersed surface. The Schur complement of this is $S = - R^T L^{-1} R$. This matrix can be obtained using the function create_RTLinvR.

Another common saddle point system is

\[A = \begin{bmatrix} L & D_s \\ G_s & R_n^T R_n \end{bmatrix}\]

where $D_s$ and $G_s$ are the surface divergence and gradient operators (surface_divergence! and surface_grad!), respectively, and $R_n$ and $R_n^T$ are regularize_normal! and normal_interpolate!. This system arises in the solution of the Poisson equation with Neumann boundary conditions on the immersed surface. The Schur complement of this is $S = R_n^T R_n - G_s L^{-1} D_s$. Each of the matrices in this are individually provided by the package, by the functions create_nRTRn and create_GLinvD, respectively. However, it is useful to know that the sum of these two matrices is exactly the matrix $-C_s L^{-1}C_s^T$, where $C_s$ and $C_s^T$ are surface curl operators surface_curl!. This complete matrix is provided by create_CLinvCT.

Another helpful matrix operator is the surface filter, given by

\[\tilde{R}^T R\]

where $\tilde{R}^T$ is a modified form of the interpolation operator, designed to return the regularized field to the surface points while maintaining the integral value of the original field [1]. We can obtain this matrix with create_surface_filter.

Matrix construction functions

ImmersedLayers.create_RTLinvRFunction
create_RTLinvR(cache::BasicILMCache[;scale=1.0])

Using the provided cache cache, construct the square matrix $-R^T L^{-1}R$, which maps data of type ScalarData to data of the same type. The operators R^T and R correspond to interpolate! and regularize! and L^{-1} is the inverse of the grid Laplacian. The optional keyword scale multiplies the matrix by the designated value.

source
ImmersedLayers.create_CLinvCTFunction
create_CLinvCT(cache::BasicILMCache[;scale=1.0])

Using the provided cache cache, construct the square matrix $-C_s L^{-1}C_s^T$, which maps data of the primary point data type of the cache to data of the same type. The operators C_s and C_s^T correspond to surface_curl! and L is the grid Laplacian. The optional keyword scale multiplies the matrix by the designated value.

source
ImmersedLayers.create_CLinvCT_scalarFunction
create_CLinvCT_scalar(cache::BasicILMCache[;scale=1.0])

Using the provided cache cache, construct the square matrix $-C_s L^{-1}C_s^T$, which maps data of the scalar point data type of the cache to data of the same type. The operators C_s and C_s^T correspond to surface_curl! and L is the grid Laplacian. The optional keyword scale multiplies the matrix by the designated value.

source
ImmersedLayers.create_CL2invCTFunction
create_CL2invCT(cache::BasicILMCache[;scale=1.0])

Using the provided cache cache, construct the square matrix $-C_s L^{-2} C_s^T$, which maps data of the primary point data type of the cache to data of the same type. The operators C_s and C_s^T correspond to surface_curl! and L is the grid Laplacian. The optional keyword scale multiplies the matrix by the designated value.

source
ImmersedLayers.create_GLinvDFunction
create_GLinvD(cache::BasicILMCache[;scale=1.0])

Using the provided cache cache, construct the square matrix $-G_s L^{-1}D_s$, which maps data of type ScalarData to data of the same type. The operators G_s and D_s correspond to surface_grad! and surface_divergence!, and L is the grid Laplacian. The optional keyword scale multiplies the matrix by the designated value.

source
ImmersedLayers.create_GLinvD_crossFunction
create_GLinvD_cross(cache::BasicILMCache[;scale=1.0])

Using the provided cache cache, construct the square matrix $-\hat{G}_s L^{-1}\hat{D}_s$, which maps data of type ScalarData to data of the same type. The operators G_s and D_s correspond to surface_grad_cross! and surface_divergence_cross!, and L is the grid Laplacian. The optional keyword scale multiplies the matrix by the designated value.

source
ImmersedLayers.create_GLinvD_symmFunction
create_GLinvD_symm(cache::BasicILMCache[;scale=1.0])

Using the provided cache cache, construct the square matrix $-G_s L^{-1}D_s$, which maps data of type ScalarData to data of the same type. The operators G_s and D_s correspond to surface_grad_symm! and surface_divergence_symm!, and L is the grid Laplacian. The optional keyword scale multiplies the matrix by the designated value.

source
ImmersedLayers.create_nRTRnFunction
create_nRTRn(cache::BasicILMCache[;scale=1.0])

Using the provided cache cache, construct the square matrix $n\cdot R_f^T R_f n \circ$, which maps data of type ScalarData to data of the same type. The operators R_f^T and R_f correspond to the interpolation and regularization matrices. The optional keyword scale multiplies the matrix by the designated value.

source
ImmersedLayers.create_surface_filterFunction
create_surface_filter(cache::BasicILMCache)

Create a surface filtering matrix operator $\tilde{R}^T R$, where $\tilde{R}^T$ represents a modified version of the interpolation operator. The resulting matrix can be applied to surface data to filter out high-frequency components.

source

This page was generated using Literate.jl.

  • 1Goza, A., et al., (2016) "Accurate computation of surface stresses and forces with immersed boundary methods," J. Comput. Phys., 321, 860–873.