Physics & special functions

Physics

VoronoiFVM.PhysicsType
struct Physics

Physics data record with the following fields:

  • flux::Function: Flux between neighboring control volumes: flux(f,u,edge) or flux(f,u,edge,data) should return in f[i] the flux of species i along the edge joining circumcenters of neighboring control volumes. u is a 2D array such that for species i, u[i,1] and u[i,2] contain the unknown values at the corresponding ends of the edge.
  • storage::Function: Storage term (term under time derivative): storage(f,u,node) or storage(f,u,node,data)

    It should return in f[i] the storage term for the i-th equation. u[i] contains the value of the i-th unknown.

  • reaction::Function: Reaction term: reaction(f,u,node) or reaction(f,u,node,data)

    It should return in f[i] the reaction term for the i-th equation. u[i] contains the value of the i-th unknown.

  • edgereaction::Function: Edge reaction term: edgereaction(f,u,edge) or edgereaction(f,u,edge,data)

    It should return in f[i] the reaction term for the i-th equation. u[i] contains the value of the i-th unknown. u is a 2D array such that for species i, u[i,1] and u[i,2] contain the unknown values at the corresponding ends of the edge.

  • source::Function: Source term: source(f,node) or source(f,node,data).

    It should return the in f[i] the value of the source term for the i-th equation.

  • bflux::Function: Flux between neighboring control volumes on the boundary. Called on edges fully adjacent to the boundary: bflux(f,u,bedge) or `bflux(f,u,bedge,data)
  • breaction::Function: Boundary reaction term: breaction(f,u,node) or breaction(f,u,node,data) Similar to reaction, but restricted to the inner or outer boundaries.
  • bsource::Function: Boundary source term: bsource(f,node) or bsource(f,node,data).

    It should return in f[i] the value of the source term for the i-th equation.

  • bstorage::Function: Boundary storage term: bstorage(f,u,node) or bstorage(f,u,node,data) Similar to storage, but restricted to the inner or outer boundaries.
  • boutflow::Function: Outflow boundary term boutflow(f,u,edge) or boutflow(f,u,edge,data) This function is called for edges (including interior ones) which have at least one ode on one of the outflow boundaries. Within this function, outflownode and isoutflownode can be used to identify that node. There is some ambiguity in the case that both nodes are outflow nodes, in that case it is assumed that the contribution is zero.
  • outflowboundaries::Vector{Int64}: List (Vector) of boundary regions which carry outflow bondary conditions. Influences when boutflow is called.
  • generic_operator::Function: Generic operator generic_operator(f,u,sys). This operator acts on the full solution u of a system. Sparsity is detected automatically unless generic_operator_sparsity is given.
  • generic_operator_sparsity::Function: Function defining the sparsity structure of the generic operator. This should return the sparsity pattern of the generic_operator.
  • data::Any: User data (parameters). This allows to pass various parameters to the callback functions.
  • num_species::Int8: Number of species including boundary species.
source

Edge and node data

VoronoiFVM.NodeType
mutable struct Node{Tc, Tp, Ti} <: VoronoiFVM.AbstractNode{Tc, Tp, Ti}

Structure holding local node information.

  • index::Any: Index in grid
  • region::Any: Inner region number
  • nspec::Any: Number of species defined in node
  • icell::Any: Number of discretization cell the node is invoked from
  • coord::Matrix: Grid coordinates
  • cellnodes::Matrix: Grid cell nodes
  • cellregions::Vector: Grid cell regions
  • time::Float64: System time
  • embedparam::Float64: Current value of embedding parameter
  • params::Vector: parameters
  • fac::Float64: Form factor
  • _idx::Any: Local loop index
source
VoronoiFVM.BNodeType
mutable struct BNode{Tv, Tc, Tp, Ti} <: VoronoiFVM.AbstractNode{Tc, Tp, Ti}

Structure holding local boundary node information.

  • index::Any: Index in grid
  • ibface::Any: BFace number it is called from
  • ibnode::Any: local node number
  • region::Any: Boundary region number
  • cellregions::Vector

  • nspec::Any: Number of species defined in node

  • coord::Matrix: Grid coordinates
  • bfacenodes::Matrix

  • bfaceregions::Vector

  • allcellregions::Vector

  • bfacecells::Adjacency

  • Dirichlet::Any

  • time::Float64: System time

  • embedparam::Float64: Current value of embedding parameter
  • params::Vector

  • dirichlet_value::Vector

  • fac::Float64

source
VoronoiFVM.EdgeType
mutable struct Edge{Tc, Tp, Ti} <: VoronoiFVM.AbstractEdge{Tc, Tp, Ti}

Structure holding local edge information.

  • index::Any: Index in grid
  • node::Vector: Index
  • region::Any: Inner region number corresponding to edge
  • nspec::Any: Number of species defined in edge
  • icell::Any: Number of discretization cell the edge is invoked from
  • coord::Matrix: Grid coordinates
  • cellx::Matrix

  • edgenodes::Matrix

  • cellregions::Vector

  • has_celledges::Bool

  • time::Float64: System time

  • embedparam::Float64: Current value of embedding parameter
  • params::Vector

  • fac::Float64: Form factor

  • _idx::Any: Local loop index
  • outflownoderegions::Union{Nothing, SparseArrays.SparseMatrixCSC{Bool, Int64}}

  • outflownode::Int64: Outflow node

source
VoronoiFVM.BEdgeType
mutable struct BEdge{Tc, Tp, Ti} <: VoronoiFVM.AbstractEdge{Tc, Tp, Ti}

Structure holding local edge information.

  • index::Any: Index in grid
  • node::Vector: Index
  • region::Any: Inner region number corresponding to edge
  • nspec::Any: Number of species defined in edge
  • icell::Any: Number of discretization cell the edge is invoked from
  • coord::Matrix: Grid coordinates
  • bedgenodes::Matrix

  • bfaceedges::Matrix

  • bfaceregions::Vector

  • time::Float64: System time

  • embedparam::Float64: Current value of embedding parameter
  • params::Vector

  • fac::Float64

source

Special functions

VoronoiFVM.fbernoulliFunction
fbernoulli(x)

Bernoulli function $B(x)=\frac{x}{e^x-1}$ for exponentially fitted upwinding.

The name fbernoulli has been chosen to avoid confusion with Bernoulli from JuliaStats/Distributions.jl

Returns a real number containing the result.

source
VoronoiFVM.fbernoulli_pmFunction
fbernoulli_pm(x)

Bernoulli function $B(x)=\frac{x}{e^x-1}$ for exponentially fitted upwind, joint evaluation for positive and negative argument

Usually, we need $B(x), B(-x)$ together, and it is cheaper to calculate them together.

Returns two real numbers containing the result for argument x and argument -x.

The error in comparison with the evaluation of the original expression with BigFloat is less than 1.0e-15

source
VoronoiFVM.inplace_linsolve!Method
inplace_linsolve!(A, b)

Non-allocating, non-pivoting inplace solution of square linear system of equations A*x=b using Doolittle's method.

After solution, A will contain the LU factorization, and b the result.

A pivoting version is available with Julia v1.9.

source
VoronoiFVM.inplace_linsolve!Method
inplace_linsolve!(A, b, ipiv)

Non-allocating, pivoting, inplace solution of square linear system of equations A*x=b using LU factorization from RecursiveFactorizations.jl.

After solution, A will contain the LU factorization, and b the result. ipiv must be an Int64 vector of the same length as b.

Info

Needs Julia v1.9 or later

source