FiniteVolumeMethod.ConstrainedConstant
Constrained

Instance of a ConditionType used for declaring that an edge has a Constrained condition. When an edge has this condition associated with it, it will be treated as any normal edge and no boundary condition will be applied. With this condition, it is assumed that you will later setup your problem as a differential-algebraic equation (DAE) and provide the appropriate constraints. See the docs for some examples.

When you provide a Constrained condition, for certain technical reasons you do still need to provide a function that corresponds to it in the function list provided to BoundaryConditions. For this function, any function will work, e.g. sin - it will not be called. The proper constraint function is to be provided after-the-fact, as explained in the docs.

FiniteVolumeMethod.DirichletConstant
Dirichlet

Instance of a ConditionType used for declaring that an edge, or a point, has a Dirichlet boundary condition. Dirichlet conditions take the form

\[u(x, y, t) = a(x, y, t, u).\]

When providing a Dirichlet condition, the function you provide takes the form

a(x, y, t, u, p)

where (x, y) is the point, t is the current time, and u is the solution at the point (x, y) at time t, as above, with an extra argument p for additional parameters.

FiniteVolumeMethod.DudtConstant
Dudt

Instance of a ConditionType used for declaring that an edge, or a point, has a Dudt-type boundary condition. Dudt-type conditions take the form

\[\dv{u(x, y, t)}{t} = a(x, y, t, u).\]

When providing a Dudt condition, the function you provide takes the form

a(x, y, t, u, p)

where (x, y) is the point, t is the current time, and u is the solution at the point (x, y) at time t, as above, with an extra argument p for additional parameters.

FiniteVolumeMethod.NeumannConstant
Neumann

Instance of a ConditionType used for declaring that an edge has a Neumann condition. Neumann conditions take the form

\[\vb q(x, y, t) \vdot \vu n_\sigma = a(x, y, t, u)\]

where $\vb q$ is the flux function and $\vu n_\sigma$ is the outward unit normal vector field on the associated edge (meaning, for example, the normal vector to an edge pq would point to the right of pq).

When providing a Neumann condition, the function you provide takes the form

a(x, y, t, u, p)

where (x, y) is the point, t is the current time, and u is the solution at the point (x, y) at time t, as above, with an extra argument p for additional parameters.

FiniteVolumeMethod.AbstractFVMTemplateType
abstract type AbstractFVMTemplate <: AbstractFVMProblem

An abstract type that defines some specific problems. These problems are those that could be defined directly using FVMProblems, but are common enough that (1) it is useful to have them defined here, and (2) it is useful to have them defined in a way that is more efficient than with a default implementation (e.g. exploiting linearity). The problems are all defined as subtypes of a common abstract type, namely, AbstractFVMTemplate (the home of this docstring), which itself is a subtype of AbstractFVMProblem.

To understand how to make use of these specific problems, either see the docstring for each problem, or see the "Solvers for Specific Problems, and Writing Your Own" section of the docs.

To see the full list of problems, do

julia> using FiniteVolumeMethod

julia> subtypes(FiniteVolumeMethod.AbstractFVMTemplate)
5-element Vector{Any}:
 DiffusionEquation
 LaplacesEquation
 LinearReactionDiffusionEquation
 MeanExitTimeProblem
 PoissonsEquation

The constructor for each problem is defined in its docstring. Note that all the problems above are exported.

These problems can all be solved using the standard solve interface from DifferentialEquations.jl, just like for FVMProblems. The only exception is for steady state problems, in which case the solve interface is still used, except the interface is from LinearSolve.jl.

FiniteVolumeMethod.BoundaryConditionsType
BoundaryConditions(mesh::FVMGeometry, functions, conditions; parameters=nothing)

This is a constructor for the BoundaryConditions struct, which holds the boundary conditions for the PDE. See also Conditions (which FVMProblem wraps this into), ConditionType, and InternalConditions.

Arguments

  • mesh::FVMGeometry

The mesh on which the PDE is defined.

  • functions::Union{<:Tuple,<:Function}

The functions that define the boundary conditions. The ith function should correspond to the part of the boundary of the mesh corresponding to the ith boundary index, as defined in DelaunayTriangulation.jl.

  • conditions::Union{<:Tuple,<:ConditionType}

The classification for the boundary condition type corresponding to each boundary index as above. See ConditionType for possible conditions - should be one of Neumann, Dudt, Dirichlet, or Constrained.

Keyword Arguments

  • parameters=ntuple(_ -> nothing, length(functions))

The parameters for the functions, with parameters[i] giving the argument p in functions[i].

Outputs

The returned value is the corresponding BoundaryConditions struct.

