Functions

comododir

Comodo.comododirFunction
comododir()

Description

This function simply returns the string for the Comodo path. This is helpful for instance to load items, such as meshes, from the assets` folder.

source

slidercontrol

Comodo.slidercontrolFunction
slidercontrol(hSlider,ax)

Adds arrow key control to sliders

Description

This function adds arrow key control to Makie sliders. The inputs are the slider handle hSlider as well as the axis ax. If this function is called the slider can be advanced a step by pressing the right arrow, and returned one step by pressing the left arrow. When one presses and holds the right or left arrow key, the slider will continue to move (as fast as graphics updating is possible on your system) up to the end or start slider position respectively. Users may also use the up or down arrow keys. These function the same as the right and left arrow keys, however, rather than stopping at the slider extrema, the sliders position will "wrap" back to the start when advancing beyond the end position, and vice versa.

source

slider2anim

Comodo.slider2animFunction
slider2anim(fig::Figure,hSlider::Slider,fileName::String; backforth=true, duration=2)

Exports movies from slider based visualisations.

Description

Converts the effect of the slider defined by the slider handle hSlider for the figure fig to an animation/movie file

source

elements2indices

Comodo.elements2indicesFunction
elements2indices(F)

Returns the indices contained in F

Description

This function obtains the unique set of indices for the vertices (nodes) used by the the simplices defined by F. The vector F may contain any type of simplices. For instance the elements in F may be of the type GeometryBasics.TriangleFace or GeometryBasics.QuadFace (or any other) for surface mesh data. However, volumetric elements of any type are permitted. In essence this function simply returns unique(reduce(vcat,F)). Hence any suitable vector containing vectors of numbers permitted by reduce(vcat,F) is supported.

source

gridpoints

Comodo.gridpointsFunction
gridpoints(x::Vector{T}, y=x, z=x) where T<:Real

Returns 3D grids of points

Description

The gridpoints function returns a vector of 3D points which span a grid in 3D space. Points are defined as per the input ranges or range vectors. The output point vector contains elements of the type Point.

source

gridpoints_equilateral

Comodo.gridpoints_equilateralFunction
gridpoints_equilateral(xSpan,ySpan,pointSpacing::T; return_faces = false, rectangular=false) where T <: Real

Returns a "grid" of 3D points that are located on the corners of an equilateral triangle tesselation.

Description

This function returns 3D point data in the form of a Vector{Point{3,Float64}}. The point distribution is for an equilateral triangle tesselation. The input consists of the span in the x-, and y-direction, i.e. xSpan and ySpan respectively, as well as the desired pointSpacing. The "spans" should be vectors or tuples defining the minimum and maximum coordinates for the grid. The true point spacing in the x-direction is computed such that a nearest whole number of steps can cover the required distance. Next this spacing is used to create the equilateral triangle point grid. Although the xSpan is closely adhered to through this method, the ySpan is not fully covered. In the y-direction the grid does start at the minimum level, but may stop short of reaching the maximum y as it may not be reachable in a whole number of steps from the minimum. Optional arguments include return_faces (default is false), which will cause the function to return triangular faces F as well as the vertices V. Secondly the option rectangular will force the grid to conform to a rectangular domain. This means the "jagged" sides are forced to be flat such that all x-coordinates on the left are at the minimum in xSpan and all on the right are at the maximum in xSpan, however, this does result in a non-uniform spacing at these edges.

source

interp_biharmonic_spline

Comodo.interp_biharmonic_splineFunction
interp_biharmonic_spline(x::Union{Vector{T}, AbstractRange{T}},y::Union{Vector{T}, AbstractRange{T}},xi::Union{Vector{T}, AbstractRange{T}}; extrapolate_method=:linear,pad_data=:linear) where T<:Real

Interpolates 1D (curve) data using biharmonic spline interpolation

Description

This function uses biharmonic spline interpolation [1], which features radial basis functions. The input is assumed to represent ordered data, i.e. consequtive unique points on a curve. The curve x-, and y-coordinates are provided through the input parameters x and y respectively. The third input xi defines the sites at which to interpolate. Each of in the input parameters can be either a vector or a range.

References

  1. David T. Sandwell, Biharmonic spline interpolation of GEOS-3 and SEASAT altimeter data, Geophysical Research Letters, 2, 139-142, 1987. doi: 10.1029/GL014i002p00139
source

interp_biharmonic

Comodo.interp_biharmonicFunction
interp_biharmonic(x,y,xi)

Interpolates n-dimensional data using biharmonic spline interpolation

Description

This function uses biharmonic interpolation [1]. The input x should define a vector consisting of m points which are n-dimensional, and the input y should be a vector consisting of m scalar data values.

References

  1. David T. Sandwell, Biharmonic spline interpolation of GEOS-3 and SEASAT altimeter data, Geophysical Research Letters, 2, 139-142, 1987. doi: 10.1029/GL014i002p00139
source

nbezier

Comodo.nbezierFunction
nbezier(P,n)

Returns a Bezier spline for the control points P whose order matches the numbe of control points provided.

Description

This function returns n points for an m-th order Bézier spline, based on the m control points contained in the input vector P. This function supports point vectors with elements of the type AbstractPoint{3} (e.g. Point{3, Float64}) or Vector{Float64}.

source

lerp

Comodo.lerpFunction
lerp(x::Union{T,Vector{T}, AbstractRange{T}},y,xi::Union{T,Vector{T}, AbstractRange{T}}) where T <: Real

Linear interpolation

Description

This linearly interpolates (lerps) the input data specified by the sites x and data y at the specified sites xi.

source

dist

Comodo.distFunction
dist(V1,V2)

Computes n-dimensional Euclidean distances

Description

Function compute an nxm distance matrix for the n inputs points in V1, and the m input points in V2. The input points may be multidimensional, in fact they can be any type supported by the euclidean function of Distances.jl. See also: https://github.com/JuliaStats/Distances.jl

source

mindist

Comodo.mindistFunction
mindist(V1,V2; getIndex=false, skipSelf = false )

Returns nearest point distances

Description

Returns the closest point distance for the input points V1 with respect to the input points V2. If the optional parameter getIndex is set to true (false by default) then this function also returns the indices of the nearest points in V2 for each point in V1. For self-distance evaluation, i.e. if the same point set is provided twice, then the optional parameter skipSelf can be set t0 true (default is false) if "self distances" (e.g. the nth point to the nth point) are to be avoided.

source

unique_dict_index

Comodo.unique_dict_indexFunction
unique_dict_index(X::Union{Array{T},Tuple{T}}; sort_entries=false) where T <: Any

Returns unique values and indices

Description

Returns the unique entries in X as well as the indices for them. The optional parameter sort_entries (default is false) can be set to true if each entry in X should be sorted, this is helpful to allow the entry [1,2] to be seen as the same as [2,1] for instance.

source

unique_dict_index_inverse

Comodo.unique_dict_index_inverseFunction
unique_dict_index_inverse(X::Union{Array{T},Tuple{T}}; sort_entries=false) where T <: Any

Returns unique values, indices, and inverse indices

Description

Returns the unique entries in X as well as the indices for them and the reverse indices to retrieve the original from the unique entries. The optional parameter sort_entries (default is false) can be set to true if each entry in X should be sorted, this is helpful to allow the entry [1,2] to be seen as the same as [2,1] for instance.

source

unique_dict_index_count

Comodo.unique_dict_index_countFunction
unique_dict_index_count(X::Union{Array{T},Tuple{T}}; sort_entries=false) where T <: Any

Returns unique values, indices, and counts

Description

Returns the unique entries in X as well as the indices for them and the counts in terms of how often they occured. The optional parameter sort_entries (default is false) can be set to true if each entry in X should be sorted, this is helpful to allow the entry [1,2] to be seen as the same as [2,1] for instance.

source

unique_dict_index_inverse_count

Comodo.unique_dict_index_inverse_countFunction
unique_dict_index_inverse_count(X::Union{Array{T},Tuple{T}}; sort_entries=false) where T <: Any

Returns unique values, indices, inverse indices, and counts

Description

Returns the unique entries in X as well as the indices for them and the reverse indices to retrieve the original from the unique entries, and also the counts in terms of how often they occured. The optional parameter sort_entries (default is false) can be set to true if each entry in X should be sorted, this is helpful to allow the entry [1,2] to be seen as the same as [2,1] for instance.

source

unique_dict_count

Comodo.unique_dict_countFunction
unique_dict_count(X::Union{Array{T},Tuple{T}}; sort_entries=false) where T <: Any

Returns unique values and counts

Description

Returns the unique entries in X as well as the counts in terms of how often they occured. The optional parameter sort_entries (default is false) can be set to true if each entry in X should be sorted, this is helpful to allow the entry [1,2] to be seen as the same as [2,1] for instance.

source

unique_dict_inverse

Comodo.unique_dict_inverseFunction
unique_dict_inverse(X::Union{Array{T},Tuple{T}}; sort_entries=false) where T <: Any

Returns unique values and inverse indices

Description

Returns the unique entries in X as well as the reverse indices to retrieve the original from the unique entries. The optional parameter sort_entries (default is false) can be set to true if each entry in X should be sorted, this is helpful to allow the entry [1,2] to be seen as the same as [2,1] for instance.

source

unique_dict

Comodo.unique_dictFunction
unique_dict(X::AbstractVector{T}) where T <: Real

Returns unique values, indices, and inverse indices. Uses an OrderedDict.

Description

Returns the unique entries in X as well as the indices for them and the reverse indices to retrieve the original from the unique entries.

source

gunique

Comodo.guniqueFunction
gunique(X; return_unique=true, return_index=false, return_inverse=false, return_counts=false, sort_entries=false)

Returns unique values and allows users to choose if they also want: sorting, indices, inverse indices, and counts.

Description

Returns the unique entries in X. Depending on the optional parameter choices the indices for the unique entries, the reverse indices to retrieve the original from the unique entries, as well as counts in terms of how often they occured, can be returned. The optional parameter sort_entries (default is false) can be set to true if each entry in X should be sorted, this is helpful to allow the entry [1,2] to be seen as the same as [2,1] for instance.

source

unique_simplices

Comodo.unique_simplicesFunction
unique_simplices(F,V=nothing)

Returns unique simplices (such as faces), independant of node order

Description

Returns the unique simplices in F as well as the indices of the unique simplices and the reverse indices to retrieve the original faces from the unique faces. Entries in F are sorted such that the node order does not matter.

source

ind2sub

Comodo.ind2subFunction
ind2sub(siz,ind)

Converts linear indices to subscript indices.

Description

Converts the linear indices in ind, for a matrix/array with size siz, to the equivalent subscript indices.

source

sub2ind

Comodo.sub2indFunction
sub2ind(siz,A)

Converts subscript indices to linear indices.

Description

Converts the subscript indices in A, for a matrix/array with size siz, to the equivalent linear indices.

source

meshedges

Comodo.meshedgesFunction
meshedges(F::Array{NgonFace{N,T},1}; unique_only=false) where N where T<:Integer

Returns a mesh's edges.

Description

This function returns the edges E for the input faces defined by F. The input F can either represent a vector of faces or a GeometryBasics.Mesh. The convention is such that for a face referring to the nodes 1-2-3-4, the edges are 1-2, 2-3, 3-4, 4-1.

source

icosahedron

Comodo.icosahedronFunction
icosahedron(r=1.0)

Creates an icosahedron mesh.

Description

Creates a GeometryBasics.Mesh for an icosahedron with radius r. The default radius, when not supplied, is 1.0.

source

octahedron

Comodo.octahedronFunction
octahedron(r=1.0)

Creates an octahedron mesh.

Description

Creates a GeometryBasics.Mesh for an octahedron with radius r. The default radius, when not supplied, is 1.0.

source

dodecahedron

Comodo.dodecahedronFunction
dodecahedron(r=1.0)

Creates a dodecahedron mesh.

Description

Creates a GeometryBasics.Mesh for an dodecahedron with radius r. The default radius, when not supplied, is 1.0.

source

cube

Comodo.cubeFunction
cube(r=1.0)

Creates a cube mesh.

Description

Creates a GeometryBasics.Mesh for an cube with radius r. The default radius, when not supplied, is 1.0.

source

tetrahedron

Comodo.tetrahedronFunction
tetrahedron(r=1.0)

Creates a tetrahedron mesh.

Description

Creates a GeometryBasics.Mesh for an tetrahedron with radius r. The default radius, when not supplied, is 1.0.

source

platonicsolid

Comodo.platonicsolidFunction
platonicsolid(n,r=1.0)

Returns a platonic solid mesh.

Description

Creates a GeometryBasics mesh description for a platonic solid of choice. The input n defines the choice.

  1. tetrahedron
  2. cube
  3. octahedron
  4. icosahedron
  5. dodecahedron

The final input parameter r defines the radius of the platonic solid (the radius of the circumsphere to the vertices). The default radius, when not supplied, is 1.0.

Arguments

n::Integer, defining platonic solid type r::Float64, defining circumsphere radius

source

tofaces

Comodo.tofacesFunction
tofaces(FM::Vector{Vector{TF}}) where TF<:Integer
tofaces(FM::Matrix{TF})  where TF<:Integer
tofaces(FM::Vector{NgonFace{m, OffsetInteger{-1, TF}}} ) where m where TF <: Integer
tofaces(FM::Vector{NgonFace{m, TF}} ) where m where TF <: Integer

Converts input to GeometryBasics compliant faces with standard integer types.

Description

The tofaces function converts "non-standard" (for Comodo) face set descriptions to "standard" ones. The following is considered such as standard: Vector{GeometryBasics.NgonFace{N,T}} where N where T<:Integer The input faces FM are converted to this format. FM can be of the following types:

  • FM::Vector{Vector{TF}} where TF<:Integer, whereby each Vector entry is

considered a face

  • FM::Matrix{TF} where TF<:Integer, whereby each row is considered a face
  • Vector{NgonFace{m, OffsetInteger{-1, TF}}} where TF<:Integer, whereby the

special integer type OffsetInteger{-1, TF} is converted to Int. If the intput is already of the right type this function leaves the input unchanged.

source

topoints

Comodo.topointsFunction
topoints(VM::Matrix{T}) where T<: Real
topoints(VM::Union{Array{Vec{N, T}, 1}, GeometryBasics.StructArray{TT,1} }) where TT <: AbstractPoint{N,T} where T <: Real where N   
topoints(VM::Vector{Vector{T}}) where T <: Real  
topoints(VM::Vector{Point{ND,TV}}) where ND where TV <: Real

Converts input to GeometryBasics compliant simple points without meta content.

Description

The topoints function converts the "non-standard" (for Comodo) input points defined by VM to the "standard" format: VM::Vector{Point{ND,TV}} where ND where TV <: Real. For matrix input each row is considered a point. For vector input each vector entry is considered a point.

source

togeometrybasics_mesh

Comodo.togeometrybasics_meshFunction
togeometrybasics_mesh

Converts the input to a GeometryBasics.Mesh

Description

This function converts the input faces F and vertices V to a GeometryBasics.Mesh. The function tofaces and topoints are used prior to conversion, to ensure standard faces and point types are used.

source

edgecrossproduct

Comodo.edgecrossproductFunction
edgecrossproduct(F,V::Vector{Point{ND,T}}) where ND where T<:Real  
edgecrossproduct(M::GeometryBasics.Mesh)

Returns the edge cross product, useful for nomal direction and area computations.

# Description

This function computes the so-called edge-cross-product for a input mesh that is either defined by the faces F and vertices V or the mesh M.

source

facenormal

Comodo.facenormalFunction
facenormal(F,V)

Returns the normal directions for each face.

Description

This function computes the per face normal directions for the input mesh defined either by the faces F and vertices V or by the GeometryBasics mesh M.

source

facearea

Comodo.faceareaFunction
facearea(F,V)

Returns the area for each face.

Description

This function computes the per face area for the input mesh defined either by the faces F and vertices V or by the GeometryBasics mesh M.

source

vertexnormal

Comodo.vertexnormalFunction
vertexnormal(F,V; weighting=:area)

Returns the surface normal at each vertex.

Description

This function computes the per vertex surface normal directions for the input mesh defined either by the faces F and vertices V or by the GeometryBasics mesh M. The optional parameter weighting sets how the face normal directions are averaged onto the vertices. If weighting=:none a plain average for the surrounding faces is used. If instead weighting=:area (default), then the average is weighted based on the face areas.

source

edgelengths

Comodo.edgelengthsFunction
edgelengths(E::LineFace,V)
edgelengths(F,V)
edgelengths(M::GeometryBasics.Mesh)

Returns edge lengths.

Description

This function computes the lengths of the edges defined by edge vector E (e.g as obtained from meshedges(F,V), where F is a face vector, and V is a vector of vertices. Alternatively the input mesh can be a GeometryBasics mesh M.

source

subtri

Comodo.subtriFunction
subtri(F,V,n; method = :linear)
subtri(F,V,n; method = :Loop)

Refines triangulations through splitting.

Description

The subtri function refines triangulated meshes iteratively. For each iteration each original input triangle is split into 4 triangles to form the refined mesh (one central one, and 3 at each corner). The following refinement methods are implemented:

method=:linear : This is the default method, and refines the triangles in a simple linear manor through splitting. Each input edge simply obtains a new mid-edge node.

method=:Loop : This method features Loop-subdivision [1,2]. Rather than linearly splitting edges and maintaining the original coordinates, as for the linear method, this method computes the new points in a special weighted sense such that the surface effectively approaches a "quartic box spline". Hence this method both refines and smoothes the geometry through spline approximation.

References

  1. Charles Loop, Smooth Subdivision Surfaces Based on Triangles, M.S. Mathematics Thesis, University of Utah. 1987.
  2. Jos Stam, Charles Loop, Quad/Triangle Subdivision, doi: 10.1111/1467-8659.t01-2-00647
source

subquad

Comodo.subquadFunction

subquad(F::Vector{NgonFace{4,TF}},V::Vector{Point{ND,TV}},n::Int; method=:linear) where TF<:Integer where ND where TV <: Real subquad(F::Vector{NgonFace{4,TF}},V::Vector{Point{ND,TV}},n::Int; method=:Catmull_Clark) where TF<:Integer where ND where TV <: Real

Refines quadrangulations through splitting.

Description

The subquad function refines quad meshes iteratively. For each iteration each original input quad is split into 4 smaller quads to form the refined mesh. The following refinement methods are implemented:

method=:linear : This is the default method, and refines the quads in a simple linear manor through splitting. Each input edge simply obtains a new mid-edge node, and each face obtains a new central node.

method=:Catmull_Clark : This method features Catmull_Clark-subdivision [1]. Rather than linearly splitting edges and maintaining the original coordinates, as for the linear method, this method computes the new points in a special weighted sense such that the surface effectively approaches a bicubic B-spline surface. Hence this method both refines and smoothes the geometry through spline approximation.

References

  1. E. Catmull and J. Clark, Recursively generated B-spline surfaces on arbitrary topological meshes, Computer-Aided Design, vol. 10, no. 6, pp. 350-355, Nov. 1978, doi: 10.1016/0010-4485(78)90110-0.
source

geosphere

Comodo.geosphereFunction
geosphere(n::Int,r::T; method=:linear) where T <: Real

Returns a geodesic sphere triangulation

Description

This function returns a geodesic sphere triangulation based on the number of refinement iterations n and the radius r. Geodesic spheres (aka Buckminster-Fuller spheres) are triangulations of a sphere that have near uniform edge lenghts. The algorithm starts with a regular icosahedron. Next this icosahedron is refined n times, while nodes are pushed to a sphere surface with radius r at each iteration. Two methods are available, i.e. :linear (default) and :Loop (see also subtri). The former features simply linear splitting while the latter features the Loop method which may produce a smoother result.

source

hexbox

Comodo.hexboxFunction
hexbox(boxDim,boxEl)

Returns a hexahedral mesh of a box

Description

This function returns a hexahedral mesh for a 3D rectangular box domain.

source

con_face_edge

Comodo.con_face_edgeFunction
con_face_edge(F,E_uni=nothing,indReverse=nothing)

Returns the edges connected to each face.

Description

This function computes the face-edge connectivity. The input faces F (and optionally also the unique edges E_uni and reverse indices indReverse to map to the non-unique edges, see also gunique) are used to create a list of edges connected to each face. If F contains N faces then the output contains N such lists. For triangles the output contains 3 edges per faces, for quads 4 per face and so on.

source

con_edge_face

Comodo.con_edge_faceFunction

conedgeface(F,E_uni=nothing,indReverse=nothing)

Returns the faces connected to each edge.

Description

This function computes the edge-face connectivity. The input faces F (and optionally also the unique edges E_uni and reverse indices indReverse to map to the non-unique edges, see also gunique) are used to create a list of faces connected to each edges. If E_uni contains N edges then the output contains N such lists. For non-boundary edges each edge should connect to 2 faces. Boundary edges connect to just 1 face.

source

con_face_face

Comodo.con_face_faceFunction
con_face_face(F,E_uni=nothing,indReverse=nothing,con_E2F=nothing,con_F2E=nothing)

Returns the edge-connected faces for each face.

Description

This function computes the face-face connectivity for each face. The input faces F are used to create a list of faces connected to each face by a shared edge. For non-boundary triangles for instance the output contains 3 edges per faces (which may be less for boundary triangles), and similarly non-boundary quads would each have 4 edge-connected faces. Additional optional inputs include: the unique edges E_uni, the reverse indices indReverse to map to the non-unique edges (see also gunique), as well as the edge-face con_E2F and face-edge con_F2E connectivity. These are all needed for computing the face-face connectivity and supplying them if already computed therefore save time.

source

con_face_face_v

Comodo.con_face_face_vFunction
con_face_face_v(F,V=nothing,con_V2F=nothing)

Returns the vertex-connected faces for each face.

Description

This function computes the face-face connectivity for each face. The input faces F are used to create a list of faces connected to each face by a shared vertex. Additional optional inputs include: the vertices V, and the vertex-face connectivity con_V2F. In terms of vertices only the number of vertices, i.e. length(V) is neede, if V is not provided it is assumed that length(V) corresponds to the largest index in F. The vertex-face connectivity if not supplied, will be computed by this function, hence computational time may be saved if it was already computed.

source

con_vertex_simplex

Comodo.con_vertex_simplexFunction
con_vertex_simplex(F,V=nothing)

Returns how vertices connect to simplices

Description

This function computes the vertex-simplex connectivity for each vertex. The input simplices F are used to create a list of simplices connected to each vertex. Additional optional inputs include: the vertices V. In terms of vertices only the number of vertices, i.e. length(V) is needed, if V is not provided it is assumed that length(V) corresponds to the largest index in F.

source

con_vertex_face

Comodo.con_vertex_faceFunction
con_vertex_face(F,V=nothing)

Returns how vertices connect to faces

Description

This function is an alias of con_vertex_simplex, and computes the vertex-face connectivity for each vertex. The input faces F are used to create a list of faces connected to each vertex. Additional optional inputs include: the vertices V. In terms of vertices only the number of vertices, i.e. length(V) is needed, if V is not provided it is assumed that length(V) corresponds to the largest index in F.

source

con_vertex_edge

Comodo.con_vertex_edgeFunction
con_vertex_edge(F,V=nothing)

Returns how vertices connect to edges

Description

This function is an alias of con_vertex_simplex, and computes the vertex-edge connectivity for each vertex. The input edges E are used to create a list of edges connected to each vertex. Additional optional inputs include: the vertices V. In terms of vertices only the number of vertices, i.e. length(V) is needed, if V is not provided it is assumed that length(V) corresponds to the largest index in E.

source

con_edge_edge

Comodo.con_edge_edgeFunction
con_edge_edge(E_uni,con_V2E=nothing)

Returns the vertex-connected edges for each edge.

Description

This function computes the edge-edge connectivity for each edge. The input edges F are used to create a list of edges connected to each edge by a shared vertex. Additional optional inputs include: con_V2E (the vertex-edge connectivity), which is instead computed when not provided.

source

con_vertex_vertex_f

Comodo.con_vertex_vertex_fFunction
con_vertex_vertex_f(F,V=nothing,con_V2F=nothing)

Returns the face-connected vertices for each vertex.

Description

This function computes the vertex-vertex connectivity for each vertex using the vertex connected faces. The input faces F are used to create a list of vertices connected to each vertex by a shared face. Additional optional inputs include: the vertices V and con_V2F (the vertex-face connectivity). In terms of vertices only the number of vertices, i.e. length(V) is needed, if V is not provided it is assumed that length(V) corresponds to the largest index in F. The vertex-face connectivity con_V2F is needed, hence is computed when not provided.

source

con_vertex_vertex

Comodo.con_vertex_vertexFunction
con_vertex_vertex(E,V=nothing,con_V2E=nothing)

Returns the edge-connected vertices for each vertex.

Description

This function computes the vertex-vertex connectivity for each vertex using the vertex connected edges. The input edges E are used to create a list of vertices connected to each vertex by a shared edge. Additional optional inputs include: the vertices V and con_V2E (the vertex-edge connectivity). In terms of vertices only the number of vertices, i.e. length(V) is needed, if V is not provided it is assumed that length(V) corresponds to the largest index in E. The vertex-edge connectivity con_V2E is needed, hence is computed when not provided.

source

meshconnectivity

Comodo.meshconnectivityFunction
meshconnectivity(F::Vector{NgonFace{N,TF}},V::Vector{Point{ND,TV}}) where N where TF<:Integer where ND where TV<:Real

Returns all mesh connectivity data

Description

This function returns the ConnectivitySet, i.e. all mesh connectivity data for the input mesh defined by the faces F and the vertices V. The ConnectivitySet contains the following connectivity descriptions:

  • face-edge
  • edge-face
  • face-face
  • face-face (wrt vertices)
  • vertex-face
  • vertex-edge
  • edge-edge
  • vertex-vertex
  • vertex-vertex (wrt faces)
source

mergevertices

Comodo.mergeverticesFunction
mergevertices(F::Vector{NgonFace{N,TF}},V::Vector{Point{ND,TV}}; roundVertices = true, numDigitsMerge=nothing) where N where TF<:Integer where ND where TV<:Real

Merges points that coincide

Description

This function take the faces F and vertices V and merges points that are sufficiently similar. Once points are merged the indices in F are corrected for the new reduced point set.

source

smoothmesh_laplacian

Comodo.smoothmesh_laplacianFunction
smoothmesh_laplacian(F,V,con_V2V=nothing; n=1, λ=0.5)

Description

This function implements weighted Laplacian mesh smoothing. At each iteration, this method replaces each point by an updated coordinate based on the mean coordinates of that point's Laplacian umbrella. The update features a lerp like weighting between the previous iterations coordinates and the mean coordinates. The code features Vs[q] = (1.0-λ).*Vs[q] .+ λ*mean(V[con_V2V[q]]) As can be seen, the weighting is controlled by the input parameter λ which is in the range (0,1). If λ=0 then no smoothing occurs. If λ=1 then pure Laplacian mean based smoothing occurs. For intermediate values a linear blending between the two occurs.

source

smoothmesh_hc

Comodo.smoothmesh_hcFunction
smoothmesh_hc(F::Vector{NgonFace{N,TF}},V::Vector{Point{ND,TV}}, n=1, α=0.1, β=0.5; con_V2V=nothing, tolDist=nothing, constrained_points=nothing) where N where TF<:Integer where ND where TV<:Real

Description

This function implements HC (Humphrey's Classes) smoothing [1]. This method uses Laplacian like smoothing but aims to compensate for shrinkage/swelling by also "pushing back" towards the original coordinates.

Reference

  1. Vollmer et al., Improved Laplacian Smoothing of Noisy Surface Meshes, 1999. doi: 10.1111/1467-8659.00334
source

quadplate

Comodo.quadplateFunction
quadplate(plateDim,plateElem; orientation=:up)

Returns a quad mesh for a plate

Description

This function creates a quadrilateral mesh (faces F and vertices V) for a plate. The dimensions in the x-, and y-direction are specified in the input vector plateDim, and the number of elements to use in each direction in the input vector plateElem.

source

quadsphere

Comodo.quadsphereFunction
quadsphere(n::Int,r::T) where T <: Real

Returns a quadrangulated sphere

Description

This function creates a quadrilateral mesh (faces F and vertices V) for a sphere with a radius defined by the input r. The input n defines the density of sphere mesh. The quad mesh is constructed using subquad subdivision of a regular cube, whereby n sets the number of splitting iterations to use. Using n=0 therefore returns a cube.

source

loflinear

Comodo.loftlinearFunction
loftlinear(V1,V2;num_steps=2,close_loop=true,face_type=:tri)

Loft a surface mesh between two input curves

Description

The loftlinear function spans a surface from input curve V1 to curve V2. The surface is formed by "lerping" curves from V1 to V2 in num_steps steps, and forming mesh faces between each curve. If close_loop==true then it is assumed the curves (and therefore the output surface mesh should be closed over, i.e. that a connection should be made between each curve end and start point. The user can request different face types for the output. The default is face_type=:tri which will form isoceles triangles (or equilateral triangles if the spacing is even) for a planar curve. The other face_type options supported are :quad (quadrilateral), and :tri_slash. For the latter, triangles are formed by slashing the quads.

Arguments:

  • V1::Vector: n-vector
  • V2::Vector: n-vector
source

pointspacingmean

Comodo.pointspacingmeanFunction
pointspacingmean(V::Vector{Point3{Float64}})
pointspacingmean(F::Array{NgonFace{N, Int}, 1},V::Vector{Point3{Float64}}) where N

The pointspacingmean function computes the mean spacing between points. The input can be just the coordinate set V, a vector of Point3 points, or also a set of edges E or faces F. If only V is provided it is assumed that V represents an ordered set of "adjacent" points, e.g. as for a curve. If a vector of edges E or a vector of faces F is also provided, then the average edge length is computed. If instead a set of facesF` is provided then edges are first computed after which the mean edge spacing is return.

