jInv.Mesh

Discretization of differential operators is an essential ingredient of PDE parameter estimation problems. Depending on the problem and computational resources, different mesh types are needed in practical applications. The Mesh submodule of jInv provides different mesh geometry under the abstract type AbstractMesh. Currently, the Mesh submodule methods and provides operators on regular and tensor meshes but is easily extensible.

Regular Meshes

List of types and methods

# jInv.Mesh.RegularMeshType.

type jInv.Mesh.RegularMesh <: AbstractTensorMesh

Regular mesh in 1D, 2D, and 3D

Fields:

    domain::Vector{Float64}  - physical domain [min(x1) max(x1) min(x2) max(x2)]
    h::Vector{Float64}       - cell size
    x0::Vector{Float64}      - origin
    dim::Int                 - dimension of mesh
    n::Vector{Int64}         - number of cells in each dimension
    nc::Int                  - total number of cells
    nf::Vector{Int64}        - number of faces in each dimension
    ne::Vector{Int64}        - number of edges in each dimension


    Persistent Operators:

    Operators should not be accessed directly. They will be built, if needed,
    when accessing them using specified method. clear!(M) will release all 
    memory.

        Div::SparseMatrixCSC    - divergence (faces -> cell-centers)
                                  Access via: getDivergenceMatrix(M)
        Grad::SparseMatrixCSC   - gradient (nodal -> edges)
                                  Access via: getNodalGradientMatrix(M)
        Curl::SparseMatrixCSC   - curl (edges -> faces)
                                  Access via: getCurlMatrix(M)
        Af::SparseMatrixCSC     - face average (faces -> cell-centers)
                                  Access via: getFaceAverageMatrix(M)
        Ae::SparseMatrixCSC     - edge average (edges -> cell-centers)
                                  Access via: getEdgeAverageMatrix(M)
        An::SparseMatrixCSC     - nodal average (nodes -> cell-centers)
                                  Access via: getNodalAverageMatrix(M)
        V::SparseMatrixCSC      - cell volumes (diagonal matrix)
                                  Access via: getVolume(M)
        F::SparseMatrixCSC      - face area (diagonal matrix)
                                  Access via: getFaceArea(M)
        L::SparseMatrixCSC      - edge length (diagonal matrix)
                                  Access via: getLength(M)
        Vi::SparseMatrixCSC     - inverse cell volumes (diagonal matrix)
                                  Access via: getVolumeInv(M)
        Fi::SparseMatrixCSC     - inverse face area (diagonal matrix)
                                  Access via: getFaceAreaInv(M)
        Li::SparseMatrixCSC     - inverse edge length (diagonal matrix)
                                  Access via: getLengthAreaInv(M)
        nLap::SparseMatrixCSC   - nodal Laplacian
                                  Access via: getNodalLaplacian(M)

Examples:
M2D  = getRegularMesh([1.2 2.4 2.2 5.0],[3,4])  
M3D  = getRegularMesh([1.2 2.4 2.2 5.0 0 1],[3,4,7])

source

# jInv.Mesh.TensorMesh3DType.

type jInv.Mesh.TensorMesh3D <: AbstractTensorMesh

Fields:

    h1::Vector{Float64}     - cell size in x1 direction
    h2::Vector{Float64}     - cell size in x2 direction
    h3::Vector{Float64}     - cell size in x3 direction
    x0::Vector{Float64}     - origin
    dim::Int                - dimension (dim=3)
    n::Vector{Int64}        - number of cells in each direction
    nc::Int                 - nc total number of cells (nc=prod(n))
    nf::Vector{Int64}       - number of faces
    ne::Vector{Int64}       - number of edges

Persistent Operators:

Operators should not be accessed directly. They will be built, if needed,
when accessing them using specified method. clear!(M) will release all 
memory.

    Div::SparseMatrixCSC    - divergence (faces -> cell-centers)
                              Access via: getDivergenceMatrix(M)
    Grad::SparseMatrixCSC   - gradient (nodal -> edges)
                              Access via: getNodalGradientMatrix(M)
    Curl::SparseMatrixCSC   - curl (edges -> faces)
                              Access via: getCurlMatrix(M)
    Af::SparseMatrixCSC     - face average (faces -> cell-centers)
                              Access via: getFaceAverageMatrix(M)
    Ae::SparseMatrixCSC     - edge average (edges -> cell-centers)
                              Access via: getEdgeAverageMatrix(M)
    An::SparseMatrixCSC     - nodal average (nodes -> cell-centers)
                              Access via: getNodalAverageMatrix(M)
    V::SparseMatrixCSC      - cell volumes (diagonal matrix)
                              Access via: getVolume(M)
    F::SparseMatrixCSC      - face area (diagonal matrix)
                              Access via: getFaceArea(M)
    L::SparseMatrixCSC      - edge length (diagonal matrix)
                              Access via: getLength(M)
    Vi::SparseMatrixCSC     - inverse cell volumes (diagonal matrix)
                              Access via: getVolumeInv(M)
    Fi::SparseMatrixCSC     - inverse face area (diagonal matrix)
                              Access via: getFaceAreaInv(M)
    Li::SparseMatrixCSC     - inverse edge length (diagonal matrix)
                              Access via: getLengthAreaInv(M)
    nLap::SparseMatrixCSC   - nodal Laplacian
                              Access via: getNodalLaplacian(M)

Example: 

h1 = rand(4); h2 = rand(6); h3 = rand(5);
M  = getTensorMesh2D(h1,h2,h3)

source

# jInv.Mesh.getBoundaryNodesMethod.

function jInv.Mesh.getBoundaryNodes(Mesh)

Returns an array tuple (boundary node indices, inner node indices)

Input: 
    Mesh::Abstract Mesh

Output:
    Tuple{Array{Int64,1},Array{Int64,1}}

source

# jInv.Mesh.getEdgeAverageMatrixMethod.

function jInv.Mesh.getEdgeAverageMatrix

Returns Edge-to-CellCenter average matrix from Mesh.Ae. 
Matrix is constructed if Mesh.Ae is empty. 

For 3D Mesh: Ae = [A1 A2 A3]

Input: 
    Mesh::Abstract Mesh

source

# jInv.Mesh.getEdgeIntegralOfPolygonalChainMethod.

    function jInv.Mesh.getEdgeIntegralOfPolygonalChain

    s = getEdgeIntegralPolygonalChain(mesh,polygon)
    s = getEdgeIntegralPolygonalChain(mesh,polygon,normalize)

    Compute the integral of a piecewise linear edge grid basis projected onto
    the edges of a polygonal chain. This function can be used to evaluate the
    source term of a line current carried by the polygonal chain in electromagnetic
    modelling.

    The piecewise linear edge grid basis consists of the functions
      phix_i (i = 1 ... numXEdges)
      phiy_i (i = 1 ... numYEdges)
      phiz_i (i = 1 ... numZEdges)
    where phix_i, phiy_i, phiz_i = 1 at the i-th x/y/z-directed edge
      and phix_i, phiy_i, phiz_i = 0 at all other edges.

    INPUT
    mesh ...... Tensor mesh
    polygon ... vertices of polygonal chain as numVertices x 3 array
    normalize . divide line integral by length of integration path (boolean, optional)

    OUTPUT
    s ......... source vector

    For a closed current loop, specify polygon such that
    polygon[1,:] == polygon[end,:]

source

# jInv.Mesh.getEdgeMassMatrixMethod.

function jInv.Mesh.getEdgeMassMatrix(mesh,sigma)

Returns mass matrix on cell edges, weighted by vector sigma.
Matrix is always constructed. Uses pre-constructed edge averaging
and cell volume matrices if available.

Input: 
    Mesh::Abstract Mesh
       sigma::Vector

Output:
    SparseMatrixCSC{Float64,Int64}

source

# jInv.Mesh.getFaceAverageMatrixMethod.

function jInv.Mesh.getFaceAverageMatrix

Returns Face-to-CellCenter average matrix from Mesh.Af. 
Matrix is constructed if Mesh.Af is empty. 

For 2D Mesh: Af = [A1 A2]
For 3D Mesh: Af = [A1 A2 A3]

Input: 
    Mesh::Abstract Mesh

source

# jInv.Mesh.getFaceMassMatrixMethod.

function jInv.Mesh.getFaceMassMatrix(M,sigma)

Returns face mass matrix, weighted by vector sigma. 
Matrix is always constructed. Uses pre-constructed face averaging
and cell volume matrices if available.