FiniteVolumeMethod.ConditionsType
Conditions{F} <: AbstractConditions

This is a struct that holds the boundary and internal conditions for the PDE. The constructor is

Conditions(mesh::FVMGeometry, bc::BoundaryConditions, ic::InternalConditions=InternalConditions())

The fields are:

Fields

  • neumann_conditions::Dict{NTuple{2,Int},Int}

A Dict that stores all Neumann edges (u, v) as keys, with keys mapping to indices idx that refer to the corresponding condition function and parameters in functions.

  • constrained_conditions::Dict{NTuple{2,Int},Int}

A Dict that stores all Constrained edges (u, v) as keys, with keys mapping to indices idx that refer to the corresponding condition function and parameters in functions.

  • dirichlet_conditions::Dict{Int,Int}

A Dict that stores all Dirichlet points u as keys, with keys mapping to indices idx that refer to the corresponding condition function and parameters in functions.

  • dudt_conditions::Dict{Int,Int}

A Dict that stores all Dudt points u as keys, with keys mapping to indices idx that refer to the corresponding condition function and parameters in functions.

  • functions::F<:Tuple

The functions that define the boundary and internal conditions. The ith function should correspond to the part of the boundary of the mesh corresponding to the ith boundary index, as defined in DelaunayTriangulation.jl. The ith function is stored as a ParametrisedFunction.

FiniteVolumeMethod.DiffusionEquationType
DiffusionEquation <: AbstractFVMTemplate

A struct for defining a problem representing a diffusion equation:

\[\pdv{u}{t} = \div\left[D(\vb x)\grad u\right]\]

inside a domain $\Omega$.

You can solve this problem using solve.

Warning

The solution to this problem will have an extra component added to it. The original solution will be inside sol[begin:end-1, :], where sol is the solution returned by solve.

Constructor

DiffusionEquation(mesh::FVMGeometry,
    BCs::BoundaryConditions,
    ICs::InternalConditions=InternalConditions();
    diffusion_function,
    diffusion_parameters=nothing,
    initial_condition,
    initial_time=0.0,
    final_time,
    kwargs...)

Arguments

  • mesh::FVMGeometry: The FVMGeometry.
  • BCs::BoundaryConditions: The BoundaryConditions. For these boundary conditions, all functions should still be of the form (x, y, t, u, p) -> Number, but the t and u arguments should be unused as they will be replaced with nothing.
  • ICs::InternalConditions=InternalConditions(): The InternalConditions. For these internal conditions, all functions should still be of the form (x, y, t, u, p) -> Number, but the t and u arguments should be unused as they will be replaced with nothing.

Keyword Arguments

  • diffusion_function: The diffusion function. Should be of the form (x, y, p) -> Number, where p = diffusion_parameters below.
  • diffusion_parameters=nothing: The argument p in diffusion_function.
  • initial_condition: The initial condition.
  • initial_time=0.0: The initial time.
  • final_time: The final time.
  • kwargs...: Any other keyword arguments are passed to the ODEProblem (from DifferentialEquations.jl) that represents the problem.

Fields

The struct has extra fields in addition to the arguments above:

  • A: This is a sparse matrix A so that du/dt = Au + b.
  • b: The b above.
  • Aop: The MatrixOperator that represents the system so that du/dt = Aop*u (with u padded with an extra component since A is now inside Aop).
  • problem: The ODEProblem that represents the problem. This is the problem that is solved when you call solve on the struct.
FiniteVolumeMethod.FVMGeometryType
FVMGeometry(tri::Triangulation)

This is a constructor for the FVMGeometry struct, which holds the mesh and associated data for the PDE.

Note

It is assumed that all vertices in tri are in the triangulation, meaning v is in tri for each v in DelaunayTriangulation.each_point_index(tri).

Fields

  • triangulation: The underlying Triangulation from DelaunayTriangulation.jl.
  • triangulation_statistics: The statistics of the triangulation.
  • cv_volumes::Vector{Float64}: A Vector of the volumes of each control volume.
  • triangle_props::Dict{NTuple{3,Int},TriangleProperties}: A Dict mapping the indices of each triangle to its [TriangleProperties].
FiniteVolumeMethod.FVMProblemType
FVMProblem(mesh, boundary_conditions[, internal_conditions];
    diffusion_function=nothing,
    diffusion_parameters=nothing,
    source_function=nothing,
    source_parameters=nothing,
    flux_function=nothing,
    flux_parameters=nothing,
    initial_condition,
    initial_time=0.0,
    final_time)

Constructs an FVMProblem. See also FVMSystem and SteadyFVMProblem.

Arguments

  • mesh::FVMGeometry