source

ray_triangle_intersect

Comodo.ray_triangle_intersectFunction
ray_triangle_intersect(F::Vector{TriangleFace{Int}},V,ray_origin,ray_vector; rayType = :ray, triSide = 1, tolEps = eps(Float64))
ray_triangle_intersect(f::TriangleFace{Int},V,ray_origin,ray_vector; rayType = :ray, triSide = 1, tolEps = eps(Float64))

Description

This function can compute triangle-ray or triangle-line intersections through the use of the "Möller-Trumbore triangle-ray intersection algorithm" [1]. The required inputs are as follows:

F an single face or a vector of faces, e.g. Vector{TriangleFace{Int}} V The triangle vertices as a vector of points, i.e. Vector{Point{3, Float64}} ray_vector The ray vector which can be Vector{Point{3, Float64}} or Vec3{Float64}

The following optional input parameters can be provided: rayType = :ray (default) or :line. This defines wether the vector is treated as a ray (extends indefinately) or as a line (finite length) triSide = 1 (default) or 0 or -1. When triSide=1 only the inward intersections are considered, e.g. when the ray or line enters the shape (ray/line is pointing against face normal) When triSide=-1 only the outward intersections are considered, e.g. when the ray or line exits the shape (ray/line is pointing allong face normal) When triSide=0 both inward and outward intersections are considered. tolEps = eps(Float64) (default)

