API

API

DifferentialOperator(x̄, bc::Tuple{BoundaryCondition, BoundaryCondition}, method::DiscretizationMethod)

Returns a discretized differential operator of length(interiornodes(x̄)) by length(interiornodes(x̄)) matrix under mixed boundary conditions from bc using a discretization method specified by method.

Examples

julia> x̄ = 0:5
0:5

julia> DifferentialOperator(x̄, (Reflecting(), Reflecting()), BackwardFirstDifference())
4×4 LinearAlgebra.Tridiagonal{Float64,Array{Float64,1}}:
  0.0   0.0    ⋅    ⋅
 -1.0   1.0   0.0   ⋅
   ⋅   -1.0   1.0  0.0
   ⋅     ⋅   -1.0  1.0

julia> DifferentialOperator(x̄, (Reflecting(), Reflecting()), ForwardFirstDifference())
4×4 LinearAlgebra.Tridiagonal{Float64,Array{Float64,1}}:
 -1.0   1.0    ⋅    ⋅
  0.0  -1.0   1.0   ⋅
   ⋅    0.0  -1.0  1.0
   ⋅     ⋅    0.0  0.0

julia> DifferentialOperator(x̄, (Reflecting(), Reflecting()), CentralSecondDifference())
4×4 LinearAlgebra.Tridiagonal{Float64,Array{Float64,1}}:
 -1.0   1.0    ⋅     ⋅
  1.0  -2.0   1.0    ⋅
   ⋅    1.0  -2.0   1.0
   ⋅     ⋅    1.0  -1.0

julia> x̄ = 0:5
0:5

julia> DifferentialOperator(x̄, (Mixed(ξ = 1.0), Mixed(ξ = 1.0)), BackwardFirstDifference())
4×4 LinearAlgebra.Tridiagonal{Float64,Array{Float64,1}}:
 Inf     0.0    ⋅    ⋅
  -1.0   1.0   0.0   ⋅
    ⋅   -1.0   1.0  0.0
    ⋅     ⋅   -1.0  1.0

julia> DifferentialOperator(x̄, (Mixed(ξ = 1.0), Mixed(ξ = 1.0)), ForwardFirstDifference())
4×4 LinearAlgebra.Tridiagonal{Float64,Array{Float64,1}}:
 -1.0   1.0    ⋅     ⋅
  0.0  -1.0   1.0    ⋅
   ⋅    0.0  -1.0   1.0
   ⋅     ⋅    0.0  -0.5

julia> DifferentialOperator(x̄, (Mixed(ξ = 1.0), Mixed(ξ = 1.0)), CentralSecondDifference())
4×4 LinearAlgebra.Tridiagonal{Float64,Array{Float64,1}}:
 -Inf     1.0    ⋅     ⋅
    1.0  -2.0   1.0    ⋅
     ⋅    1.0  -2.0   1.0
     ⋅     ⋅    1.0  -1.5
source
ExtensionDifferentialOperator(x̄, method::DiscretizationMethod)

Returns a discretized differential operator of length(interiornodes(x̄)) by length(x̄) matrix under no boundary condition using a discretization method specified by method.

Examples

julia> x̄ = 0:5
0:5

julia> ExtensionDifferentialOperator(x̄, BackwardFirstDifference())
4×6 SparseArrays.SparseMatrixCSC{Float64,Int64} with 8 stored entries:
  [1, 1]  =  -1.0
  [1, 2]  =  1.0
  [2, 2]  =  -1.0
  [2, 3]  =  1.0
  [3, 3]  =  -1.0
  [3, 4]  =  1.0
  [4, 4]  =  -1.0
  [4, 5]  =  1.0

julia> ExtensionDifferentialOperator(x̄, ForwardFirstDifference())
4×6 SparseArrays.SparseMatrixCSC{Float64,Int64} with 8 stored entries:
  [1, 2]  =  -1.0
  [1, 3]  =  1.0
  [2, 3]  =  -1.0
  [2, 4]  =  1.0
  [3, 4]  =  -1.0
  [3, 5]  =  1.0
  [4, 5]  =  -1.0
  [4, 6]  =  1.0

julia> ExtensionDifferentialOperator(x̄, CentralSecondDifference())
4×6 SparseArrays.SparseMatrixCSC{Float64,Int64} with 12 stored entries:
  [1, 1]  =  1.0
  [1, 2]  =  -2.0
  [2, 2]  =  1.0
  [1, 3]  =  1.0
  [2, 3]  =  -2.0
  [3, 3]  =  1.0
  [2, 4]  =  1.0
  [3, 4]  =  -2.0
  [4, 4]  =  1.0
  [3, 5]  =  1.0
  [4, 5]  =  -2.0
  [4, 6]  =  1.0

