Postprocessing
Plotting
Plotting can be performed using the package GridVisualize.jl. This package extends the API with a couple of methods:
GridVisualize.gridplot
— Functiongridplot(sys::VoronoiFVM.AbstractSystem; kwargs...) -> Any
Plot grid behind system
GridVisualize.gridplot!
— Functiongridplot!(
vis,
sys::VoronoiFVM.AbstractSystem;
kwargs...
) -> Any
Plot grid behind system
GridVisualize.scalarplot
— Functionscalarplot(
sys::VoronoiFVM.AbstractSystem,
sol::AbstractMatrix;
species,
scale,
kwargs...
) -> Any
Plot one species from solution
scalarplot(
sys::VoronoiFVM.AbstractSystem,
sol::VoronoiFVM.AbstractTransientSolution;
species,
scale,
tscale,
kwargs...
) -> Any
Plot one species from transient solution
GridVisualize.scalarplot!
— Functionscalarplot!(
vis,
sys::VoronoiFVM.AbstractSystem,
sol::AbstractMatrix;
species,
scale,
kwargs...
) -> Any
Plot one species from solution
scalarplot!(
vis,
sys::VoronoiFVM.AbstractSystem,
sol::VoronoiFVM.AbstractTransientSolution;
species,
scale,
tscale,
tlabel,
kwargs...
) -> Any
Plot one species from transient solution
VoronoiFVM.plothistory
— Functionplothistory(tsol;
plots=[:timestepsizes,
:timestepupdates,
:newtonsteps,
:newtonupdates],
size=(700,600),
logmin=1.0e-20,
Plotter=GridVisualize.default_plotter(),
kwargs...)
Plot solution history stored in tsol
. The plots
argument allows to choose which plots are shown.
Grid verification
VoronoiFVM.nondelaunay
— Functionnondelaunay(grid;tol=1.0e-14, verbose=true)
Return non-Delaunay edges. Returns a vector of tuples: Each tuple consists of (node1, node2, edge factor, region)
If the vector has length 0, the grid is boundary conforming Delaunay with respect to each cell region. This means that up to tol
, all edge form factors are nonnegative.
Norms & volumes
LinearAlgebra.norm
— Functionnorm(system, u)
norm(system, u, p)
Calculate Euklidean norm of the degree of freedom vector.
VoronoiFVM.lpnorm
— Functionlpnorm(sys, u, p)
lpnorm(sys, u, p, species_weights)
Calculate weighted discrete $L^p$ norm of a solution vector.
VoronoiFVM.l2norm
— Functionl2norm(sys, u)
l2norm(sys, u, species_weights)
Calculate weigthed discrete $L^2(\Omega)$ norm of a solution vector.
VoronoiFVM.w1pseminorm
— Functionw1pseminorm(sys, u, p)
w1pseminorm(sys, u, p, species_weights)
Calculate weighted discrete $W^{1,p}(\Omega)$ seminorm of a solution vector.
VoronoiFVM.h1seminorm
— Functionh1seminorm(sys, u)
h1seminorm(sys, u, species_weights)
Calculate weighted discrete $H^1(\Omega)$ seminorm of a solution vector.
VoronoiFVM.w1pnorm
— Functionw1pnorm(sys, u, p)
w1pnorm(sys, u, p, species_weights)
Calculate weighted discrete $W^{1,p}(\Omega)$ norm of a solution vector.
VoronoiFVM.h1norm
— Functionh1norm(sys, u)
h1norm(sys, u, species_weights)
Calculate weighted discrete $H^1(\Omega)$ norm of a solution vector.
VoronoiFVM.lpw1pseminorm
— Functionlpw1pseminorm(sys, u, p)
lpw1pseminorm(sys, u, p, species_weights)
Calculate weighted discrete $L^p([0,T];W^{1,p}(\Omega))$ seminorm of a transient solution.
VoronoiFVM.l2h1seminorm
— Functionl2h1seminorm(sys, u)
l2h1seminorm(sys, u, species_weights)
Calculate weighted discrete $L^2([0,T];H^1(\Omega))$ seminorm of a transient solution.
VoronoiFVM.lpw1pnorm
— Functionlpw1pnorm(sys, u, p)
lpw1pnorm(sys, u, p, species_weights)
Calculate weighted discrete $L^p([0,T];W^{1,p}(\Omega))$ norm of a transient solution.
VoronoiFVM.l2h1norm
— Functionl2h1norm(sys, u)
l2h1norm(sys, u, species_weights)
Calculate weighted discrete $L^2([0,T];H^1(\Omega))$ norm of a transient solution.
VoronoiFVM.nodevolumes
— Functionnodevolumes(system)
Calculate volumes of Voronoi cells.
Solution integrals
VoronoiFVM.integrate
— Methodintegrate(system, tf, U; kwargs...)
Calculate test function integral for steady state solution.
VoronoiFVM.integrate
— Methodintegrate(system,F,U; boundary=false)
Integrate node function (same signature as reaction or storage) F
of solution vector region-wise over domain or boundary. The result is an nspec x nregion
matrix.
VoronoiFVM.integrate
— Methodintegrate(system,U; boundary=false)
Integrate solution vector region-wise over domain or boundary. The result is an nspec x nregion
matrix.
VoronoiFVM.integrate
— Methodintegrate(system,F,U; boundary=false)
Integrate node function (same signature as reaction or storage) F
of solution vector region-wise over domain or boundary. The result is an nspec x nregion
matrix.
VoronoiFVM.edgeintegrate
— Functionedgeintegrate(system,F,U; boundary=false)
Integrate edge function (same signature as flux function) F
of solution vector region-wise over domain or boundary. The result is an nspec x nregion
matrix.
Nodal flux reconstruction
VoronoiFVM.nodeflux
— Functionnodeflux(system, U; data)
Reconstruction of edge flux as vector function on the nodes of the triangulation. The result can be seen as a piecewiesw linear vector function in the FEM space spanned by the discretization nodes exhibiting the flux density.
The reconstruction is based on the "magic formula" R. Eymard, T. Gallouet, R. Herbin, IMA Journal of Numerical Analysis (2006) 26, 326−353, Lemma 2.4 .
The return value is a dim x nspec x nnodes
vector. The flux of species i can e.g. plotted via GridVisualize.vectorplot.
Example:
ispec=3
vis=GridVisualizer(Plotter=Plotter)
scalarplot!(vis,grid,solution[ispec,:],clear=true,colormap=:summer)
vectorplot!(vis,grid,nf[:,ispec,:],clear=false)
reveal(vis)
CAVEAT: there is a possible unsolved problem with the values at domain corners in the code. Please see any potential boundary artifacts as a manifestation of this issue and report them.
Boundary flux calculation
VoronoiFVM.TestFunctionFactory
— Typemutable struct TestFunctionFactory{Tu, Tv}
Data structure containing DenseSystem used to calculate test functions for boundary flux calculations.
Type parameters:
Tu
: value type of test functionsTv
: Default value type of systemsystem::VoronoiFVM.AbstractSystem{Tv} where Tv
: Original system
state::VoronoiFVM.SystemState{Tu, Tp, TMatrix, TSolArray} where {Tu, Tp, TMatrix<:AbstractMatrix{Tu}, TSolArray<:AbstractMatrix{Tu}}
: Test function system state
control::SolverControl
: Solver control
VoronoiFVM.TestFunctionFactory
— MethodTestFunctionFactory(
system::VoronoiFVM.AbstractSystem{Tv};
control
) -> TestFunctionFactory
Constructor for TestFunctionFactory from System
VoronoiFVM.integrate
— Methodintegrate(system, tf, U; kwargs...)
Calculate test function integral for steady state solution.
VoronoiFVM.integrate
— Methodintegrate(system, tf, U, Uold, tstep; params, data)
Calculate test function integral for transient solution.
VoronoiFVM.integrate_stdy
— Methodintegrate_stdy(system, tf, U; data)
Steady state part of test function integral.
VoronoiFVM.integrate_tran
— Methodintegrate_tran(system, tf, U; data)
Calculate transient part of test function integral.
VoronoiFVM.testfunction
— Methodtestfunction(
factory::TestFunctionFactory{Tv},
bc0,
bc1
) -> Any
Create testfunction which has Dirichlet zero boundary conditions for boundary regions in bc0 and Dirichlet one boundary conditions for boundary regions in bc1.
Impedance calculatiom
Impedance calculation can be seen as a postprocessing step after the solution of the unexcited stationary system.
VoronoiFVM.AbstractImpedanceSystem
— Typeabstract type AbstractImpedanceSystem{Tv<:Number}
Abstract type for impedance system.
VoronoiFVM.ImpedanceSystem
— Typemutable struct ImpedanceSystem{Tv} <: VoronoiFVM.AbstractImpedanceSystem{Tv}
Concrete type for impedance system.
sysnzval::AbstractArray{Complex{Tv}, 1} where Tv
: Nonzero pattern of time domain system matrix
storderiv::AbstractMatrix
: Derivative of storage term
matrix::AbstractArray{Complex{Tv}, 2} where Tv
: Complex matrix of impedance system
F::AbstractArray{Complex{Tv}, 2} where Tv
: Right hand side of impedance system
U0::AbstractMatrix
: Stationary state
VoronoiFVM.ImpedanceSystem
— MethodImpedanceSystem(system, U0, excited_spec, excited_bc)
Construct impedance system from time domain system sys
and steady state solution U0
under the assumption of a periodic perturbation of species excited_spec
at boundary excited_bc
.
VoronoiFVM.ImpedanceSystem
— MethodImpedanceSystem(system, U0; λ0)
Construct impedance system from time domain system sys
and steady state solution U0
CommonSolve.solve!
— Methodsolve!(UZ, impedance_system, ω)
Solve the impedance system for given frequency ω
.
VoronoiFVM.freqdomain_impedance
— Methodfreqdomain_impedance(
impedance_system,
ω,
U0,
excited_spec,
excited_bc,
excited_bcval,
dmeas_stdy,
dmeas_tran
)
Calculate reciprocal value of impedance.
- excitedspec,excitedbc,excited_bcval are ignored.
This is deprecated: use impedance
.
VoronoiFVM.impedance
— Methodimpedance(impedance_system,ω, U0 ,
excited_spec, excited_bc, excited_bcval,
dmeas_stdy,
dmeas_tran
)
Calculate impedance.
- ω: frequency
- U0: steady state slution
- dmeas_stdy: Derivative of steady state part of measurement functional
- dmeas_tran Derivative of transient part of the measurement functional
VoronoiFVM.measurement_derivative
— Methodmeasurement_derivative(system, measurement_functional, U0)
Calculate the derivative of the scalar measurement functional at steady state U0
Usually, this functional is a test function integral. Initially, we assume that its value depends on all unknowns of the system.
VoronoiFVM.unknowns
— Methodunknowns(impedance_system)
Create a vector of unknowns of the impedance system