References

  1. Möller, Tomas; Trumbore, Ben (1997). Fast, Minimum Storage Ray-Triangle Intersection. Journal of Graphics Tools. 2: 21-28. doi: 10.1080/10867651.1997.10487468.
source

mesh_curvature_polynomial

Comodo.mesh_curvature_polynomialFunction
mesh_curvature_polynomial(F::Vector{TriangleFace{Int}},V::Vector{Point3{Float64}})
mesh_curvature_polynomial(M::GeometryBasics.Mesh)

Description

This function computes the mesh curvature at each vertex for the input mesh defined by the face F and the vertices V. A local polynomial is fitted to each point's "Laplacian umbrella" (point neighbourhood), and the curvature of this fitted form is derived. Instead of the mesh faces and vertices one may instead specify the GeometryBasics.Mesh M as the input.

The reference below [1] provides more detail on the algorithm. In addition, this implementation was created with the help of this helpful document, which features a nice overview of the theory/steps involved in this algorithm.

References

  1. F. Cazals and M. Pouget, Estimating differential quantities using polynomial fitting of osculating jets, Computer Aided Geometric Design, vol. 22, no. 2, pp. 121-146, Feb. 2005, doi: 10.1016/j.cagd.2004.09.004
source

separate_vertices

Comodo.separate_verticesFunction
separate_vertices(F::Array{NgonFace{N, Int}, 1},V::Array{Point{M, T}, 1}) where N where M where T<:Real
separate_vertices(M::GeometryBasics.Mesh)