The mesh on which the PDE is defined, given as a FVMGeometry.

  • boundary_conditions::BoundaryConditions

The boundary conditions for the PDE, given as a BoundaryConditions.

  • internal_conditions::InternalConditions=InternalConditions()

The internal conditions for the PDE, given as an InternalConditions. This argument is optional.

Keyword Arguments

  • diffusion_function=nothing

If isnothing(flux_function), then this can be provided to give the diffusion-source formulation. See also construct_flux_function. Should be of the form D(x, y, t, u, p).

  • diffusion_parameters=nothing

The argument p for diffusion_function.

  • source_function=(x, y, t, u, p) -> zero(u)

The source term, given in the form S(x, y, t, u, p).

  • source_parameters=nothing

The argument p for source_function.

  • flux_function=nothing

The flux function, given in the form q(x, y, t, α, β, γ, p) ↦ (qx, qy), where (qx, qy) is the flux vector, (α, β, γ) are the shape function coefficients for estimating u ≈ αx+βy+γ. If isnothing(flux_function), then diffusion_function is used instead to construct the function.

  • flux_parameters=nothing

The argument p for flux_function.

  • initial_condition

The initial condition, with initial_condition[i] the initial value at the ith node of the mesh.

  • initial_time=0.0

The initial time.

  • final_time

The final time.

Outputs

The returned value is the corresponding FVMProblem struct. You can then solve the problem using solve(::Union{FVMProblem,FVMSystem}, ::Any; kwargs...) from DifferentialEquations.jl.

FiniteVolumeMethod.FVMSystemType
FVMSystem(prob1, prob2, ..., probN)

Constructs a representation for a system of PDEs, where each probi is a FVMProblem for the ith component of the system.

For these FVMProblems, the functions involved, such as the condition functions, should all be defined so that the u argument assumes the form u = (u₁, u₂, ..., uN) (both Tuples and Vectors will be passed), where uᵢ is the solution for the ith component of the system. For the flux functions, which for a FVMProblem takes the form

q(x, y, t, α, β, γ, p) ↦ (qx, qy),

the same form is used, except α, β, γ are all Tuples so that α[i]*x + β[i]*y + γ[i] is the approximation to uᵢ.

Providing default flux functions

Due to this difference in flux functions, and the need to provide α, β, and γ to the flux function, for FVMSystems you need to provide a flux function rather than a diffusion function. If you do provide a diffusion function, it will error when you try to solve the problem.

This problem is solved in the same way as a FVMProblem, except the problem is defined such that the solution returns a matrix at each time, where the (j, i)th component corresponds to the solution at the ith node for the jth component.

See also FVMProblem and SteadyFVMProblem.

Note

To construct a steady-state problem for an FVMSystem, you need to apply SteadyFVMProblem to the system rather than first applying it to each individual FVMProblem in the system.

FiniteVolumeMethod.InternalConditionsType
InternalConditions(functions=();
    dirichlet_nodes::Dict{Int,Int}=Dict{Int,Int}(),
    dudt_nodes::Dict{Int,Int}=Dict{Int,Int}(),
    parameters::Tuple=ntuple(_ -> nothing, length(functions)))

This is a constructor for the InternalConditions struct, which holds the internal conditions for the PDE. See also Conditions (which FVMProblem wraps this into), ConditionType, and BoundaryConditions.

Arguments

  • functions::Union{<:Tuple,<:Function}=()

The functions that define the internal conditions. These are the functions referred to in edge_conditions and point_conditions.

Keyword Arguments

  • dirichlet_nodes::Dict{Int,Int}=Dict{Int,Int}()

A Dict that stores all Dirichlet points u as keys, with keys mapping to indices idx that refer to the corresponding condition function and parameters in functions and parameters.

  • dudt_nodes::Dict{Int,Int}=Dict{Int,Int}()

A Dict that stores all Dudt points u as keys, with keys mapping to indices idx that refer to the corresponding condition function and parameters in functions and parameters.

  • parameters::Tuple=ntuple(_ -> nothing, length(functions))

The parameters for the functions, with parameters[i] giving the argument p in functions[i].

Outputs

The returned value is the corresponding InternalConditions struct.

Note

When the internal conditions get merged with the boundary conditions, any internal conditions that are placed onto the boundary will be replaced with the boundary condition at that point on the boundary.

FiniteVolumeMethod.LaplacesEquationType
LaplacesEquation

A struct for defining a problem representing a (generalised) Laplace's equation:

\[\div[D(\vb x)\grad u] = 0\]

inside a domain $\Omega$. See also PoissonsEquation.

You can solve this problem using solve.

Constructor

