API Documentation
TwoPointFluxFVM.SubGrid
— Type.struct SubGrid{Tc} <: AbstractGrid
Subgrid of parent grid (mainly for visualization purposes). Intended to hold support of species which are not defined everywhere.
TwoPointFluxFVM.SubGrid
— Method.function SubGrid(parent::Grid,
subregions::AbstractArray;
transform::Function=copytransform!,
boundary=false)
Create subgrid of list of regions.
TwoPointFluxFVM.add_boundary_species
— Method.function add_boundary_species(this::System, ispec::Integer, regions::AbstractArray)
Add species to a list of boundary regions. Species numbers for bulk and boundary species have to be distinct.
TwoPointFluxFVM.add_species
— Method.function add_species(this::System,ispec::Integer, regions::AbstractArray)
Add species to a list of bulk regions. Species numbers for bulk and boundary species have to be distinct.
TwoPointFluxFVM.cellmask!
— Method.function cellmask!(grid::Grid,
maskmin::AbstractArray, # lower left corner
maskmax::AbstractArray, # upper right corner
ireg::Integer; # new region number for elements under mask
eps=1.0e-10) # tolerance.
Edit region numbers of grid cells via rectangular mask.
TwoPointFluxFVM.dof
— Method.function dof(a::SysArray,ispec, inode)
Get number of degree of freedom. Return 0 if species is not defined in node.
TwoPointFluxFVM.fbernoulli
— Method.function fbernoulli(x::Real)
Bernoulli function implementation for exponentially fitted finite volumes.
The name fbernoulli has been chosen to avoid confusion with Bernoulli from JuliaStats/Distributions.jl
Returns a real number containing the result.
TwoPointFluxFVM.fbernoulli_pm
— Method.function fbernoulli_pm(x::Real)
Bernoulli function implementation for exponentially fitted finite volumes, joint evaluation for positive and negative argument
Usually, we need B(x), B(-x) togehter, and it is cheaper to calculate them together.
Returns two real numbers containing the result for argument x
and argument -x
.
TwoPointFluxFVM.getdof
— Method.function getdof(a::SysArray,i::Integer)
Return value for degree of freedom.
TwoPointFluxFVM.integrate
— Method.function integrate(this::System,F::Function,U)
Integrate solution vector over domain. Returns an Array{Int64,1}
containing the integral for each species.
TwoPointFluxFVM.num_bfaces
— Method.num_bfaces(grid::Grid)
Number of boundary faces in grid.
TwoPointFluxFVM.num_cells
— Method.num_cells(grid)
Number of cells in grid
TwoPointFluxFVM.num_nodes
— Method.num_nodes(grid)
Number of nodes in grid
TwoPointFluxFVM.num_nodes
— Method.num_nodes(a::SysArray)
Number of nodes (size of second dimension) of solution array.
TwoPointFluxFVM.setdof!
— Method.function setdof!(a::SysArray,v,i::Integer)
Set value for degree of freedom.
TwoPointFluxFVM.solve
— Method.function solve(
this::System, # Finite volume system
oldsol::Array{Tv,1}; # old time step solution resp. initial value
control=NewtonControl(), # Solver control information (optional)
tstep::Tv=Inf # Time step size. Inf means stationary solution. (optional)
)
Solution method for instance of System.
Perform solution of stationary system (if tstep==Inf
) or implicit Euler time step system.
TwoPointFluxFVM.unknowns
— Method.function unknowns(system)
Create a solution vector for system.
TwoPointFluxFVM.value
— Function.value(x)
Extract value from dual number. Use to debug physics callbacks. Re-exported from ForwardDiff.jl
TwoPointFluxFVM.AbstractGrid
— Type.abstract type AbstractGrid
Abstract type for grid like datastructures.
TwoPointFluxFVM.Edge
— Type.mutable struct Edge
Structure holding local edge information.
Fields:
index::Int32
region::Int32
nodeK::Int32
nodeL::Int32
coordK::Array{Float64,1}
coordL::Array{Float64,1}
TwoPointFluxFVM.Grid
— Type.mutable struct Grid
Structure holding grid data.
TwoPointFluxFVM.Grid
— Method.Grid(X::Array{Tc,1})
Constructor for 1D grid.
Construct 1D grid from an array of node cordinates. It creates two boundary regions with index 1 at the left end and index 2 at the right end.
Primal grid holding unknowns: marked by o
, dual grid marking control volumes: marked by |
.
o-----o-----o-----o-----o-----o-----o-----o-----o
|--|-----|-----|-----|-----|-----|-----|-----|--|
TwoPointFluxFVM.Grid
— Method.Grid(X::Array{Tc,1},X::Array{Tc,1})
Constructor for 2D grid from coordinate arrays. Boundary region numbers count counterclockwise:
location | number |
---|---|
south | 1 |
east | 2 |
north | 3 |
west | 4 |
TwoPointFluxFVM.NewtonControl
— Type.Control parameter structure for Newton method.
Fields:
tol_absolute::Float64 # Tolerance (in terms of norm of Newton update)
tol_relative::Float64 # Tolerance (relative to the first residual)
damp_initial::Float64 # Initial damping parameter
damp_growth::Float64 # Damping parameter growth factor
max_iterations::Int32 # Maximum number of iterations
max_lureuse::Int32 # Maximum number of reuses of lu factorization
tol_linear::Float64 # Tolerance of iterative linear solver
verbose::Bool # Verbosity
TwoPointFluxFVM.NewtonControl
— Method.NewtonControl()
Default constructor
TwoPointFluxFVM.Node
— Type.mutable struct Node
Structure holding local node information. Fields:
index::Int32
region::Int32
coord::Array{Real,1}
struct SubgridSysArrayView{Tv} <: AbstractArray{Tv,2}
Struct holding information for solution array view on subgrid
TwoPointFluxFVM.SysArray
— Type.struct SysArray{Tv} <: AbstractArray{Tv,2}
node_dof::SparseMatrixCSC{Tv,Int16}
end
Struct holding solution information for system. Solution is stored in a sparse matrix structure.
This class plays well with the abstract array interface
TwoPointFluxFVM.System
— Type.mutable struct System{Tv}
Main structure holding data for system solution.
TwoPointFluxFVM.System
— Method.function System(grid::Grid, physics::Any, maxspec::Integer)
Constructor for System. physics
provides some user data, maxspec
is the maximum number of species.
Base.copy
— Method.copy(this::SysArray)
Create a copy of solution array
Base.eltype
— Method.eltype(grid)
Return element type of grid coordinates.
Base.getindex
— Method.getindex(aview::SubgridSysArrayView,ispec::Integer,inode::Integer)
Accessor method for subgrid array view.
Base.getindex
— Method. getindex!(a::SysArray, ispec, inode)
Accessor for solution array.
Base.setindex!
— Method.setindex!(aview::SubgridSysArrayView,v,ispec::Integer,inode::Integer)
Accessor method for subgrid array view.
Base.setindex!
— Method. setindex!(a::SysArray, v, ispec, inode)
Accessor for solution array.
Base.size
— Method.size(a::SubgridSysArrayView)
Return size of solution array view.
Base.size
— Method.size(a::SysArray)
Return size of solution array.
Base.view
— Method.view(a::SysArray{Tv},sg::SubGrid)
Create a view of the solution array on a subgrid.
TwoPointFluxFVM.bfacefactors!
— Method.bfacefactors!(grid::Grid,icell,nodefac)
Calculate node volume and voronoi surface contributions for boundary face.
TwoPointFluxFVM.bfacenode
— Method.bfacenode(grid::Grid,inode,ibface)
Index of boundary face node.
TwoPointFluxFVM.celledgenode
— Method.celledgenode(grid::Grid,inode,iedge,icell)
Index of cell edge node.
TwoPointFluxFVM.cellfactors!
— Method.cellfactors!(grid::Grid,icell,nodefac,edgefac)
Calculate node volume and voronoi surface contributions for cell.
TwoPointFluxFVM.cellnode
— Method.cellnode(grid,i,icell)
Return index of i-th local node in cell icell
TwoPointFluxFVM.dim_grid
— Method.dim_grid(grid)
Topological dimension of grid
TwoPointFluxFVM.dim_space
— Method.dim_space(grid)
Space dimension of grid
TwoPointFluxFVM.edgelength
— Method.function edgelength(edge::Edge)
Calculate the length of an edge.
TwoPointFluxFVM.geomspace
— Function.function geomspace(a::Real, b::Real, ha::Real, hb::Real, tol=1.0e-10)
(Try to) create a subdivision of interval (a,b) stored in the returned array X such that
X[1]==a, X[end]==b
(X[2]-X[1])<=ha+tol*(b-a)
(X[end]-X[end-1])<=hb+tol*(b-a)
- There is a number q such that
X[i+1]-X[i] == q*(X[i]-X[i-1])
- X is the array with the minimal possible number of points with the above property
Caveat: the algorithm behind this is well tested but unproven.
Returns an Array containing the points of the subdivision.
TwoPointFluxFVM.inidirichlet!
— Method.function inidirichlet!(this::System,U)
Initialize dirichlet boundary values for solution.
TwoPointFluxFVM.is_boundary_species
— Method.function is_boundary_species(this::System, ispec::Integer)
Check if species number corresponds to boundary species.
TwoPointFluxFVM.is_bulk_species
— Method.function is_bulk_species(this::System, ispec::Integer)
Check if species number corresponds bulk species.
TwoPointFluxFVM.nodecoord
— Method.nodecoord(grid,inode)
Return view of coordinates of node inode
.
TwoPointFluxFVM.num_bfaceregions
— Method.num_bfaceregions(grid::Grid)
Number of boundary face regions in grid.
TwoPointFluxFVM.num_cellregions
— Method.num_cellregions(grid::Grid)
Number of cell regions in grid.
TwoPointFluxFVM.num_dof
— Method.num_dof(this::System)
Number of degrees of freedom for system.
TwoPointFluxFVM.num_edges_per_cell
— Method.num_edges_per_cell(grid::Grid)
Number of edges per grid cell.
TwoPointFluxFVM.num_nodes_per_bface
— Method.num_nodes_per_bface(grid::Grid)
Number of nodes per boundary face
TwoPointFluxFVM.num_nodes_per_cell
— Method.num_nodes_per_cell(grid)
Return number of nodes per cell in grid.
TwoPointFluxFVM.num_species
— Method.num_species(a::SysArray)
Number of species (size of first dimension) of solution array.
TwoPointFluxFVM.num_species
— Method.num_species(this::System)
Number of species in system
TwoPointFluxFVM.reg_bface
— Method.reg_bface(grid, ibface)
Boundary region number for boundary face
TwoPointFluxFVM.reg_cell
— Method.reg_cell(grid,icell)
Bulk region number for cell
TwoPointFluxFVM.values
— Method.values(a::SysArray)=a.node_dof
Array of values in solution array.
TwoPointFluxFVM.Dirichlet
— Constant.Constant to be used as boundary condition factor to mark Dirichlet boundary conditons.
TwoPointFluxFVM.fbernoulli_eps
— Constant.Constant for switch between Taylor series and full implementation