This function takes the input mesh defined by the faces F and vertices V and separates any shared vertices. It does this by giving each face its own set of unshared vertices. Note that any unused points are not returned in the output point array Vn. Indices for the mapping are not created here but can simply be obtained using reduce(vcat,F).

source

evenly_sample

Comodo.evenly_sampleFunction
evenly_sample(V::Vector{Point{ND,TV}}, n::Int; rtol = 1e-8, niter = 1) where ND where TV<:Real

Evenly samples curves.

Description

This function aims to evenly resample the input curve defined by the ND points V using n points. The function returns the resampled points as well as the spline interpolator S used. The output points can also be retriebed by using: S.(range(0.0, 1.0, n)). Note that the even sampling is defined in terms of the curve length for a 4th order natural B-spline that interpolates the input data. Hence if significant curvature exists for the B-spline between two adjacent data points then the spacing between points in the output may be non-uniform (despite the allong B-spline distance being uniform).

source

invert_faces

Comodo.invert_facesFunction
invert_faces(F::Vector{NgonFace{N, TF}, 1}) where N where TF<:Integer

Flips face orientations.

Description

This function inverts the faces in F, such that the face normal will be flipped, by reversing the node order for each face.

source

kabsch_rot

sweeploft

Comodo.sweeploftFunction
F,V = sweeploft(Vc,V1,V2; face_type=:quad, num_twist = 0, close_loop=true)