LaplacesEquation(mesh::FVMGeometry,
    BCs::BoundaryConditions,
    ICs::InternalConditions=InternalConditions();
    diffusion_function=(x,y,p)->1.0,
    diffusion_parameters=nothing,
    kwargs...)

Arguments

  • mesh::FVMGeometry: The FVMGeometry.
  • BCs::BoundaryConditions: The BoundaryConditions. For these boundary conditions, all functions should still be of the form (x, y, t, u, p) -> Number, but the t and u arguments should be unused as they will be replaced with nothing.
  • ICs::InternalConditions=InternalConditions(): The InternalConditions. For these internal conditions, all functions should still be of the form (x, y, t, u, p) -> Number, but the t and u arguments should be unused as they will be replaced with nothing.

Keyword Arguments

  • diffusion_function=(x,y,p)->1.0: The diffusion function. Should be of the form (x, y, p) -> Number, where p = diffusion_parameters below.
  • diffusion_parameters=nothing: The argument p in diffusion_function.
  • kwargs...: Any other keyword arguments are passed to the LinearProblem (from LinearSolve.jl) that represents the problem.

Fields

The struct has extra fields in addition to the arguments above:

  • A: This is a sparse matrix A so that Au = b.
  • b: The b above.
  • problem: The LinearProblem that represents the problem. This is the problem that is solved when you call solve on the struct.
FiniteVolumeMethod.LinearReactionDiffusionEquationType
LinearReactionDiffusionEquation

A struct for defining a problem representing a linear reaction-diffusion equation:

\[\pdv{u}{t} = \div\left[D(\vb x)\grad u\right] + f(\vb x)u\]

inside a domain $\Omega$.

You can solve this problem using solve.

Warning

The solution to this problem will have an extra component added to it. The original solution will be inside sol[begin:end-1, :], where sol is the solution returned by solve.

Constructor

LinearReactionDiffusionEquation(mesh::FVMGeometry,
    BCs::BoundaryConditions,
    ICs::InternalConditions=InternalConditions();
    diffusion_function,
    diffusion_parameters=nothing,
    source_function,
    source_parameters=nothing,
    initial_condition,
    initial_time=0.0,
    final_time,
    kwargs...)

Arguments

  • mesh::FVMGeometry: The FVMGeometry.
  • BCs::BoundaryConditions: The BoundaryConditions. For these boundary conditions, all functions should still be of the form (x, y, t, u, p) -> Number, but the t and u arguments should be unused as they will be replaced with nothing.
  • ICs::InternalConditions=InternalConditions(): The InternalConditions. For these internal conditions, all functions should still be of the form (x, y, t, u, p) -> Number, but the t and u arguments should be unused as they will be replaced with nothing.

Keyword Arguments

  • diffusion_function: The diffusion function. Should be of the form (x, y, p) -> Number, where p = diffusion_parameters below.
  • diffusion_parameters=nothing: The argument p in diffusion_function.
  • source_function: The source function. Should be of the form (x, y, p) -> Number, where p = source_parameters below.
  • source_parameters=nothing: The argument p in source_function.
  • initial_condition: The initial condition.
  • initial_time=0.0: The initial time.
  • final_time: The final time.
  • kwargs...: Any other keyword arguments are passed to the ODEProblem (from DifferentialEquations.jl) that represents the problem.

Fields

The struct has extra fields in addition to the arguments above:

  • A: This is a sparse matrix A so that du/dt = Au + b.
  • b: The b above.
  • Aop: The MatrixOperator that represents the system so that du/dt = Aop*u (with u padded with an extra component since A is now inside Aop).
  • problem: The ODEProblem that represents the problem. This is the problem that is solved when you call solve on the struct.
FiniteVolumeMethod.MeanExitTimeProblemType
MeanExitTimeProblem

A struct for defining a problem representing a mean exit time problem:

\[\div\left[D(\vb x)\grad T\right] =-1\]

inside a domain $\Omega$. This problem is a special case of PoissonsEquation, but is defined separately since it is common enough to warrant its own definition; MeanExitTimeProblem is constructed using PoissonsEquation.

You can solve this problem using solve.

Constructor

MeanExitTimeProblem(mesh::FVMGeometry,
    BCs::BoundaryConditions,
    ICs::InternalConditions=InternalConditions();
    diffusion_function,
    diffusion_parameters=nothing,
    kwargs...)

Arguments

The functions for BCs and ICs are not used. Whenever a Neumann condition is encountered, or a Dirichlet condition, it is assumed that the conditon is homogeneous. If any of the conditions are Dudt or Constrained types, then an error is thrown.

