API Documentation

API Documentation

struct SubGrid{Tc} <: AbstractGrid

Subgrid of parent grid (mainly for visualization purposes). Intended to hold support of species which are not defined everywhere.

function SubGrid(parent::Grid, 
                 subregions::AbstractArray; 
                 transform::Function=copytransform!,
                 boundary=false)

Create subgrid of list of regions.

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.

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.

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.

function dof(a::SysArray,ispec, inode)

Get number of degree of freedom. Return 0 if species is not defined in node.

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.

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.

function getdof(a::SysArray,i::Integer)

Return value for degree of freedom.

function integrate(this::System,F::Function,U)

Integrate solution vector over domain. Returns an Array{Int64,1} containing the integral for each species.

num_bfaces(grid::Grid)

Number of boundary faces in grid.

num_cells(grid)

Number of cells in grid

num_nodes(grid)

Number of nodes in grid

num_nodes(a::SysArray)

Number of nodes (size of second dimension) of solution array.

function setdof!(a::SysArray,v,i::Integer)

Set value for degree of freedom.

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.

function unknowns(system)

Create a solution vector for system.

TwoPointFluxFVM.valueFunction.
value(x)

Extract value from dual number. Use to debug physics callbacks. Re-exported from ForwardDiff.jl

abstract type AbstractGrid

Abstract type for grid like datastructures.

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}
mutable struct Grid

Structure holding grid data.

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
 |--|-----|-----|-----|-----|-----|-----|-----|--|
Grid(X::Array{Tc,1},X::Array{Tc,1})

Constructor for 2D grid from coordinate arrays. Boundary region numbers count counterclockwise:

locationnumber
south1
east2
north3
west4

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
NewtonControl()

Default constructor

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

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

mutable struct System{Tv}

Main structure holding data for system solution.

function  System(grid::Grid, physics::Any, maxspec::Integer)

Constructor for System. physics provides some user data, maxspec is the maximum number of species.

Base.copyMethod.
copy(this::SysArray)

Create a copy of solution array

Base.eltypeMethod.
eltype(grid)

Return element type of grid coordinates.

Base.getindexMethod.
getindex(aview::SubgridSysArrayView,ispec::Integer,inode::Integer)

Accessor method for subgrid array view.

Base.getindexMethod.
 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.sizeMethod.
size(a::SubgridSysArrayView)

Return size of solution array view.

Base.sizeMethod.
size(a::SysArray)

Return size of solution array.

Base.viewMethod.
view(a::SysArray{Tv},sg::SubGrid)

Create a view of the solution array on a subgrid.

bfacefactors!(grid::Grid,icell,nodefac)

Calculate node volume and voronoi surface contributions for boundary face.

bfacenode(grid::Grid,inode,ibface)

Index of boundary face node.

celledgenode(grid::Grid,inode,iedge,icell)

Index of cell edge node.

cellfactors!(grid::Grid,icell,nodefac,edgefac)

Calculate node volume and voronoi surface contributions for cell.

cellnode(grid,i,icell)

Return index of i-th local node in cell icell

dim_grid(grid)

Topological dimension of grid

dim_space(grid)

Space dimension of grid

function edgelength(edge::Edge)

Calculate the length of an edge.

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.

function inidirichlet!(this::System,U)

Initialize dirichlet boundary values for solution.

function is_boundary_species(this::System, ispec::Integer)

Check if species number corresponds to boundary species.

function is_bulk_species(this::System, ispec::Integer)

Check if species number corresponds bulk species.

nodecoord(grid,inode)

Return view of coordinates of node inode.

num_bfaceregions(grid::Grid)

Number of boundary face regions in grid.

num_cellregions(grid::Grid)

Number of cell regions in grid.

num_dof(this::System)

Number of degrees of freedom for system.

num_edges_per_cell(grid::Grid)

Number of edges per grid cell.

num_nodes_per_bface(grid::Grid)

Number of nodes per boundary face

num_nodes_per_cell(grid)

Return number of nodes per cell in grid.

num_species(a::SysArray)

Number of species (size of first dimension) of solution array.

num_species(this::System)

Number of species in system

reg_bface(grid, ibface)

Boundary region number for boundary face

reg_cell(grid,icell)

Bulk region number for cell

values(a::SysArray)=a.node_dof

Array of values in solution array.

Constant to be used as boundary condition factor to mark Dirichlet boundary conditons.

Constant for switch between Taylor series and full implementation