Periodic functions

To make a function f periodic with respect to a (partially) periodic boundary b::Boundary or geometry VG::VoronoiGeometry use the following

f2 = PeriodicFunction(f::Function,b::Boundary)
f2 = PeriodicFunction(f::Function,VG::VoronoiGeometry)

Step functions

Use one of the following methods to create a step function on the Voronoi grid:

f = StepFunction(VG::VoronoiGeometry, u<:AbstractVector; tree::Union{VoronoiKDTree,KDTree})
f = StepFunction(VG::VoronoiGeometry, u::Function; tree::Union{VoronoiKDTree,KDTree})
f = StepFunction(VG::VoronoiGeometry; tree::Union{VoronoiKDTree,KDTree})
f = StepFunction(nodes::VoronoiNodes, u<:AbstractVector; tree::Union{VoronoiKDTree,KDTree}=KDTree(nodes))

This yields a step function f that is constant on every cell of the VoronoiGeometry VG or on the Voronoi tessellation given by nodes. If u is an abstract vector, the value f(x)=u[i] is assigned if - according to tree - the nearest neighbor of x is the i-th node of VG or nodes. If no value for u is provided, StepFunction will retrieved bulk-integral data stored in VG. If VG has no bulk-data, the step-function will return nothing.

tree can be a KDTree from NearestNeighbors.jl or a VoronoiKDTree. It is highly recommended to use the last one as it accounts for periodicity.

Finally, consider the following advanced code:

    # create a composed function for integration
    f = FunctionComposer(reference_argument = [0.0,0.0], super_type = Float64, alpha = x->norm(x)*x, beta = x->sum(abs,x) )
    # create a VoronoiGeometry and integrate :alpha, :beta
    VG = VoronoiGeometry(VoronoiNodes(rand(2,40)), cuboid(2,periodic=[1]), integrator=HighVoronoi.VI_MONTECARLO, integrand=f.functions)
    # make a step function from integrated values:
    f_all = StepFunction(VG)
    # retrieve the alpha and beta- components as a single (real) valued stepfunctions
    alpha_step = x-> HighVoronoi.decompose(f, f_all(x),scalar=true)[:alpha]
    beta_step = x-> HighVoronoi.decompose(f, f_all(x),scalar=true)[:beta]
    # generate some sample output
    println(alpha_step([0.5,0.5]))
    println(beta_step([0.5,0.5]))

VoronoiKDTree

vt = VoronoiKDTree(VG::VoronoiGeometry; restrict_to_periodic=true)

This will create a KDTree that accounts for periodicity of VG. restrict_to_periodic=true implies that the "official" nodes are used only. It is highly recommended not to change this option if you are not knowing what you are doing.

Diameters of cells

f = DiameterFunction(VG::VoronoiGeometry; tree = VoronoiKDTree(VG))

This yields $f(x):=(r,R)$ where r is the inner and R the outer radius of the Voronoi cell that contains x. This is by its nature a step function.

Functions on interfaces

It can be usefull to consider the integrated values over the interfaces of the Voronoi tessellation as a function. This is achieved by InterfaceFunction:

f = InterfaceFunction(VD::VoronoiData,range,symbol=nothing;scalar=true)

This takes the VD::VoronoiData and creates a function that locally takes the value VD.interface_integral[i][k] over the respective interface. The value $f(x)$ of the function is chosen according to the two nearest neighbors, hence there is ambiguity in points with more than 2 nearest neighbors.

  • range: This can be a FunctionComposer-opject, in which case symbol has to be provided. It can also be a:b or [a1,a2,...,aN] to take a subarray of the values. It can also be :all in which case the full vector of values is taken.
  • scalar: If true, then vectors with only one index will be returned as scalar values.
f = InterfaceFunction(VG::VoronoiGeometry,range,symbol=nothing;scalar=true)

Calculates the VoronoiData and calls the first instance of the method.

f = InterfaceFunction(VG::VoronoiGeometry)

Sets range to the full data range. Similar to the above example for StepFunctions one may consider the following setting:

    # create a composed function for integration
    f = FunctionComposer(reference_argument = [0.0,0.0], super_type = Float64, alpha = x->norm(x)*x, beta = x->sum(abs,x) )
    # create a VoronoiGeometry and integrate :alpha, :beta
    VG = VoronoiGeometry(VoronoiNodes(rand(2,40)), cuboid(2,periodic=[1]), integrator=HighVoronoi.VI_MONTECARLO, integrand=f.functions)
    # make a step function from integrated values:
    f_all = InterfaceFunction(VG)
    # retrieve the alpha and beta- components as a single (real) valued stepfunctions
    alpha_i = x-> HighVoronoi.decompose(f, f_all(x),scalar=true)[:alpha]
    beta_i = x-> HighVoronoi.decompose(f, f_all(x),scalar=true)[:beta]
    # generate some sample output
    println(alpha_i([0.5,0.5]))
    println(beta_i([0.5,0.5]))