Keyword Arguments

  • diffusion_function: The diffusion function. Should be of the form (x, y, p) -> Number, where p = diffusion_parameters below.
  • diffusion_parameters=nothing: The argument p in diffusion_function.
  • kwargs...: Any other keyword arguments are passed to the LinearProblem (from LinearSolve.jl) that represents the problem.

Fields

The struct has extra fields in addition to the arguments above:

  • A: This is a sparse matrix A so that AT = b.
  • b: The b above.
  • problem: The LinearProblem that represents the problem. This is the problem that is solved when you call solve on the struct.
FiniteVolumeMethod.ParametrisedFunctionType
ParametrisedFunction{F<:Function,P} <: Function

This is a struct that wraps a function f and some parameters p into a single object.

Fields

  • fnc::F

The function that is wrapped.

  • parameters::P

The parameters that are wrapped.

FiniteVolumeMethod.PoissonsEquationType
PoissonsEquation

A struct for defining a problem representing a (generalised) Poisson's equation:

\[\div[D(\vb x)\grad u] = f(\vb x)\]

inside a domain $\Omega$. See also LaplacesEquation, a special case of this problem with $f(\vb x) = 0$.

You can solve this problem using solve.

Constructor

PoissonsEquation(mesh::FVMGeometry,
    BCs::BoundaryConditions,
    ICs::InternalConditions=InternalConditions();
    diffusion_function=(x,y,p)->1.0,
    diffusion_parameters=nothing,
    source_function, 
    source_parameters=nothing,
    kwargs...)

Arguments

  • mesh::FVMGeometry: The FVMGeometry.
  • BCs::BoundaryConditions: The BoundaryConditions. For these boundary conditions, all functions should still be of the form (x, y, t, u, p) -> Number, but the t and u arguments should be unused as they will be replaced with nothing.
  • ICs::InternalConditions=InternalConditions(): The InternalConditions. For these internal conditions, all functions should still be of the form (x, y, t, u, p) -> Number, but the t and u arguments should be unused as they will be replaced with nothing.

Keyword Arguments

  • diffusion_function=(x,y,p)->1.0: The diffusion function. Should be of the form (x, y, p) -> Number, where p = diffusion_parameters below.
  • diffusion_parameters=nothing: The argument p in diffusion_function.
  • source_function: The source function. Should be of the form (x, y, p) -> Number, where p = source_parameters below.
  • source_parameters=nothing: The argument p in source_function.
  • kwargs...: Any other keyword arguments are passed to the LinearProblem (from LinearSolve.jl) that represents the problem.

Fields

The struct has extra fields in addition to the arguments above:

  • A: This is a sparse matrix A so that Au = b.
  • b: The b above.
  • problem: The LinearProblem that represents the problem. This is the problem that is solved when you call solve on the struct.
FiniteVolumeMethod.TrianglePropertiesType
TriangleProperties(shape_function_coefficients, cv_edge_midpoints, cv_edge_normals, cv_edge_lengths)

This is a struct for holding the properties of a control volume's intersection with a triangle.

Fields

  • shape_function_coefficients::NTuple{9,Float64}: The shape function coefficients for the triangle.
  • cv_edge_midpoints::NTuple{3,NTuple{2,Float64}}: The midpoints of the control volume edges. If the triangle is (i, j, k), then the edges are given in the order (i, j), (j, k), and (k, i), where 'edge' refers to the edge joining e.g. the midpoint of the edge (i, j) to the centroid of the triangle.
  • cv_edge_normals::NTuple{3,NTuple{2,Float64}}: The normal vectors to the control volume edges, in the same order as in cv_edge_midpoints.
  • cv_edge_lengths::NTuple{3,Float64}: The lengths of the control volume edges, in the same order as in cv_edge_midpoints.
Notes

The shape function coefficients are defined so that, if s are the coefficients and the triangle is T = (i, j, k), with function values u[i], u[j], and u[k] at the vertices i, j, and k, respectively, then αxₙ + βyₙ + γₙ = u[n] for n = i, j, k, where xₙ and yₙ are the x- and y-coordinates of the nth vertex, respectively, α = s₁u₁ + s₂u₂ + s₃u₃, β = s₄u₁ + s₅u₂ + s₆u₃, and γ = s₇u₁ + s₈u₂ + s₉u₃.

CommonSolve.solveMethod
solve(prob::AbstractFVMTemplate, args...; kwargs...)

Solve the problem prob using the standard solve interface from DifferentialEquations.jl. For steady state problems, the interface is from LinearSolve.jl.

CommonSolve.solveMethod
solve(prob::SteadyFVMProblem, alg; 
    specialization=SciMLBase.AutoSpecialize, 
    jac_prototype=jacobian_sparsity(prob),
    parallel::Val{<:Bool}=Val(true),
    kwargs...)