julia> ExtensionDifferentialOperator(x̄, JumpProcess(x̄, -1.0))
4×6 BandedMatrices.BandedMatrix{Float64,Array{Float64,2},Base.OneTo{Int64}}:
 0.0  0.0    ⋅     ⋅     ⋅    ⋅
  ⋅   1.0  -1.0    ⋅     ⋅    ⋅
  ⋅    ⋅    1.0  -1.0    ⋅    ⋅
  ⋅    ⋅     ⋅    1.0  -1.0   ⋅
source
L₀affine(x̄, b, bc::Tuple{BoundaryCondition, BoundaryCondition})

Returns a vector for affine operator given an extended grid `x̄` and boundary conditions `bc`

to solve the system L₀ * v(̄x) = b on interior nodes.

Returns a `length(interiornodes(x̄))`-length vector such that solving

L₀ * v(x̄) = b on interior nodes x under the boundary condition given by bc is identical with L₀bc * v(x) = b + L₀affine(x̄, bc) The first element of bc is applied to the lower bound, and second element of bc to the upper.

Examples

julia> x̄ = 0:5
0:5

julia> b = 1:3
1:3

julia> L₀affine(x̄, b, (Reflecting(), Reflecting()))
4-element Array{Float64,1}:
 0.0
 0.0
 0.0
 0.0

julia> L₀affine(x̄, b, (NonhomogeneousAbsorbing(1.0), Reflecting()))
4-element Array{Float64,1}:
 0.0
 0.0
 0.0
 0.0

julia> L₀affine(x̄, b, (NonhomogeneousAbsorbing(1.0), NonhomogeneousAbsorbing(2.0)))
4-element Array{Float64,1}:
 0.0
 0.0
 0.0
 0.0
source
L₁₊(x̄)

Returns a discretized first-order differential operator of length(interiornodes(x̄)) by length(x̄) matrix using forward difference under no boundary condition.

Examples

julia> x̄ = 0:5
0:5

julia> Array(L₁₊(x̄))
4×6 Array{Float64,2}:
 0.0  -1.0   1.0   0.0   0.0  0.0
 0.0   0.0  -1.0   1.0   0.0  0.0
 0.0   0.0   0.0  -1.0   1.0  0.0
 0.0   0.0   0.0   0.0  -1.0  1.0
source
L₁₊affine(x̄, b, bc::Tuple{BoundaryCondition, BoundaryCondition})

Returns a vector for affine operator given an extended grid `x̄` and boundary conditions `bc`

to solve the system L₁₊ * v(̄x) = b on interior nodes.

Returns a `length(interiornodes(x̄))`-length vector such that solving

L₁₊ * v(x̄) = b on interior nodes x under the boundary condition given by bc is identical with L₁₊bc * v(x) = b + L₁₊affine(x̄, bc) The first element of bc is applied to the lower bound, and second element of bc to the upper.

Examples

julia> x̄ = 0:5
0:5

julia> b = 1:3
1:3

julia> L₁₊affine(x̄, b, (Reflecting(), Reflecting()))
4-element Array{Float64,1}:
 0.0
 0.0
 0.0
 0.0

julia> L₁₊affine(x̄, b, (NonhomogeneousAbsorbing(1.0), Reflecting()))
4-element Array{Float64,1}:
 0.0
 0.0
 0.0
 0.0

julia> L₁₊affine(x̄, b, (NonhomogeneousAbsorbing(1.0), NonhomogeneousAbsorbing(2.0)))
4-element Array{Float64,1}:
  0.0
  0.0
  0.0
 -2.0
source
L₁₊bc(x̄, bc::Tuple{BoundaryCondition, BoundaryCondition})

Returns a discretized first-order differential operator of length(interiornodes(x̄)) by length(interiornodes(x̄)) matrix using forward difference under boundary conditions specified by bc.

The first element of bc is applied to the lower bound, and second element of bc to the upper.

Examples

julia> x̄ = 0:5
0:5

julia> L₁₊bc(x̄, (Reflecting(), Reflecting()))
4×4 LinearAlgebra.Tridiagonal{Float64,Array{Float64,1}}:
 -1.0   1.0    ⋅    ⋅
  0.0  -1.0   1.0   ⋅
   ⋅    0.0  -1.0  1.0
   ⋅     ⋅    0.0  0.0
source
L₁₋(x̄)

Returns a discretized first-order differential operator of length(interiornodes(x̄)) by length(x̄) matrix using backward difference under no boundary condition.