Description

This function implements swept lofting. The start curve V1 is pulled allong the guide curve Vc while also gradually (linearly) morphing into the end curve V2. The optional parameter face_type (default :quad) defines the type of mesh faces uses. The same face types as loftlinear and extrudecurve are supported, i.e. :quad, :tri_slash, tri, or quad2tri. The optional parameter num_twist (default is 0) can be used to add an integer number (negative or positive) of full twists to the loft. Finally the optional parameter close_loop (default is true) determines if the section curves are deemed closed or open ended.

source

revolvecurve

Comodo.revolvecurveFunction
revolvecurve(Vc::Vector{Point{ND,TV}}; extent = 2.0*pi, direction=:positive, n=Vec{3, Float64}(0.0,0.0,1.0),num_steps=nothing,close_loop=true,face_type=:quad)  where ND where TV<:Real

Revolves curves to build surfaces

Description

This function rotates the curve Vc by the angle extent, in the direction defined by direction (:positive, :negative, :both), around the vector n, to build the output mesh defined by the faces F and vertices V.

source

batman

Comodo.batmanFunction
batman(n::Int)

Description

The batman function creates points on the curve for the Batman logo. The curve is useful for testing surface meshing algorithms since it contains sharp transitions and pointy features. The user requests n points on the curve. The default forces exactly n points which may result in an assymetric curve. To instead force symmetry the user can set the optional parameter symmetric=true. In this case the output will be symmetric allong the y-axis, however the number of points on the curve may have increased (if the input n is not even). The second optional input is the direction of the curve, i.e. if it is clockwise, dir=:cw or anti-clockwise dir=:acw. The implementation is based on a "parameterised Batman equation" 1. The following modifications where made, the curve is here centered around [0,0,0], scaled to be 2 in width, resampled evenly, and the default curve direction is anti-clockwise.