Solves the given FVMProblem or FVMSystem prob with the algorithm alg.

Arguments

  • prob: The problem to be solved.
  • alg: The algorithm to be used to solve the problem. This can be any of the algorithms in NonlinearSolve.jl.

Keyword Arguments

  • specialization=SciMLBase.AutoSpecialize: The type of specialization to be used. See https://docs.sciml.ai/DiffEqDocs/stable/features/low_dep/#Controlling-Function-Specialization-and-Precompilation.
  • jac_prototype=jacobian_sparsity(prob): The prototype for the Jacobian matrix, constructed by default from jacobian_sparsity.
  • parallel::Val{<:Bool}=Val(true): Whether to use multithreading. Use Val(false) to disable multithreading.
  • kwargs...: Any other keyword arguments to be passed to the solver.

Outputs

The returned value sol depends on whether the underlying problem is a FVMProblem or an FVMSystem, but in each case it is an ODESolution type that can be accessed like the solutions in DifferentialEquations.jl:

In this case, sol is such that the ith component of sol refers to the ith node of the underlying mesh.

In this case, the (j, i)th component of sol refers to the ith node of the underlying mesh for the jth component of the system.

CommonSolve.solveMethod
solve(prob::Union{FVMProblem,FVMSystem}, alg; 
    specialization=SciMLBase.AutoSpecialize, 
    jac_prototype=jacobian_sparsity(prob), 
    parallel::Val{<:Bool}=Val(true),
    kwargs...)

Solves the given FVMProblem or FVMSystem prob with the algorithm alg.

Missing vertices

When the underlying triangulation, tri, has points in get_points(tri) that are not vertices of the triangulation itself, the associated values of the solution at these points will not be updated by the solver, and will remain at their initial values.

Arguments

  • prob: The problem to be solved.
  • alg: The algorithm to be used to solve the problem. This can be any of the algorithms in DifferentialEquations.jl.

Keyword Arguments

  • specialization=SciMLBase.AutoSpecialize: The type of specialization to be used. See https://docs.sciml.ai/DiffEqDocs/stable/features/low_dep/#Controlling-Function-Specialization-and-Precompilation.
  • jac_prototype=jacobian_sparsity(prob): The prototype for the Jacobian matrix, constructed by default from jacobian_sparsity.
  • parallel::Val{<:Bool}=Val(true): Whether to use multithreading. Use Val(false) to disable multithreading.
  • kwargs...: Any other keyword arguments to be passed to the solver.

Outputs

The returned value sol depends on the type of the problem.

In this case, sol::ODESolution is such that the ith component of sol refers to the ith node of the underlying mesh.

In this case, the (j, i)th component of sol::ODESolution refers to the ith node of the underlying mesh for the jth component of the system.

FiniteVolumeMethod.apply_dirichlet_conditions!Method
apply_dirichlet_conditions!(initial_condition, mesh, conditions)

Applies the Dirichlet conditions specified in conditions to the initial_condition. The boundary conditions are assumed to take the form a(x, y, t, u, p) -> Number, but t and u are passed as nothing. Note that this assumes that the associated system (A, b) is such that A[i, :] is all zero, and b[i] is zero, where i is a node with a Dirichlet condition.

FiniteVolumeMethod.apply_dudt_conditions!Method
apply_dudt_conditions!(b, mesh, conditions)

Applies the Dudt conditions specified in conditions to the b vector. The boundary conditions are assumed to take the form a(x, y, t, u, p) -> Number, but t and u are passed as nothing. Note that this assumes that the associated system (A, b) is such that A[i, :] is all zero, so that replacing b[i] with the boundary condition will set duᵢ/dt = b[i].

FiniteVolumeMethod.apply_steady_dirichlet_conditions!Method
apply_steady_dirichlet_conditions!(A, b, mesh, conditions)

Applies the Dirichlet conditions specified in conditions to the initial_condition. The boundary conditions are assumed to take the form a(x, y, t, u, p) -> Number, but t and u are passed as nothing. Note that this assumes that the associated system (A, b) is such that A[i, :] is all zero, and b[i] is zero, where i is a node with a Dirichlet condition.

For a steady problem Au = b, applies the Dirichlet boundary conditions specified by conditions so that A[i, i] = 1 and b[i] is the condition, where i is a boundary node. Note that this assumes that all of A[i, :] is zero before setting A[i, i] = 1.

FiniteVolumeMethod.boundary_edge_contributions!Method
boundary_edge_contributions!(A, b, mesh, conditions, diffusion_function, diffusion_parameters)

Add the contributions from each boundary edge to the matrix A, based on the equation