Examples

julia> x̄ = 1:3
1:3

julia> Array(L₁₋(x̄))
1×3 Array{Float64,2}:
 -1.0  1.0  0.0
source
L₁₋affine(x̄, b, bc::Tuple{BoundaryCondition, BoundaryCondition})

Returns a vector for affine operator given an extended grid `x̄` and boundary conditions `bc`

to solve the system L₁₋ * v(̄x) = b on interior nodes.

Returns a `length(interiornodes(x̄))`-length vector such that solving

L₁₋ * v(x̄) = b on interior nodes x under the boundary condition given by bc is identical with L₁₋bc * v(x) = b + L₁₋affine(x̄, bc). The first element of bc is applied to the lower bound, and second element of bc to the upper.

Examples

julia> x̄ = 0:5
0:5

julia> b = 1:3
1:3

julia> L₁₋affine(x̄, b, (Reflecting(), Reflecting()))
4-element Array{Float64,1}:
 0.0
 0.0
 0.0
 0.0

julia> L₁₋affine(x̄, b, (NonhomogeneousAbsorbing(1.0), Reflecting()))
4-element Array{Float64,1}:
 1.0
 0.0
 0.0
 0.0

julia> L₁₋affine(x̄, b, (NonhomogeneousAbsorbing(1.0), NonhomogeneousAbsorbing(2.0)))
4-element Array{Float64,1}:
 1.0
 0.0
 0.0
 0.0
source
L₁₋bc(x̄, bc::Tuple{BoundaryCondition, BoundaryCondition})

Returns a discretized first-order differential operator of length(interiornodes(x̄)) by length(interiornodes(x̄)) matrix using backward difference under boundary conditions specified by bc.

The first element of bc is applied to the lower bound, and second element of bc to the upper.

Examples

julia> x̄ = 0:5
0:5

julia> L₁₋bc(x̄, (Reflecting(), Reflecting()))
4×4 LinearAlgebra.Tridiagonal{Float64,Array{Float64,1}}:
  0.0   0.0    ⋅    ⋅
 -1.0   1.0   0.0   ⋅
   ⋅   -1.0   1.0  0.0
   ⋅     ⋅   -1.0  1.0
source
L₂(x̄)

Returns a discretized second-order differential operator of length(interiornodes(x̄)) by length(x̄) matrix using central difference under no boundary condition.

Examples

julia> x̄ = 0:5
0:5 

julia> Array(L₂(x̄))
4×6 Array{Float64,2}:
 1.0  -2.0   1.0   0.0   0.0  0.0
 0.0   1.0  -2.0   1.0   0.0  0.0
 0.0   0.0   1.0  -2.0   1.0  0.0
 0.0   0.0   0.0   1.0  -2.0  1.0
source
L₂affine(x̄, bc::Tuple{BoundaryCondition, BoundaryCondition})

Returns a vector for affine operator given an extended grid `x̄` and boundary conditions `bc`

to solve the system L₂ * v(̄x) = b on interior nodes.

Returns a `length(interiornodes(x̄))`-length vector such that solving

L₂ * v(x̄) = b on interior nodes x under the boundary condition given by bc is identical with L₂bc * v(x) = b + L₂affine(x̄, bc) The first element of bc is applied to the lower bound, and second element of bc to the upper.

Examples

julia> x̄ = 0:5
0:5

julia> b = 1:3
1:3

julia> L₂affine(x̄, b, (Reflecting(), Reflecting()))
4-element Array{Float64,1}:
 0.0
 0.0
 0.0
 0.0

julia> L₂affine(x̄, b, (NonhomogeneousAbsorbing(1.0), Reflecting()))
4-element Array{Float64,1}:
 -1.0
  0.0
  0.0
  0.0

julia> L₂affine(x̄, b, (NonhomogeneousAbsorbing(1.0), NonhomogeneousAbsorbing(2.0)))
4-element Array{Float64,1}:
 -1.0
  0.0
  0.0
 -2.0
source
L₂bc(x̄, bc::Tuple{BoundaryCondition, BoundaryCondition})

Returns a discretized second-order differential operator of length(interiornodes(x̄)) by length(interiornodes(x̄)) matrix using central difference under boundary conditions specified by bc.

The first element of bc is applied to the lower bound, and second element of bc to the upper.

Examples

julia> x̄ = 0:5
0:5

julia> L₂bc(x̄, (Reflecting(), Reflecting()))
4×4 LinearAlgebra.Tridiagonal{Float64,Array{Float64,1}}:
 -1.0   1.0    ⋅     ⋅
  1.0  -2.0   1.0    ⋅
   ⋅    1.0  -2.0   1.0
   ⋅     ⋅    1.0  -1.0