References

  1. https://www.desmos.com/calculator/ajnzwedvql
source

tridisc

Comodo.tridiscFunction
tridisc(r=1.0,n=0; ngon=6, method = :linear, orientation=:up)

Description

Generates the faces F and vertices V for a triangulated disc (circle). The algorithm starts with a triangulated hexagon (obtained if n=0) and uses iterative subtriangulation, and uses iterative subdivision (and pushing of boundary points to circular boundary) to obtain the final mesh. The subdivision method is an optional input, and is either :Loop (default) or :linear. Lastly the optional input orientation, which can be :up or :down sets the face normal direction.

source

regiontrimesh

Comodo.regiontrimeshFunction
regiontrimesh(VT,R,P)

Description

Generates a multi-region triangle mesh for the input regions. The boundary curves for all regions are containedin the tuple VT. Each region to be meshed is next defined using a tuple R containing indices into the curve typle VT. If an entry in R contains only one index then the entire curve domain is meshed. If R contains multiple indices then the first index is assumed to be for the outer boundary curve, while all subsequent indices are for boundaries defining holes in this region.

source

scalesimplex

Comodo.scalesimplexFunction
scalesimplex(F,V,s)

Scales faces (or general simplices) wrt their centre.

Description

This function scales each simplex (e.g. a face) wrt their centre (mean of coordinates). This function is useful in generating lattice structures from elements as well as to create visualisations whereby "looking into" the mesh is needed.

