Physics & special functions
Physics callbacks
VoronoiFVM.AbstractPhysics
— Typeabstract type AbstractPhysics
Abstract type for physics.
VoronoiFVM.Physics
— Typestruct Physics
Physics data record with the following fields:
flux::Function
: Flux between neighboring control volumes:flux(f,u,edge,data)
should return inf[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]
andu[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,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,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,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]
andu[i,2]
contain the unknown values at the corresponding ends of the edge.
source::Function
: Source term: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,data)
breaction::Function
: Boundary reaction term:breaction(f,u,node,data)
Similar to reaction, but restricted to the inner or outer boundaries.
bsource::Function
: Boundary source term: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,data)
Similar to storage, but restricted to the inner or outer boundaries.
boutflow::Function
: Outflow boundary termboutflow(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
andisoutflownode
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 whenboutflow
is called.
generic_operator::Function
: Generic operatorgeneric_operator(f,u,sys)
. This operator acts on the full solutionu
of a system. Sparsity is detected automatically unlessgeneric_operator_sparsity
is given.
generic_operator_sparsity::Function
: Function defining the sparsity structure of the generic operator. This should return the sparsity pattern of thegeneric_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. Meaningless & deprecated.
VoronoiFVM.Physics
— MethodPhysics(;num_species=0,
data=nothing,
flux,
reaction,
edgereaction,
storage,
source,
breaction,
bstorage,
boutflow,
outflowboundaries,
generic,
generic_sparsity
)
Constructor for physics data. For the meaning of the optional keyword arguments, see VoronoiFVM.System(grid::ExtendableGrid; kwargs...)
.
VoronoiFVM.generic_operator_sparsity!
— Functiongeneric_operator_sparsity!(system, sparsematrix)
Set generic operator sparsity, in the case where a generic operator has been defined in physics.
Base.show
— Methodshow(io, physics)
Show physics object
VoronoiFVM.AbstractData
— Typeabstract type AbstractData{Tv}
Abstract type for user data.
It is possible but not necessary to make user data a subtype of AbstractData and get a prettyprinting show method.
Base.show
— Methodshow(
io::IO,
_::MIME{Symbol("text/plain")},
this::VoronoiFVM.AbstractData
)
Pretty print AbstractData
Base.copy!
— Methodcopy!(to::AbstractData, from::AbstractData)
Copy AbstractData
.
Handling boundary conditions
Boundary conditions are handled in the bcondition
callback passed to the system constructor. For being called in this callback, the following functions are available
VoronoiFVM.boundary_dirichlet!
— Method boundary_dirichlet!(y,u,bnode,ispec,ireg,val)
Set Dirichlet boundary condition for species ispec at boundary ibc.
VoronoiFVM.boundary_dirichlet!
— Method boundary_dirichlet!(y,u,bnode; kwargs...)
Keyword argument version:
species
: species number. Default: 1region
: boundary region number. By default, all boundary regions.value
: value
VoronoiFVM.boundary_neumann!
— Method boundary_neumann!(y,u,bnode,ispec,ireg,val)
Set Neumann boundary condition for species ispec at boundary ibc.
VoronoiFVM.boundary_neumann!
— Method boundary_neumann!(y,u,bnode; kwargs...)
Keyword argument version:
species
: species number. Default: 1region
: boundary region number. By default, all boundary regions.value
: value
VoronoiFVM.boundary_robin!
— Method boundary_robin!(y,u,bnode,ispec,ireg,fac,val)
Set Robin boundary condition for species ispec at boundary ibc.
VoronoiFVM.boundary_robin!
— Method boundary_robin!(y,u,bnode, args...; kwargs...)
Keyword argument version:
species
: species number. Default: 1region
: boundary region number. By default, all boundary regions.factor
: factorvalue
: value
VoronoiFVM.ramp
— Function ramp(t; kwargs...)
Ramp function for specifying time dependent boundary conditions
Keyword arguments:
dt
: Tuple: start and end time of ramp. Default:(0,0.1)
du
: Tuple: values at start and end time. Default:(0,0)
Outflow boundary conditions
These are characterized by the boutflow
physics callback and and the outflowboundaries
keyword argument in the system resp. physics constructor. See also the corresponding notebook
VoronoiFVM.hasoutflownode
— Functionhasoutflownode(edge)
Check if one node of the edge is situated on a boundary region listed in outflowboundaries
, see [struct Physics
].
VoronoiFVM.isoutflownode
— Functionisoutflownode(edge,inode)
Check if inode (1 or 2) is an outflow node.
isoutflownode(edge,inode,irefgion)
Check if inode (1 or 2) is an outflow node on boundary region iregion
.
VoronoiFVM.outflownode
— Functionoutflownode(edge)
Return outflow node of edge (1 or 2).
VoronoiFVM.calc_divergences
— Functioncalc_divergences(sys, evelo, bfvelo)
Calculate for each Voronoi cell associated with a node of sys.grid the divergence of the velocity field used to obtain evelo
and bfvelo
via edgevelocities
and bfacevelocities
by means of summing all evelos
and bfvelos
per Voronoi cell.
Coupling to flow
VoronoiFVM.edgevelocities
— Functionedgevelocities(grid, velofunc; kwargs...)
Project velocity onto grid edges. That is, we compute the path integrals of the given velofunc
along the Voronoi cell edges as provided by integrate
.
edgevelocities(grid, vel; kwargs...)
Compute VoronoiFVM.edgevelocities
for a finite element flow field computed by ExtendableFEM
.
VoronoiFVM.bfacevelocities
— Functionbfacevelocities(grid, velofunc; kwargs...)
Similar to edgevelocities
, but for boundary faces.
bfacevelocities(grid, vel; kwargs...)
Compute VoronoiFVM.bfacevelocities
for a finite element flow field computed by ExtendableFEM
.
VoronoiFVM.bfacenodefactors
— Functionbfacenodefactors(sys)
Node form factors per boundary face
Edge and node data
VoronoiFVM.Node
— Typemutable 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 (deprecated)
fac::Float64
: Form factor
_idx::Any
: Local loop index
VoronoiFVM.BNode
— Typemutable 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
: Parameters (deprecated)
dirichlet_value::Vector
fac::Float64
VoronoiFVM.Edge
— Typemutable 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
: Parameters (deprecated)
fac::Float64
: Form factor
_idx::Any
: Local loop index
outflownoderegions::Union{Nothing, SparseArrays.SparseMatrixCSC{Bool, Int64}}
outflownode::Int64
: Outflow node
VoronoiFVM.BEdge
— Typemutable 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
: Parameters (deprecated)
fac::Float64
VoronoiFVM.time
— Functiontime(edge_or_node)
Return actual simulation time stored in node or edge
VoronoiFVM.region
— Functionregion(edgeornode)
Return region number node or edge is belonging to
VoronoiFVM.parameters
— Functionparameters(edge_or_node)
Return abstract vector of parameters passed via vector of unknonws. This allows differentiation with respect to these parameters.
VoronoiFVM.embedparam
— Functionembedparam(edge_or_node)
Return embedding parameter stored in node or edge
VoronoiFVM.project
— Functionproject(edge, vector)
Project d-vector onto d-dimensional vector, i.e. calculate the dot product of vector
with the difference of the edge end coordinates.
VoronoiFVM.edgelength
— Functionedgelength(edge)
Return length of edge
VoronoiFVM.meas
— Functionmeas(edge::VoronoiFVM.AbstractEdge) -> Any
Calculate the length of an edge.
Special functions
VoronoiFVM.fbernoulli
— Functionfbernoulli(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.
While x/expm1(x)
appears to be sufficiently accurate for all x!=0
, it's derivative calculated via ForwardDiff.jl
is not, so we use the polynomial approximation in the interval (-bernoulli_small_threshold, bernoulli_small_threshold)
.
Also, see the discussion in #117.
VoronoiFVM.fbernoulli_pm
— Functionfbernoulli_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
.
For the general approach to the implementation, see the discussion in fbernoulli
.
VoronoiFVM.inplace_linsolve!
— Methodinplace_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.
VoronoiFVM.inplace_linsolve!
— Methodinplace_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
.