source
Lₙ(x̄, method)

Returns a discretized jump process operator of length(interiornodes(x̄)) by length(x̄) matrix specified by method

Examples

julia> x̄ = 0:5
0:5 

julia> Lₙ(x̄, JumpProcess(x̄, -1.0))
4×6 BandedMatrices.BandedMatrix{Float64,Array{Float64,2},Base.OneTo{Int64}}:
 0.0  0.0    ⋅     ⋅     ⋅    ⋅
  ⋅   1.0  -1.0    ⋅     ⋅    ⋅
  ⋅    ⋅    1.0  -1.0    ⋅    ⋅
  ⋅    ⋅     ⋅    1.0  -1.0   ⋅
source
Lₙbc(x̄, method, (Absorbing(), Absorbing()))

Returns a discretized jump process operator of length(interiornodes(x̄)) by length(interiornodes(x̄)) matrix specified by method under boundary conditions specified by bc.

Examples

julia> x̄ = 0:5
0:5 

julia> Lₙbc(x̄, (Absorbing(), Absorbing()), JumpProcess(x̄, -1.0))
4×4 BandedMatrices.BandedMatrix{Float64,Array{Float64,2},Base.OneTo{Int64}}:
 0.0   0.0    ⋅     ⋅
 1.0  -1.0   0.0    ⋅
  ⋅    1.0  -1.0   0.0
  ⋅     ⋅    1.0  -1.0
source
extrapolatetoboundary(v, x̄, bc::Tuple{BoundaryCondition, BoundaryCondition})

Returns a length(x̄)-vector whose 2:(length(x̄)-1) elements are v, the first and last element are extrapolated v on the boundaries of according to boundary conditions bc given.

The first element of bc is applied to the lower bound, and second element of bc to the upper.

Examples

julia> x̄ = -2:2
-2:2

julia> x = interiornodes(x̄)
-1:1

julia> v = (x -> x^2).(x)
3-element Array{Int64,1}:
 1
 0
 1

julia> extrapolatetoboundary(v, x̄, (Absorbing(), Absorbing()))
5-element Array{Int64,1}:
 0
 1
 0
 1
 0

julia> extrapolatetoboundary(v, x̄, (Absorbing(), Reflecting()))
5-element Array{Int64,1}:
 0
 1
 0
 1
 1

julia> extrapolatetoboundary(v, x̄, (Mixed(ξ = 3.0), Reflecting()))
5-element Array{Float64,1}:
 -0.5
  1.0
  0.0
  1.0
  1.0
  
source
interiornodes(x̄, bc)

Returns an interior grid corresponding to the boundary condition bc given extended grid .

julia> x̄ = 0:5
0:5

julia> interiornodes(x̄, (Reflecting(), Reflecting()))
1:4

julia> x̄ = [1.0; 1.5; 1.7]
3-element Array{Float64,1}:
 1.0
 1.5
 1.7

julia> interiornodes(x̄, (Mixed(ξ = 1.0), Mixed(ξ = 1.0)))
1-element Array{Float64,1}:
 1.5
source
interiornodes(x̄)

Returns an interior grid of length length(x̄)-2 given extended grid .

julia> x̄ = 0:5
0:5

julia> interiornodes(x̄)
1:4

julia> x̄ = [1.0; 1.5; 1.7]
3-element Array{Float64,1}:
 1.0
 1.5
 1.7

julia> interiornodes(x̄)
1-element Array{Float64,1}:
 1.5
source
jointoperator_bc(operators, Q::Array)

Returns a discretized operator that solves systems of differential equations defined by operators with transitions by Q where operators is an N-length collection of discretized operators with boundary conditions applied and Q is N by N intensity matrix.

Examples

julia> x̄ = 0:5
0:5

julia> L1bc = L₁₋bc(x̄, (Reflecting(), Reflecting()))
4×4 LinearAlgebra.Tridiagonal{Float64,Array{Float64,1}}:
  0.0   0.0    ⋅    ⋅
 -1.0   1.0   0.0   ⋅
   ⋅   -1.0   1.0  0.0
   ⋅     ⋅   -1.0  1.0

julia> L2bc = L₂bc(x̄, (Reflecting(), Reflecting()))
4×4 LinearAlgebra.Tridiagonal{Float64,Array{Float64,1}}:
 -1.0   1.0    ⋅     ⋅
  1.0  -2.0   1.0    ⋅
   ⋅    1.0  -2.0   1.0
   ⋅     ⋅    1.0  -1.0

