API reference
Main.LagrangianVoronoi.PolygonMulti
— TypePolygonMulti(; x::RealVector, kwargs...)
Predefined Voronoi Polygon for Navier-Stokes equations with implicit viscosity, optimized for multiphase problems with high density ratios, including surface tension. Mandatory keyword variable x
is the position of the generating seed.
Main.LagrangianVoronoi.PolygonNS
— TypePolygonNS(; x::RealVector, kwargs...)
Predefined Voronoi Polygon for Navier-Stokes equations. Mandatory keyword variable x
is the position of the generating seed.
Main.LagrangianVoronoi.PolygonNSF
— TypePolygonNSF(; x::RealVector, kwargs...)
Predefined Voronoi Polygon for Navier-Stokes-Fourier equations. Mandatory keyword variable x
is the position of the generating seed.
Main.LagrangianVoronoi.@Euler_vars
— MacroEuler_vars()
A convenience macro which defines the Eulerian fluid variables. See examples/doubleshear.jl
for an example of usage.
Main.LagrangianVoronoi.bdary_friction!
— Methodbdary_friction!(grid::VoronoiGrid, vDirichlet::Function, dt::Float64)
Applies a viscous drag at boundaries. This is used for Dirichlet condition for velocity when the the velocity is tangent of the surface. It cannot be used for inflows or outflows. The function vDirichlet
should accept a single argument x::RealVector
and return a RealVector
corresponding to the velocity of the boundary.
Main.LagrangianVoronoi.find_D!
— Methodfind_D!(grid::VoronoiGrid)
Compute the velocity deformation tensor $D = \frac{1}{2}(\nabla v + \nabla v^T)$ and assign it to every polygon.
Main.LagrangianVoronoi.viscous_step!
— Methodviscous_step!(grid::VoronoiGrid, dt::Float64; artificial_viscosity::Bool = true)
Applies one forward viscous step of size dt
to all Voronoi polygons. It assumes that the velocity deformation tensor is already computed using the find_D!
function. It is also assumed that every polygon p
has its dynamic viscosity p.mu
assigned by the initial condition or otherwise. A keyword parameter controls whether Stone-Norman viscosity tensor is used to damp oscillations near shocks (yes by default).
Main.LagrangianVoronoi.FastVector
— TypeFastVector{T} <: AbstractVector{T}
A vector with preallocated memory than can only increase in size. Compared to normal vector, this improves the performance of mesh generation by about 50%. If this is already implemented somewhere, please let me know.
Main.LagrangianVoronoi.fourier_step!
— Methodfourier_step!(grid::VoronoiGrid, dt::Float64)
Update the energy of Voronoi polygon by Fourier heat conduction. This assumes that every polygon p
has its value of heat conductivity p.k
assigned by the initial condition or otherwise.
Main.LagrangianVoronoi.heat_from_bdary!
— Methodheat_from_bdary!(grid::VoronoiGrid, dt::Float64, T_bc::Function)
Update the internal energy of a Voronoi cell by considering heat flux from boundary. Function T_bc
specifies the temperature at that boundary and should return NaN
to indicate adiabatic walls. This currently works only for ideal gas.
Main.LagrangianVoronoi.ideal_temperature!
— Methodideal_temperature!(grid::VoronoiGrid)
Assign temperature to every Voronoi polygon based on ideal gas law.
Main.LagrangianVoronoi.MAT0
— ConstantMAT0
A static 2x2 zero matrix.
Main.LagrangianVoronoi.MAT1
— ConstantMAT1
A static 2x2 identity matrix.
Main.LagrangianVoronoi.VEC0
— ConstantVEC0
Static zero vector. Equivalent to zero(RealVector)
.
Main.LagrangianVoronoi.VECNULL
— ConstantVECNULL
Undefined vector.
Main.LagrangianVoronoi.VECX
— ConstantVECX
Static cartesian basis vector in the X direction. Equivalent to RealVector(1,0)
Main.LagrangianVoronoi.VECY
— ConstantVECY
Static cartesian basis vector in the Y direction. Equivalent to RealVector(0,1)
.
Main.LagrangianVoronoi.Edge
— TypeEdge(v1::RealVector, v2::RealVector; label::Int = 0)
Edge defined by two endpoints and a label which indicates the neighbor's index. Zero or negative indices are used for domain boundaries. Edges are stack-allocated and immutable.
Main.LagrangianVoronoi.RealMatrix
— TypeRealMatrix(x11::Float64, x21::Float64, x12::Float64, x22::Float64)
Static Float64 matrix with 2x2 elements. Be warned that Julia has column-major ordering of elements!
Main.LagrangianVoronoi.RealVector
— TypeRealVector(x1::Float64, x2::Float64)
Static Float64 vector with 2 elements.
Main.LagrangianVoronoi.Rectangle
— TypeRectangle(xmin::RealVector, xmax::RealVector)
A rectangle aligned with the coordinate system defined by bottomleft and topright corner.
Main.LagrangianVoronoi.UnitRectangle
— MethodUnitRectangle()
The unit rectange (0,1)x(0,1).
Main.LagrangianVoronoi.area
— Methodarea(r::Rectangle)::Float64
Unoriented area of a rectangle.
Main.LagrangianVoronoi.cross2
— Methodfunction cross2(a::RealVector, b::RealVector)::Float64
A cross product in 2d. Returns a scalar equal to the z component of the corresponding 3d cross product
Main.LagrangianVoronoi.isinside
— Methodisinside(r::Rectangle, x::RealVector)::Bool
Return true
iff x
lies inside r
.
Main.LagrangianVoronoi.isnullvector
— Methodisnullvector(x::RealVector)::Bool
Returns true iff the vector is VECNULL.
Main.LagrangianVoronoi.len
— Methodlen(e::Edge)::Float64
Return the length of an edge.
Main.LagrangianVoronoi.midpoint
— Methodmidpoint(e::Edge)::RealVector
The edge midpoint.
Main.LagrangianVoronoi.midpoint
— Methodmidpoint(x::RealVector, y::RealVector)::RealVector
Midpoint inbetween two points.
Main.LagrangianVoronoi.norm_squared
— Methodnorm_squared(x::RealVector)::Float64
The squared norm of a vector. This is much faster than norm(x)^2
.
Main.LagrangianVoronoi.outer
— Methodouter(x::RealVector, y::RealVector)::RealMatrix
The outer product of two vectors.
Main.LagrangianVoronoi.verts
— Methodverts(r::Rectangle)::NTuple{4, RealVector}
Return the four vertices of a Rectangle in counter-clockwise order, starting with the bottomleft corner.
Main.LagrangianVoronoi.export_grid
— Methodexport_grid(grid::VoronoiGrid, filename::String, vars::Symbol...)
Export grid to a vtk file. Append datasets specified by vars
. Currently only numbers and vectors can be exported.
Main.LagrangianVoronoi.export_points
— Methodexport_points(grid::VoronoiGrid, filename::String, vars::Symbol...)
Export points to a vtk file as a cloud of particles. Append datasets specified by vars
.
Main.LagrangianVoronoi.boundaries
— Methodboundaries(p::VoronoiPolygon)
Create an iterator through all edges of 'p' that lie on the boundary of the domain. This can be used to implement boundary forces.
Main.LagrangianVoronoi.neighbors
— Methodneighbors(p::VoronoiPolygon, grid::VoronoiGrid)
Create an iterator through all Voronoi polygons that neighbor the polygon p
. Use it in a for loop to iterate through all triplets (q,e,y)
where
q
= the neighboring polygone
= the edge connectingp
andq
y
= the position ofq
as a neighbor ofp
(equivalent toq.x
for non-periodic grids)
Main.LagrangianVoronoi.move!
— Methodmove!(grid::VoronoiGrid, dt::Float64)
Update the positions of all polygons, moving each by dt*p.v
. The function ensures that points will never escape the computational domain (this would lead to undefined behavior). If the grid is peridoic, points will be wrapped around automatically. The grid is remeshed after this update.
Main.LagrangianVoronoi.CubicExpansion
— TypeCubicExpansion
Static vector to store coefficients of a cubic Taylor expansion.
Main.LagrangianVoronoi.LinearExpansion
— TypeLinearExpansion
Static vector to store coefficients of a linear Taylor expansion.
Main.LagrangianVoronoi.QuadraticExpansion
— TypeQuadraticExpansion
Static vector to store coefficients of a quadratic Taylor expansion.
Main.LagrangianVoronoi.integral
— Methodintegral(p::VoronoiPolygon, val::Float64, taylor::SVector)::Float64
Compute an integral of a polynomial over polygon p
defined by value val
at p.x
and the Taylor expansion taylor
.
Main.LagrangianVoronoi.movingls
— Methodmovingls(::Type{T}, grid::VoronoiGrid, p::VoronoiPolygon, fun; h, kernel)
Finds the Taylor expansion of a given function using moving least squares. The first argument T
should be one of following:
- LinearExpansion
- QuadraticExpansion
- CubicExpansion
Polygon p
is where the Taylor expansion of a function fun
is computed. A keyword argument h
is the moving radius. This needs to be sufficiently big, otherwise the Taylor expension may be undefined. Kernel is the weighting function used for defining the (weighted) least squared problem.
Main.LagrangianVoronoi.point_value
— Methodpoint_value(grid::VoronoiGrid, x::RealVector, fun::Function)
Obtain the interpolation of a function fun
defined on polygons at an arbitrary point x
. This is a slow code and should never be used in perfomance-critical areas.
Main.LagrangianVoronoi.poly_eval
— Methodpoly_eval(val::Float64, taylor::T, dx::RealVector)::Float64 where {T <: SVector}
Use this function to interpolate a value from the known Taylor expansion.
Main.LagrangianVoronoi.VoronoiPolygon
— TypeVoronoiPolygon
Abstract supertype for Voronoi Polygons. We call them Voronoi Polygons and not Voronoi cells to avoid confusion with cells of a cell list. No functions need to be defined for subtypes of Voronoi Polygon but we require that it has these fields
x::RealVector
edges::FastVector{Edge}
where x
is the generating seed and edges is a list of all edges in no particular order. The edges
can be initialized as edges = emptypolygon()
and are generated by the remesh!
function. The subtypes may contain additional fields depending on the physics involved. See src/celldefs.jl
for examples.
Main.LagrangianVoronoi.area
— Methodarea(p::VoronoiPolygon)::Float64
Return the unoriented area of a polygon.
Main.LagrangianVoronoi.centroid
— Methodcentroid(p::VoronoiPolygon)::RealVector
Return the centroid of a Voronoi polygon.
Main.LagrangianVoronoi.is_inside
— Methodis_inside(p::VoronoiPolygon, x::RealVector)::Bool
Determine whether a point x
lies inside the polygon p
.
Main.LagrangianVoronoi.isboundary
— Methodisboundary(e::Edge)::Bool
Return true iff edge e lies on domain boundary.
Main.LagrangianVoronoi.isboundary
— Methodisboundary(p::VoronoiPolygon)::Bool
Return true iff the Voronoi Polygon has at least one edge on domain boundary.
Main.LagrangianVoronoi.lr_ratio
— Methodlr_ratio(dx::RealVector, e::Edge)::Float64
The ratio between the length of a vector dx
and length of the edge e
. This may sound like a dubious function but is used so very often in this method. It is slightly faster than norm(dx)/len(e)
because it computes only one square root.
Main.LagrangianVoronoi.normal_vector
— Methodnormal_vector(e::Edge)::RealVector
The vector is rotated by 90% counter-clockwise with respect to the edge e
pointing from e.v1
to e.v2
.
Main.LagrangianVoronoi.surface_element
— Methodsurface_element(p::VoronoiPolygon)::Float64
The net surface element of a Voronoi Polygon. The resulting vector points outward and is zero for interior polygons.
Main.LagrangianVoronoi.tri_area
— Methodtri_area(a::RealVector, b::RealVector, c::RealVector)::Float64
Return the unsigned area of a triangle given by three points. The points can be in arbitrary order.
Main.LagrangianVoronoi.populate_circ!
— Methodpopulate_circ!(grid::VoronoiGrid{T}; charfun, center, ic!)
Populate computational domain with polygons arranged in concentric circles. Keyword parameters:
charfun
: the characteristic function; only those areas wherecharfun(x) = true
are populatedcenter
: the center of each circleic!
: the initial condition;ic!(p)
is called on every polygon after the mesh is generated
Main.LagrangianVoronoi.populate_hex!
— Methodpopulate_hex!(grid::VoronoiGrid{T}; charfun, ic!)
Populate computational domain with polygons arranged in hexagonal grid. Use this initializing method when you are not sure. Keyword parameters:
charfun
: the characteristic function; only those areas wherecharfun(x) = true
are populatedic!
: the initial condition;ic!(p)
is called on every polygon after the mesh is generated
Main.LagrangianVoronoi.populate_lloyd!
— Methodpopulate_lloyd!(grid::VoronoiGrid{T}; charfun, niterations, ic!)
Populate computational domain with polygons arranged in concentric circles. Keyword parameters:
charfun
: the characteristic function; only those areas wherecharfun(x) = true
are populatedcenter
: the center of each circleic!
: the initial condition;ic!(p)
is called on every polygon after the mesh is generated
Main.LagrangianVoronoi.populate_rand!
— Methodpopulate_rand!(grid::VoronoiGrid{T}; charfun, ic!)
Populate computational domain with polygons arranged randomly. This initializing method is not recommended because the mesh will have very low quality. Keyword parameters:
charfun
: the characteristic function; only those areas wherecharfun(x) = true
are populatedic!
: the initial condition;ic!(p)
is called on every polygon after the mesh is generated
Main.LagrangianVoronoi.populate_rect!
— Methodpopulate_rect!(grid::VoronoiGrid{T}; charfun, ic!)
Populate computational domain with Cartesian grid. Keyword parameters:
charfun
: the characteristic function; only those areas wherecharfun(x) = true
are populatedic!
: the initial condition;ic!(p)
is called on every polygon after the mesh is generated
Main.LagrangianVoronoi.populate_vogel!
— Methodpopulate_vogel!(grid::VoronoiGrid{T}; charfun, center, ic!)
Populate computational domain with polygons arranged on Vogel spiral. Keyword parameters:
charfun
: the characteristic function; only those areas wherecharfun(x) = true
are populatedcenter
: the center of the spiralic!
: the initial condition;ic!(p)
is called on every polygon after the mesh is generated
Main.LagrangianVoronoi.PressureOperator
— TypePressureOperator(grid::VoronoiGrid)
Construct the operator for pressure system. It is symmetric and positive-definite.
Main.LagrangianVoronoi.PressureSolver
— TypePressureSolver(grid::VoronoiGrid; verbose::Bool)
Construct the solver for pressure system.
Main.LagrangianVoronoi.eint
— Methodeint(p::VoronoiPolygon)::Float64
Return the internal energy of a Voronoi Polygon.
Main.LagrangianVoronoi.find_pressure!
— Functionfind_pressure!(solver::PressureSolver, dt::Float64, niter::Int64)
Find the estimated value of pressure field at the next time step. Integer niter
is the number fixed point iteratrions.
Main.LagrangianVoronoi.gravity_step!
— Methodgravity_step!(grid::VoronoiGrid, g::RealVector, dt::Float64)
Update velocity field and specific energy by a gravitational force.
Main.LagrangianVoronoi.ideal_eos!
— Functionideal_eos!(grid::VoronoiGrid, gamma = 1.4; Pmin)
Compute pressure and sound speed from internal energy and density using ideal gas equation of state. Number gamma
is adiabatic index and Pmin
can specify least possible value pressure. In the semi-implicit scheme, this value of pressure is used as an initial condition for an implicit solver.
Main.LagrangianVoronoi.pressure_step!
— Methodpressure_step!(grid::VoronoiGrid, dt::Float64)
Update the velocity and energy by the pressure field. This assumes that pressure was already determined.
Main.LagrangianVoronoi.stiffened_eos!
— Functionstiffened_eos!(grid::VoronoiGrid, gamma = 1.4; P0)
Compute pressure and sound speed from internal energy and density using ideal gas equation of state. Number gamma
is adiabatic index and P0
is the stiffness constant. Stiffened equation of state fares better for flows with very low Mach number.
Main.LagrangianVoronoi.MultiphaseProjector
— TypeMultiphaseProjector(grid::VoronoiGrid)
A linear operator which appears in the multiphase projection system. It is symmetric and positive-definite.
Main.LagrangianVoronoi.MultiphaseSolver
— TypeMultiphaseSolver(grid::VoronoiGrid, quality_threshold::Float64 = 0.25)
Construct a solver for the multiphase projection step. A keyword argument quality_threshold
determines when dv
should not be projected because the quality of the cell is too poor. In these problematic cases, stability of the simulation takes priority over physical accuracy.
Main.LagrangianVoronoi.find_dv!
— Functionfind_dv!(grid::VoronoiGrid, dt::Float64, alpha::Float64 = 1.0)
Determine dv
, the repair velocity. This is the first step of mesh relaxation procedure and is needed to maintain sufficiently high mesh quality. The repair velocity is chosen to be proportional to Frobenius norm of velocity deformation tensor D
and inversely proportional to the mesh quality squared. There is also multiplicative dimensionless parameter alpha
with default value alpha = 1
. Higher alpha
causes the Voronoi cells to be rounder (= better) but also less genuinly Lagrangian.
Main.LagrangianVoronoi.multiphase_projection!
— Methodmultiphase_projection!(solver::MultiphaseSolver)
Solve the linear system to find an orthogonal projection of dv
to a constraint space which gurantees conservation of area for every fluid phase. You only need this for multiphase flows. You can also run multiphase problems without it but expect some artificial density ridges between phases of different densities.
Main.LagrangianVoronoi.relaxation_step!
— Methodrelaxation_step!(grid::VoronoiGrid, dt::Float64; rusanov::Bool = true)
Update the variables by taking into account the repair velocity dv
. The repair velocity should be known before the relalaxation_step!
is called. See find_dv!
. The mass, momentum and energy of each cell is updated by solving one step of an advection problem. A keyword boolean argument rusanov
controls the use of a Rusanov approximate Riemann solver. The Rusanov flux is recommended to prevent the generation of new extrema (like negative density) and decrease of entropy. This is a second step of the mesh relaxation procedure (or maybe third, if your simulation involves a multiphase projector).
Main.LagrangianVoronoi.SimulationWorkspace
— TypeSimulationWorkspace
Abstract type for containing the computational grid, solvers and global variables. This methods should be defined for an instance sim
of SimulationWorkspace:
- step!(sim, time::Float64)
- postproc!(sim, time::Float64) [optional]
Main.LagrangianVoronoi.run!
— Methodrun!(sim::SimulationWorkspace, dt::Float64, t_end::Float64, step!::Function; kwargs...)
A utility function for time-marching and exporting the results. Specify the time step dt
, final time t_end
and a function $step!(sim::SimulationWorkspace, time::Float64).$ Keyword arguments:
nframes::Bool
(how many times should I export the data?)path::String
(where to export the data?)save_points::Bool
(should I export the point data?)save_grid:::Bool
(should I export the grid data?)save_csv::Bool
(should I export global variables as a csv data?)vtp_vars
(which local variables should I export?)csv_vars
(which global variables should I export?)postproc!::Function
(If provided,postproc!(sim, dt)
is called each time just before the data is saved.
Use this for some (expensive) post proccessing that affects the output files, not the simulation itself.)
Note this assumes constant time step. See examples/sedov.jl
for an example with adaptive time step.
Main.LagrangianVoronoi.ThreadedVec
— TypeThreadedVec
A mulithreaded vector which uses @batch for algebraic operations.
Main.LagrangianVoronoi.VoronoiGrid
— TypeVoronoiGrid{T}(boundary_rect::Rectangle, dr::Float64, kwargs...)
This struct contains all information about geometry, the Voronoi mesh and a cell list. Type variable T
specifies the Voronoi Polygon, boundary_rect
defines the computational domain and dr
is the default resolution (a particle will typically occupy an area of size dr^2
).
Keyword arguments:
xperiodic::Bool
: Is our domain periodic in the horizontal direction?yperiodic::Bool
: Is our domain periodic in the vertical direction?h
the size of cells in the cell_listr_max
the maximum possible distance between neighboring cells
Main.LagrangianVoronoi.get_arrow
— Methodget_arrow(x::RealVector, y::RealVector, grid::VoronoiGrid)::RealVector
Return a vector from y
to x
. Equal to x - y
on non-periodic domains. For periodic rectangle, the shortest possible vector is returned.
Main.LagrangianVoronoi.nearest_polygon
— Methodnearest_polygon(grid::VoronoiGrid{T}, x::RealVector)::T
Find the nearest polygon to the reference point x
.
Main.LagrangianVoronoi.remesh!
— Methodremesh!(grid::VoronoiGrid)
Generate all Voronoi cells by performing all Voronoi cuts. This is called automatically when the VoronoiGrid is populated and each the Voronoi generators move in space.