Functions from Data

If you want to generate a function from various integrated data in your own way, you can call

HighVoronoi.FunctionFromDataFunction
FunctionFromData(args...)

comes in to variations:

FunctionFromData(vg::VoronoiGeometry,tree=VoronoiKDTree(vg),composer=nothing; function_generator)

generates a function

x->function_generator( data=VoronoiData(vg), composer=composer, _Cell=nearest_neighbor_of(x) )

from vg, tree and a function

function_generator(;data ,composer , _Cell )

which takes data::VoronoiData generated from vg, composer from above and _Cell::Int for the number of the current node and returns a Float64 or a Vector{Float64} or anything else if you do not plan to hand it over to the routines of HighVoronoi. You can access every entry of VoronoiData to generate the value you want to be associated with the Voronoi cell belonging to vd.nodes[_Cell].

FunctionFromData(vd::VoronoiData,tree::VoronoiKDTree,composer=nothing; function_generator)

basically does the same but takes a vd::VoronoiData and tree mandatorily and passes vd to the function_generator.

The FunctionComposer: Passing function arguments

Always glue functions with a FunctionComposer

The FunctionComposer is internally used to glue together real valued functions. Therefore, if a user wants to glue together functions and afterwards work with "glued" information generated from HighVoronoi, using FunctionComposer is the way unify internal and external calculations.

The FunctionComposer is the element implemented in HighVoronoi to concatenate several Float or Vector{Float} valued functions into one single Vector{Float}-valued function using vcat(...). It is built using a call of the following method.

HighVoronoi.FunctionComposerMethod
FunctionComposer(;reference_argument, super_type, _functions...)

The composer takes the following arguments:

  • _functions: This is a list of named funcions.
  • super_type: suppose your functions return values of type T and Vector{T} you should set super_type=T
  • reference_argument: Your functions take values of type Float and are well defined in 0.0? Then you can put e.g. 0.0 here. If your function accepts StaticArray{3,Float64} put e.g. SVector{3,Float64}([0.0,1.2,3.4])

A typical example would be

f = FunctionComposer(reference_argument = [0.0,0.0,0.0], super_type = Float64, alpha = x->norm(x)*x, beta = x->sum(abs,x) )

or:

myfunctions=(alpha = x->norm(x)*x, beta = x->sum(abs,x))
f = FunctionComposer(reference_argument = [0.0,0.0,0.0], super_type = Float64; myfunctions...  )

The latter has the advantage that you can define your set of functions once and for all and use it again and again ensuring you always have the same order in the arguments. This brings us to an important point:

Don't mess with the order of arguments

FunctionComposer takes the order of functions as given in the argument. That is if you make function calls

f1 = FunctionComposer(reference_argument = [0.0,0.0,0.0], super_type = Float64, alpha = exp, beta = sin  )
f2 = FunctionComposer(reference_argument = [0.0,0.0,0.0], super_type = Float64; beta = sin, alpha = exp  )    

the algorithm will create two different functions x->[exp(x),sin(x)] and x->[sin(x),exp(x)] and it will NOT be able to clear up the mess this creates....

Retrieving the full (combined) function

The full function is stored in the variable FunctionComposer.functions.

myfunctions=(alpha = x->norm(x)*x, beta = x->sum(abs,x))
f = FunctionComposer(reference_argument = [0.0,0.0,0.0], super_type = Float64; myfunctions...  )

myvalue = f.functions([1.2,3.4,5.6])

Decomposing the Composer

To retrieve single information from an array like myvalue in the last example, you can simply use the internal function HighVoronoi.decompose(...):

myfunctions=(alpha = x->norm(x)*x, beta = x->sum(abs,x))
f = FunctionComposer(reference_argument = [0.0,0.0,0.0], super_type = Float64; myfunctions...  )

myvalue = f.functions([1.2,3.4,5.6])

values = HighVoronoi.decompose(f, myvalue)

println(values[:alpha], values[:beta])

If you whish $1d$-vectors to be returned as scalars, try out this one:

values2 = HighVoronoi.decompose(f, myvalue, scalar=true)

println(values2[:alpha], values2[:beta])