source

subcurve

Comodo.subcurveFunction
subcurve(V,n)

Adds n points between each curve point.

Description

This function adds n points between each current point on the curve V.

source

dualclad

Comodo.dualcladFunction
dualclad(F::Vector{NgonFace{N, TF}},V::Vector{Point{ND,TV}},s::T; connectivity=:face) where N where TF<:Integer where ND where TV<:Real where T<:Real

Returns a surface conforming dual lattice

Description

source

tet2hex

Comodo.tet2hexFunction
tet2hex(E::Vector{Tet4{T}},V::Vector{Point{ND,TV}}) where T<:Integer where ND where TV<:Real

Converts tetrahedra to hexahedra

Description

This function converts the input tetrahedra defined by the element set E and the vertex set V to a set of hexahedral elements Eh with vertices Vh. The conversion involves a splitting of each tetrahedron into 4 hexahedra.

source

element2faces

Comodo.element2facesFunction
element2faces(E::Vector{Element{N,T}}) where N where T

Returns element faces

Description

This function computes the faces for the input elements defined by E. The elements should be Vectors consisting of Tet4, Hex8 elements.

source

subhex

Comodo.subhexFunction
subhex(E::Vector{Hex8{T}},V::Vector{Point{ND,TV}},n::Int; direction=0) where T<:Integer where ND where TV<:Real

Split hexahedral elements

Description

This function splits the hexahedral elements defined by the elements E and vertices V. Splitting is done n times as requested. By default the splitting occurs in all direction (corresponding to the default direction=0). If instead direction is set to 1, 2, or 3, then the splitting only occur in the first, second or third local element direction respectively. Note that this direction depends on node order used. For a hexahedron where by nodes 1:4 are for the bottom, and nodes 5:8 are for the top of the element then the directions 1, 2, and 3 correspond to the x-, y-, and z-direction respectively.

source

rhombicdodecahedron

Comodo.rhombicdodecahedronFunction
rhombicdodecahedron(r = 1.0)

Creates mesh for rhombicdodecahedron

Description

This function creates the faces F and vertices V for a rhombicdodecahedron. The radius of the shape is set using the input radius r.

source

tri2quad

Comodo.tri2quadFunction
tri2quad(F,V; method=:split)