Input: 
    M::Abstract Mesh
       sigma::Vector


Output:
    SparseMatrixCSC{Float64,Int64}

source

# jInv.Mesh.getInterpolationMatrixMethod.

function jInv.Mesh.getInterpolationMatrix

computes bi/trilinear interpolation matrix P for cell-centered data from Mesh1 to Mesh2. If I1 is a cell-centerd discretization of some function on Mesh1 then, its interpolant on Mesh2 is given by

I2 = P*I1

Required Input:

M1::AbstractTensorMesh 
M2::AbstractTensorMesh

Example:

In mesh-decoupling, we use different meshes for the inverse solution
and the different forward problems. 

Mesh2Mesh = getInterpolationMatrix(Minv,Mfor)

source

# jInv.Mesh.getNodalAverageMatrixMethod.

function jInv.Mesh.getNodalAverageMatrix

Returns Nodal-to-CellCenter average matrix from Mesh.An. 
Matrix is constructed if Mesh.An is empty. 

Input: 
    Mesh::Abstract Mesh

source

# jInv.Mesh.getNodalMassMatrixMethod.

function jInv.Mesh.getNodalMassMatrix(M,sigma)

Returns nodal mass matrix, weighted by vector sigma.
Matrix is always constructed. Uses pre-constructed nodal averaging
and cell volume matrices if available.

Input: 
    M::Abstract Mesh
       sigma::Vector

Output:
    SparseMatrixCSC{Float64,Int64}

source

# jInv.Mesh.getRegularMeshMethod.

function jInv.Mesh.getRegularMesh

Constructs regular mesh

Input: 
    domain - physical domain rectangular
    n      - number of cells in each dimension

Examples:
M2D  = getRegularMesh([1.2 2.4 2.2 5.0],[3,4])  
M3D  = getRegularMesh([1.2 2.4 2.2 5.0 0 1],[3,4,7])

source

# jInv.Mesh.getStraightLineCurrentIntegralMethod.

    function jInv.Mesh.getStraightLineCurrentIntegral

    [sx,sy,sz] = getStraightLineCurrentIntegral(hx,hy,hz,ax,ay,az,bx,by,bz)

    Compute integral int(W . J dx^3) in brick of size hx x hy x hz
    where W denotes the 12 local bilinear edge basis functions
    and where J prescribes a unit line current
    between points (ax,ay,az) and (bx,by,bz).

source

# jInv.Mesh.getTensorMesh3DFunction.

function jInv.Mesh.getTensorMesh3D

constructs TensorMesh3D

Required Input:

   h1::Array - cell-sizes in x1 direction
   h2::Array - cell-sizes in x2 direction
   h3::Array - cell-sizes in x3 direction

Optional Input
   x0::Array - origin (default = zeros(3))

source

# jInv.Mesh.getdEdgeMassMatrixMethod.

function jInv.Mesh.getdEdgeMassMatrix(M,v)

Returns directional derivative of edge mass matrix w.r.t. sigma, i.e.,

        d_sigma (M(sigma)*v)

Matrix is always constructed. Uses pre-constructed edge averaging
and cell volume matrices if available.

Input: 
    Mesh::Abstract - Mesh
       v::Vector   - edge vector defining derivative

Output:
    SparseMatrixCSC{Float64,Int64}

source

# jInv.Mesh.getdFaceMassMatrixMethod.

function jInv.Mesh.getdFaceMassMatrix(M,v)

Returns directional derivative of face mass matrix w.r.t sigma, i.e.

    d_sigma (M(sigma)*v)

Matrix is always constructed. Uses pre-constructed face averaging
and cell volume matrices if available.

Input: 
    M::Abstract  - Mesh
       v::Vector -  face vector defining derivative

Output:
    SparseMatrixCSC{Float64,Int64}

source

# jInv.Mesh.getdNodalMassMatrixMethod.

function jInv.Mesh.getdNodalMassMatrix(M,v)

Returns directional derivative of nodal mass matrix w.r.t. sigma, i.e.,

    d_sigma (M(sigma)*v)

Matrix is always constructed. Uses pre-constructed nodal averaging
and cell volume matrices if available.

Input: 
    M::Abstract  - Mesh
       v::Vector - nodal vector defining derivative

Output:
    SparseMatrixCSC{Float64,Int64}

source