julia> Q = [-0.5 0.5; 0.3 -0.3]
2×2 Array{Float64,2}:
 -0.5   0.5
  0.3  -0.3

julia> jointoperator_bc((L1bc, L2bc), Q)
2×2-blocked 8×8 BlockBandedMatrices.BandedBlockBandedMatrix{Float64,BlockArrays.PseudoBlockArray{Float64,2,Array{Float64,2},BlockArrays.BlockSizes{2,Tuple{Array{Int64,1},Array{Int64,1}}}}}:
 -0.5   0.0    ⋅    ⋅   │   0.5   0.0    ⋅     ⋅
 -1.0   0.5   0.0   ⋅   │   0.0   0.5   0.0    ⋅
   ⋅   -1.0   0.5  0.0  │    ⋅    0.0   0.5   0.0
   ⋅     ⋅   -1.0  0.5  │    ⋅     ⋅    0.0   0.5
 ───────────────────────┼────────────────────────
  0.3   0.0    ⋅    ⋅   │  -1.3   1.0    ⋅     ⋅
  0.0   0.3   0.0   ⋅   │   1.0  -2.3   1.0    ⋅
   ⋅    0.0   0.3  0.0  │    ⋅    1.0  -2.3   1.0
   ⋅     ⋅    0.0  0.3  │    ⋅     ⋅    1.0  -1.3
source

JumpProcess(x̄, jumps::AbstractArray, truncate = (:interior, :interior))

Returns a DiscretizationMethod object that can be used to construct a discretized operator jump process.

jumps is a (length(x̄)-2)-vector whose ith element is an integer that represents a jump size in index and direction (by sign) from ith element of interiornodes(x̄) and the first and second elements of truncate represent truncation location for the lower bound and upper bound when the jump is out of the truncated boundary (:interior for the first/last element of interiornodes(̄x) and :exterior for the first/last element of ). The default parameter is (:interior, :interior).

Examples

julia> x̄ = 0:5
0:5

julia> JumpProcess(x̄, [-1; -1; -1; -1], (:interior, :interior))
JumpProcess{Array{Int64,1}}([0, -1, -1, -1])
source

JumpProcess(x̄, jumpf::Function, truncate = (:interior, :interior))

Returns a DiscretizationMethod object that can be used to construct a discretized operator jump process.

jumpf is a function that takes an element of interiornodes(x̄) and returns a scalar that represents the corresponding nominal jump size and direction (by sign) and the first and second elements of truncate represent truncation location for the lower bound and upper bound when the jump is out of the truncated boundary (:interior for the first/last element of interiornodes(̄x) and :exterior for the first/last element of ). The default parameter is (:interior, :interior). The code uses nearest-neighbour rule to determine the indices of destinations according to `jumpf.

Examples

julia> x̄ = 0:5
0:5

julia> jumpf(x) = -1.4
jumpf (generic function with 1 method)

julia> JumpProcess(x̄, jumpf)
JumpProcess{Array{Int64,1}}([0, -1, -1, -1])
source

JumpProcess(x̄, uniformjumpsize::Real, truncate = (:interior, :interior))

Returns a DiscretizationMethod object that can be used to construct a discretized operator jump process.

uniform_jump_size is a scalar in Real that represents a nominal jump size and direction (by sign) from all elements of interiornodes(x̄) and the first and second elements of truncate represent truncation location for the lower bound and upper bound when the jump is out of the truncated boundary (:interior for the first/last element of interiornodes(̄x) and :exterior for the first/last element of ). The default parameter is (:interior, :interior). The code uses nearest-neighbour rule to determine the indices of destinations according to `jumpf.

Examples

julia> x̄ = 0:5
0:5

julia> JumpProcess(x̄, -1.4)
JumpProcess{Array{Int64,1}}([0, -1, -1, -1])
source

JumpProcess(x̄, uniform_jump::AbstractArray, truncate = (:interior, :interior))

Returns a DiscretizationMethod object that can be used to construct a discretized operator jump process.

uniform_jump is a scalar Int64 that represents a jump size in index and direction (by sign) from all elements of interiornodes(x̄) and the first and second elements of truncate represent truncation location for the lower bound and upper bound when the jump is out of the truncated boundary (:interior for the first/last element of interiornodes(̄x) and :exterior for the first/last element of ). The default parameter is (:interior, :interior).

Examples

julia> x̄ = 0:5
0:5

julia> JumpProcess(x̄, -1)
JumpProcess{Array{Int64,1}}([0, -1, -1, -1])
source