Converts triangles to quads

Description

This function converts the input triangular mesh, defined by the faces F and vertices V, to a quadrangulation. The method for this conversion is set using the attribute method which can be set to :split, splitting each triangle into 3 quads by introducing a new central node, or :rhombic. whereby each triangle edge is used to construct a rhombic quadrilateral face.

source

tetgenmesh

Comodo.tetgenmeshFunction
tetgenmesh(F::Array{NgonFace{N,TF}, 1},V::Vector{Point{3,TV}}; facetmarkerlist=nothing, V_regions=nothing,region_vol=nothing,V_holes=nothing,stringOpt="paAqYQ")  where N where TF<:Integer where TV<:Real

Creates a tetrahedral mesh

Description

This function uses the TetGen.jl library to mesh the input geometry defined by the faces F and the vertices V using tetrahedral elements. Several optional input parameters are available:

  • facetmarkerlist, a vector of integers with the same length as F and defines a face label for each face.
  • V_regions, a vector of points inside regions which require tetrahedral meshing.
  • region_vol, a vector of scalar values to denote the desired tetrahedral volume for each region.
  • V_holes, a vector of points inside holes (voids) that should remain empty.
  • stringOpt, the TetGen command string to use. See also the TetGen documentation.
source

surfacevolume

Comodo.surfacevolumeFunction
surfacevolume(F::Vector{NgonFace{N,TF}},V::Vector{Point{ND,TV}}) where N where TF<:Integer where ND where TV<:Real

Computes closed surface volume

Description

This function computes the volume of a closed surface defined by the faces F and the vertices V.

source

tetvolume

Comodo.tetvolumeFunction
tetvolume(E::Vector{Tet4{T}},V::Vector{Point{ND,TV}}) where T<:Integer where ND where TV<:Real

Computes tetrahedral volumes

Description

This function computes the volume for each tetrahedron defined by the input E, a vector of Tet4 elements, and V the point coordinates.

source

extrudefaces

Comodo.extrudefacesFunction
extrudefaces(F::Vector{NgonFace{NF,TF}},V::Vector{Point{ND,TV}}; extent=1.0, direction=:positive, num_steps=2, N::Union{Vector{Point{ND,TN}},Vector{Vec{ND, TN}},Nothing}=nothing) where NF where TF<:Integer where ND where TV<:Real where TN<:Real

Extrudes/thickens faces to form elements

Description

The inputs surface mesh, defined by the faces F and vertices V is extruded to create volumetric elements. Quadrilateral and triangular input faces are supported. These extrude into hexahedral and pentahedral elements respectively. The following input parameters are defined:

  • extent<:Real (default = 1.0) the length of the extrusion
  • direction is a symbol that is either :positive (default), :negative, or :both.
  • N The extrusion vectors. The default is nothing in which case the local

vertex normals are used.

  • num_steps (default is 2) is the number of nodes in the extrude direction, the

number of elements in the extrude direction is therefore num_steps-1.

source

filletcurve

Comodo.filletcurveFunction
filletcurve(V::Vector{Point{NV,TV}}; rMax::Union{Vector{T},T,Nothing}=nothing, constrain_method = :max, n=25, close_loop = false, eps_level = 1e-6) where TV<:Real where NV where T<:Real

Fillets/rounds curves

Description

The function takes in a curve defined by the points V and applies filleting (or rounding) to each "corner" (i.e. a point between two neighbouring points). The maximum radius rMax is used as the largest possible radius to use. If this, radius is not possible (e.g. if input points are too close), then a lower Radius is used.

source

squircle

Comodo.squircleFunction
squircle(r::T,n::Int,τ=0.5; atol=1e-6, dir=:acw) where T <: Real

Creates the squircle curve

Description

This function returns n points on a squircle. The squircle curve is defined using the radius r, and the parameter τ. The latter controls the morphing between a circle (τ=0) and square (τ=1).

source

circlerange

Comodo.circlerangeFunction
circlerange(n::Int; dir=:acw, deg=false)

Creates circular angles

Description

This function returns n angles for an even circular distribution of points. The optional input dir can be set to :acw (default) resulting in an anti-clockwise set of angles, and can be set to :cw for a clockwise set of angles. Angles are returned in radians since deg is false by default. Using deg=true results in angles in degrees.

source

edgefaceangles

Comodo.edgefaceanglesFunction
edgefaceangles(F::Vector{NgonFace{NF,TF}},V::Vector{Point{ND,TV}}; deg=false) where NF where TF<:Integer where ND where TV<:Real

Computed angles between faces

Description

This function computes the angle between two faces for each unique edge in the mesh specified by F and the vertices V. If the input mesh consists of n unique edges then the output features n angles. For boundary edges, where no pair of faces exists, the angle returned is NaN. The default behaviour results in angles being computed in radians. However by specifying deg=true the user can request degrees instead. Finally two additional outputs are created, namely the unique edge vector E_uni as well as the edge-to-face connectivity vector ,con_E2F.

source

faceanglesegment

Comodo.faceanglesegmentFunction
faceanglesegment(F::Vector{NgonFace{NF,TF}},V::Vector{Point{ND,TV; deg=false, angleThreshold = pi/8, indStart = 1)  where NF where TF<:Integer where ND where TV<:Real

Segments surfaces using face angles

Description

This function takes in a surface mesh defined by the faces F and the vertices V, and segments the surface mesh based on face angles. The output consists of a "feature label" vector G (a Vector{Int}, with the same length of F) whereby adjacent faces whosee angle is smaller or equal to the angleThreshold (default is pi/8) receive the same label. Hence this function allows one to find all faces with a similar orientation, for instance all top or side faces of some mesh. The function uses radians by default. However, buy specifying the optional parameter deg=true the user request that angles and the angleThreshold are in degrees.

source