\[\dv{u_i}{t} = \frac{1}{V_i}\sum_{\sigma \in \mathcal E_i} D(\vb x_\sigma)\left[\left(s_{k, 11}n_\sigma^x + s_{k, 21}n_\sigma^y\right)u_{k1} + \left(s_{k, 12}n_\sigma^x + s_{k, 22}n_\sigma^y\right)u_{k2} + \left(s_{k, 13}n_\sigma^x + s_{k, 23}n_\sigma^y\right)u_{k3}\right]L_\sigma + S_i, \]

as explained in the docs. Will not update any rows corresponding to Dirichlet or Dudt nodes.

FiniteVolumeMethod.compute_fluxMethod
compute_flux(prob::AbstractFVMProblem, i, j, u, t)

Given an edge with indices (i, j), a solution vector u, a current time t, and a problem prob, computes the flux ∇⋅q(x, y, t, α, β, γ, p) ⋅ n, where n is the normal vector to the edge, q is the flux function from prob, (x, y) is the midpoint of the edge, (α, β, γ) are the shape function coefficients, and p are the flux parameters from prob. If prob is an FVMSystem, the returned value is a Tuple for each individual flux. The normal vector n is a clockwise rotation of the edge, meaning pointing right of the edge (i, j).

FiniteVolumeMethod.construct_flux_functionMethod
construct_flux_function(q, D, Dp)

If isnothing(q), then this returns the flux function based on the diffusion function D and diffusion parameters Dp, so that the new function is

(x, y, t, α, β, γ, _) -> -D(x, y, t, α*x + β*y + γ, Dp)[α, β]

Otherwise, just returns q again.

FiniteVolumeMethod.create_rhs_bMethod
create_rhs_b(mesh, conditions, source_function, source_parameters)

Create the vector b defined by

b = [source_function(x, y, source_parameters) for (x, y) in DelaunayTriangulation.each_point(mesh.triangulation)],

and b[i] = 0 whenever i is a Dirichlet node.

FiniteVolumeMethod.eval_condition_fncMethod
eval_condition_fnc(conds, fidx, x, y, t, u)

Evaluate the function that corresponds to the condition at fidx at the point (x, y) at time t with solution u.

FiniteVolumeMethod.fix_missing_vertices!Method
fix_missing_vertices!(A, b, mesh)

Given a system (A, b) and a mesh, sets A[i, i] = 1 and b[i] = 0 for any vertices i that are a point in mesh, but not an actual vertex in the mesh.

FiniteVolumeMethod.fvm_eqs!Method
fvm_eqs!(du, u, p, t)

Computes the finite volume method equations for the current time t and solution u. This function is public API.

The argument p depends on whether the problem is being solved in parallel or not. If it is solved serially, than the fields are:

  • p.prob: The prob <: AbstractFVMProblem.
  • p.parallel: Val(false).

If the problem is solved in parallel, then the fields are:

  • p.prob: The prob <: AbstractFVMProblem.
  • p.parallel: Val(true).
  • p.duplicated_du: A Matrix of the same size as du that is used to store the contributions to du from each thread.
  • p.solid_triangles: A Vector of the solid triangles in the triangulation.
  • p.solid_vertices: A Vector of the solid vertices in the triangulation.
  • p.chunked_solid_triangles: A Vector of tuples of the form (range, chunk_idx) where range is a range of indices into p.solid_triangles and chunk_idx is the index of the chunk.
  • p.boundary_edges: A Vector of the boundary edges in the triangulation.
  • p.chunked_boundary_edges: A Vector of tuples of the form (range, chunk_idx) where range is a range of indices into p.boundary_edges and chunk_idx is the index of the chunk.

These fields are public API, although note that they are not intended to be modified by the user, and we may freely add in new fields over new versions.

FiniteVolumeMethod.get_boundary_cv_componentsMethod
get_boundary_cv_components(tri::Triangulation, i, j)

Get the quantities for both control volume edges lying a boundary edge (i, j).

Outputs

  • nx: The x-component of the edge's normal vector.
  • ny: The y-component of the edge's normal vector.
  • mᵢx: The x-coordinate of the midpoint of the ith vertex and the edge's midpoint.
  • mᵢy: The y-coordinate of the midpoint of the ith vertex and the edge's midpoint.
  • mⱼx: The x-coordinate of the midpoint of the jth vertex and the edge's midpoint.
  • mⱼy: The y-coordinate of the midpoint of the jth vertex and the edge's midpoint.
  • ℓᵢ: Half the length of the boundary edge, which is the length of the control volume edge.
  • T: The triangle containing the boundary edge.
  • props: The TriangleProperties for T.
FiniteVolumeMethod.get_cv_componentsMethod
get_cv_components(props, edge_index)

Get the quantities for a control volume edge interior to the associated triangulation, relative to the edge_indexth edge of the triangle corresponding to props.

Outputs

  • x: The x-coordinate of the edge's midpoint.
  • y: The y-coordinate of the edge's midpoint.
  • nx: The x-component of the edge's normal vector.
  • ny: The y-component of the edge's normal vector.
  • : The length of the edge.
FiniteVolumeMethod.get_dirichlet_callbackMethod
get_dirichlet_callback(prob[, f=update_dirichlet_nodes!]; kwargs...)

Get the callback for updating Dirichlet nodes. The kwargs... argument is ignored, except to detect if a user has already provided a callback, in which case the callback gets merged into a CallbackSet with the Dirichlet callback. If the problem prob has no Dirichlet nodes, the returned callback does nothing and is never called.

You can provide f to change the function that updates the Dirichlet nodes.

FiniteVolumeMethod.neumann_boundary_edge_contributions!Method
neumann_boundary_edge_contributions!(b, mesh, conditions, diffusion_function, diffusion_parameters)

Add the contributions from each Neumann boundary edge to the vector b, based on the equation

\[\dv{u_i}{t} = \frac{1}{V_i}\sum_{\sigma \in \mathcal E_i} D(\vb x_\sigma)\left[\grad u(\vb x_\sigma) \vdot \vu n\right]L_\sigma + S_i,\]

as explained in the docs. Will not update any rows corresponding to Dirichlet or Dudt nodes. This function will pass nothing in place of the arguments u and t in the boundary condition functions.

FiniteVolumeMethod.neumann_boundary_edge_contributions!Method
neumann_boundary_edge_contributions!(F, mesh, conditions, diffusion_function, diffusion_parameters, u, t)

Add the contributions from each Neumann boundary edge to the vector F, based on the equation

\[\dv{u_i}{t} = \frac{1}{V_i}\sum_{\sigma \in \mathcal E_i} D(\vb x_\sigma)\left[\grad u(\vb x_\sigma) \vdot \vu n\right]L_\sigma + S_i,\]

as explained in the docs. Will not update any rows corresponding to Dirichlet or Dudt nodes.

FiniteVolumeMethod.non_neumann_boundary_edge_contributions!Method
non_neumann_boundary_edge_contributions!(A, mesh, conditions, diffusion_function, diffusion_parameters)

Add the contributions from each non-Neumann boundary edge to the matrix A, based on the equation

\[\dv{u_i}{t} = \frac{1}{V_i}\sum_{\sigma \in \mathcal E_i} D(\vb x_\sigma)\left[\left(s_{k, 11}n_\sigma^x + s_{k, 21}n_\sigma^y\right)u_{k1} + \left(s_{k, 12}n_\sigma^x + s_{k, 22}n_\sigma^y\right)u_{k2} + \left(s_{k, 13}n_\sigma^x + s_{k, 23}n_\sigma^y\right)u_{k3}\right]L_\sigma + S_i, \]

as explained in the docs. Will not update any rows corresponding to Dirichlet or Dudt nodes.

FiniteVolumeMethod.pl_interpolateMethod
pl_interpolate(prob, T, u, x, y)

Given a prob <: AbstractFVMProblem, a triangle T containing a point (x, y), and a set of function values u at the corresponding nodes of prob, interpolates the solution at the point (x, y) using piecewise linear interpolation.

FiniteVolumeMethod.triangle_contributions!Method
triangle_contributions!(A, mesh, conditions, diffusion_function, diffusion_parameters)

Add the contributions from each triangle to the matrix A, based on the equation

\[\dv{u_i}{t} = \frac{1}{V_i}\sum_{\sigma \in \mathcal E_i} D(\vb x_\sigma)\left[\left(s_{k, 11}n_\sigma^x + s_{k, 21}n_\sigma^y\right)u_{k1} + \left(s_{k, 12}n_\sigma^x + s_{k, 22}n_\sigma^y\right)u_{k2} + \left(s_{k, 13}n_\sigma^x + s_{k, 23}n_\sigma^y\right)u_{k3}\right]L_\sigma + S_i, \]

as explained in the docs. Will not update any rows corresponding to Dirichlet or Dudt nodes.

FiniteVolumeMethod.two_point_interpolantMethod
two_point_interpolant(mesh, u, i, j, mx, my)

Given a mesh <: FVMGeometry, a set of function values u at the nodes of mesh, and a point (mx, my) on the line segment between the nodes i and j, interpolates the solution at the point (mx, my) using two-point interpolation.