LinearAlgebraicRepresentation.LARConstant
LAR = Union{ Tuple{Points, Cells},Tuple{Points, Cells, Cells} }

Alias declation of LAR-specific data structure. LAR is a pair (Geometry, Topology), where Geometry is stored as Points, and Topology is stored as Cells.

LinearAlgebraicRepresentation.CellsType
Cells = Array{Array{Int,1},1}

Alias declation of LAR-specific data structure. Dense Array to store the indices of vertices of P-cells of a cellular complex. The linear space of P-chains is generated by Cells as a basis. Simplicial P-chains have $P+1$ vertex indices for cell element in Cells array. Cuboidal P-chains have $2^P$ vertex indices for cell element in Cells array. Other types of chain spaces may have different numbers of vertex indices for cell element in Cells array.

LinearAlgebraicRepresentation.ChainType
Chain = SparseArrays.SparseVector{Int8,Int}

Alias declation of LAR-specific data structure. Binary SparseVector to store the coordinates of a chain of N-cells. It is nnz=1 with value=1 for the coordinates of an elementary N-chain, constituted by a single N-chain.

LinearAlgebraicRepresentation.ChainComplexType
ChainComplex = Array{ChainOp,1}

Alias declation of LAR-specific data structure. It is a 1-dimensional Array of ChainOp that provides storage for either the chain of boundaries (from D to 0) or the transposed chain of coboundaries (from 0 to D), with D the dimension of the embedding space, which may be either $R^2$ or $R^3$.

LinearAlgebraicRepresentation.ChainOpType
ChainOp = SparseArrays.SparseMatrixCSC{Int8,Int}

Alias declation of LAR-specific data structure. SparseMatrix in Compressed Sparse Column format, contains the coordinate representation of an operator between linear spaces of P-chains. Operators $P-Boundary : P-Chain -> (P-1)-Chain$ and $P-Coboundary : P-Chain -> (P+1)-Chain$ are typically stored as ChainOp with elements in ${-1,0,1}$ or in ${0,1}$, for signed and unsigned operators, respectively.

LinearAlgebraicRepresentation.LARmodelType
LARmodel = Tuple{Points,Array{Cells,1}}

Alias declation of LAR-specific data structure. LARmodel is a pair (Geometry, Topology), where Geometry is stored as Points, and Topology is stored as Array of Cells. The number of Cells values may vary from 1 to N+1.

LinearAlgebraicRepresentation.PointsType
Points = Array{Number,2}

Alias declation of LAR-specific data structure. Dense Array{Number,2,1}$M x N$ to store the position of vertices (0-cells) of a cellular complex. The number of rows $M$ is the dimension of the embedding space. The number of columns $N$ is the number of vertices.

LinearAlgebraicRepresentation.StructType
Struct

A container of geometrical objects is defined by applying the function Struct to the array of contained objects. Each value is defined in local coordinates and may be transformed by affine transformation tensors.

The value returned from the application of Struct to an Array{Any, 1} of LAR values, matrices, and Struct values is a value of Struct type. The coordinate system of this value is the one associated with the first object of the Struct arguments. Also, the resulting geometrical value is often associated with a variable name.

The generation of containers may continue hierarchically by suitably applying Struct. Notice that each LAR object in a Struct container is transformed by each matrix before it within the container, going from right to left. The action of a transformation (tensor) extends to each object on its right within its own container. Whereas, the action of a tensor does not extend outside its container, according to the semantics of PHIGS structures.

Example

julia> L = LinearAlgebraicRepresentation;

julia> assembly = L.Struct([L.sphere()(), L.t(3,0,-1), L.cylinder()()])
# return
Lar.Struct(Any[([0.0 -0.173648 … -0.336824 -0.17101; 0.0 0.0 … 0.0593912 0.0301537;
-1.0 -0.984808 … 0.939693 0.984808], Array{Int64,1}[[2, 3, 1], [4, 2, 3], [4, 3, 5], [4,
5, 6], [7, 5, 6], [7, 8, 6], [7, 9, 8], … , [1.0 0.0 0.0 3.0; 0.0 1.0 0.0 0.0; 0.0 0.0 1.0
-1.0; 0.0 0.0 0.0 1.0], ([0.5 0.5 … 0.492404 0.492404; 0.0 0.0 … -0.0868241 -0.0868241;
0.0 2.0 … 0.0 2.0], Array{Int64,1}[[4, 2, 3, 1], [4, 3, 5, 6], [7, 5, 8, 6], [7, 9, 10,
8], [9, 10, 11, 12], [13, 14, 11, 12], … , [68, 66, 67, 65], [68, 69, 67, 70], [71, 69,
72, 70], [71, 2, 72, 1]])], Array{Float64,2}[[-1.0; -1.0; -1.0], [3.5; 1.0; 1.0]],
"14417445522513259426", 3, "feature")

julia> assembly.name = "simple example"
# return
"simple example"

julia> assembly
# return
Lar.Struct(Any[([0.0 -0.173648 … -0.336824 -0.17101; 0.0 0.0 … 0.0593912 0.0301537;
-1.0 -0.984808 … 0.939693 0.984808], Array{Int64,1}[[2, 3, 1], [4, 2, 3], [4, 3, 5], [4,
5, 6], [7, 5, 6], [7, 8, 6], … , [71, 2, 72, 1]])], Array{Float64,2}[[-1.0; -1.0; -1.0],
[3.5; 1.0; 1.0]], "simple example", 3, "feature")

julia> using Plasm

julia> Plasm.view(assembly)
LinearAlgebraicRepresentation.surfaceType
sphere(radius=1., angle1=pi, angle2=2*pi)(shape=[18, 36])

Compute a cellular 2-complex, approximation of the two-dimensional closed surface, embedded in a three-dimensional Euclidean space. Geographical coordinates are user to compute the 0-cells of the complex.

Example

julia> GL.VIEW([
	GL.GLGrid( Lar.sphere()()..., GL.COLORS[1],0.75 ),
	GL.GLFrame
]);
source
LinearAlgebraicRepresentation.surfaceType
surface(P::Lar.LAR, signedInt::Bool=false)::Float64

surface integral on polyhedron P.

Example # unit 3D tetrahedron

julia> V = [0.0 1.0 0.0 0.0; 0.0 0.0 1.0 0.0; 0.0 0.0 0.0 1.0]
3×4 Array{Float64,2}:
 0.0  1.0  0.0  0.0
 0.0  0.0  1.0  0.0
 0.0  0.0  0.0  1.0

julia> FV = [[1, 2, 4], [1, 3, 2], [4, 3, 1], [2, 3, 4]]
4-element Array{Array{Int64,1},1}:
 [1, 2, 4]
 [1, 3, 2]
 [4, 3, 1]
 [2, 3, 4]

julia> P = V,FV
([0.0 1.0 0.0 0.0; 0.0 0.0 1.0 0.0; 0.0 0.0 0.0 1.0], 
Array{Int64,1}[[1, 2, 4], [1, 3, 2], [4, 3, 1], [2, 3, 4]])

julia> Lar.volume(P)
0.16666666666666674
LinearAlgebraicRepresentation.DFV_visitMethod
DFV_visit( VV::cells, out::Array, count::Int, visited::Array, parent::Array, d::Array, low::Array, stack::Array, u::Int )::Array

Hopcroft-Tarjan algorithm. Depth-First-Visit recursive algorithm.

LinearAlgebraicRepresentation.IIFunction
II(P::Lar.LAR, alpha::Int, beta::Int, gamma::Int, signedInt=false)

Basic integration function on 2D plane.

Example unit 3D triangle

julia> V = [0.0 1.0 0.0; 0.0 0.0 1.0; 0.0 0.0 0.0]
3×3 Array{Float64,2}:
 0.0  1.0  0.0
 0.0  0.0  1.0
 0.0  0.0  0.0

julia> FV = [[1,2,3]]
1-element Array{Array{Int64,1},1}:
 [1, 2, 3]

julia> P = V,FV
([0.0 1.0 0.0; 0.0 0.0 1.0; 0.0 0.0 0.0], Array{Int64,1}[[1, 2, 3]])

julia> Lar.II(P, 0,0,0)
0.5
LinearAlgebraicRepresentation.IIIMethod
III(P::Lar.LAR, alpha::Int, beta::Int, gamma::Int)::Float64

Basic integration function on 3D space.

Example # unit 3D tetrahedron

julia> V = [0.0 1.0 0.0 0.0; 0.0 0.0 1.0 0.0; 0.0 0.0 0.0 1.0]
3×4 Array{Float64,2}:
 0.0  1.0  0.0  0.0
 0.0  0.0  1.0  0.0
 0.0  0.0  0.0  1.0

julia> FV = [[1, 2, 4], [1, 3, 2], [4, 3, 1], [2, 3, 4]]
4-element Array{Array{Int64,1},1}:
 [1, 2, 4]
 [1, 3, 2]
 [4, 3, 1]
 [2, 3, 4]

julia> P = V,FV
([0.0 1.0 0.0 0.0; 0.0 0.0 1.0 0.0; 0.0 0.0 0.0 1.0], 
Array{Int64,1}[[1, 2, 4], [1, 3, 2], [4, 3, 1], [2, 3, 4]])

julia> Lar.III(P, 0,0,0)
0.16666666666666674
LinearAlgebraicRepresentation.INSRMethod
INSR(f::Function)(seq::Array{Any,1})::Any

FL primitive combinator to transform a binary function to an n-ary one.

julia> mod1D = Lar.grid(repeat([.1,-.1],outer=5)...)
([0.0 0.1 … 0.9 1.0], Array{Int64,1}[[1, 2], [3, 4], [5, 6], [7, 8], [9, 10]])

julia> using ViewerGL; GL = ViewerGL

julia> GL.VIEW([ GL.GLFrame2, GL.GLGrid(mod1D..., GL.COLORS[1],1) ])

julia> mod3D = Lar.INSR(Lar.larModelProduct)([mod1D,mod1D,mod1D])
([0.0 0.0 … 1.0 1.0; 0.0 0.0 … 1.0 1.0; 0.0 0.1 … 0.9 1.0],
Array{Int64,1}[[1, 2, 12, 13, 122, 123, 133, 134], [3, 4, 14, 15, 124, 125, 135, 136],
… [1063, 1064, 1074, 1075, 1184, 1185, 1195, 1196], [1065, 1066, 1076, 1077, 1186, 1187, 1197, 1198]])

julia> GL.VIEW([ GL.GLFrame2, GL.GLPol(mod3D..., GL.COLORS[1],1) ])
LinearAlgebraicRepresentation.ballFunction
ball(radius=1, angle1=pi, angle2=2*pi)(shape=[18, 36,4])

Generate a cell decomposition of a solid 3-sphere in $R^3$. The variable shape provides the domain decomposition. Empty cells are removed after the Cartesian -> Polar coordinate mapping.

Example

julia> GL.VIEW([
	GL.GLPol( Lar.ball()()..., GL.COLORS[1],0.5 ),
	GL.GLFrame ]);
LinearAlgebraicRepresentation.bboxMethod
bbox(vertices::Points)

The axis aligned bounding box of the provided set of n-dim vertices.

The box is returned as the couple of Points of the two opposite corners of the box.

LinearAlgebraicRepresentation.bbox_containsMethod
bbox_contains(container, contained)

Check if the axis aligned bounding box container contains contained.

Each input box must be passed as the couple of Points standing on the opposite corners of the box.

LinearAlgebraicRepresentation.biconnectedComponentMethod
biconnectedComponent(model)

Main procedure for computation of biconnected components of a graph represented by a LAR unsigned pair (Points, Cells).

Example

julia> V =  [3.0  6.0  6.0  6.0  6.0  0.0  3.0  10.0  10.0  10.0  0.0  3.0  0.0;
 			 2.0  2.0  0.0  4.0  8.0  4.0  4.0   4.0   0.0   8.0  8.0  0.0  0.0]

julia> EV = [[1, 2], [2, 3], [2, 4], [5, 10], [3, 12], [5, 11], [6, 7], [8, 9], 
			 [6, 13], [12, 13], [1, 7], [1, 12], [4, 7], [4, 8], [4, 5], [8, 10]]
			 
julia> VV = [[k] for k=1:size(V,2)]

julia> Plasm.view( Plasm.numbering(3)((V,[VV, EV])) )
model = (V,EV)
V,EVs = biconnectedComponent(model)
EW = convert(Lar.Cells, cat(EVs))
Plasm.view( Plasm.numbering(3)((V,[VV,EW])) )

Another example:

V = [0. 1  0  0 -1;
	 0  0  1 -1  0]
EV = [[1,2],[2,3],[1,4],[4,5],[5,1],[1,3]  ,[2,4]]
model = (V,EV)
Plasm.view( Plasm.numbering(1)((V,[[[k] for k=1:size(V,2)],EV])) )
V,EVs = Lar.biconnectedComponent(model)

cscEV = Lar.characteristicMatrix(EV)
biconcomp = Lar.Arrangement.biconnected_components(cscEV)

Matrix(Lar.characteristicMatrix(EV))

V,EVs = biconnectedComponent(model)
EW = convert(Lar.Cells, cat(EVs))
Plasm.view( Plasm.numbering(.3)((V,[VV,EW])) )
LinearAlgebraicRepresentation.boolopsMethod
boolops( assembly::Lar.Struct, op::Symbol )

User interface to 2d Boolean ops, where op symbol $in$ {:|, :&, :-}. Return a pair (V,EV).

Example

assembly = Lar.Struct([wall, openings])
V,EV = boolops(assembly, :-)
LinearAlgebraicRepresentation.boundary_1Method
boundary_1( EV::Cells )::ChainOp

Computation of sparse signed boundary operator $C_1 -> C_0$.

Example

julia> V,(VV,EV,FV,CV) = cuboid([1.,1.,1.], true);

julia> EV
12-element Array{Array{Int64,1},1}:
[1, 2]
[3, 4]
...
[2, 6]
[3, 7]
[4, 8]

julia> boundary_1( EV::Cells )
8×12 SparseMatrixCSC{Int8,Int64} with 24 stored entries:
[1 ,  1]  =  -1
[2 ,  1]  =  1
[3 ,  2]  =  -1
...       ...
[7 , 11]  =  1
[4 , 12]  =  -1
[8 , 12]  =  1

julia> Matrix(boundary_1(EV::Cells))
8×12 Array{Int8,2}:
-1   0   0   0  -1   0   0   0  -1   0   0   0
1   0   0   0   0  -1   0   0   0  -1   0   0
0  -1   0   0   1   0   0   0   0   0  -1   0
0   1   0   0   0   1   0   0   0   0   0  -1
0   0  -1   0   0   0  -1   0   1   0   0   0
0   0   1   0   0   0   0  -1   0   1   0   0
0   0   0  -1   0   0   1   0   0   0   1   0
0   0   0   1   0   0   0   1   0   0   0   1
LinearAlgebraicRepresentation.buildFVMethod
buildFV(EV::Cells, face::Cell)

The list of vertex indices that expresses the given face.

The returned list is made of the vertex indices ordered following the traversal order to keep a coherent face orientation. The edges are need to understand the topology of the face.

In this method the input face must be expressed as a Cell(=SparseVector{Int8, Int}) and the edges as Cells.

LinearAlgebraicRepresentation.buildFVMethod
buildFV(copEV::ChainOp, face::Array{Int, 1})

The list of vertex indices that expresses the given face.

The returned list is made of the vertex indices ordered following the traversal order to keep a coherent face orientation. The edges are need to understand the topology of the face.

In this method the input face must be expressed as a list of vertex indices and the edges as ChainOp.

LinearAlgebraicRepresentation.buildFVMethod
buildFV(copEV::ChainOp, face::Cell)

The list of vertex indices that expresses the given face.

The returned list is made of the vertex indices ordered following the traversal order to keep a coherent face orientation. The edges are need to understand the topology of the face.

In this method the input face must be expressed as a Cell(=SparseVector{Int8, Int}) and the edges as ChainOp.

LinearAlgebraicRepresentation.cartMethod
cart(args::Array{Array{Any,1},1})::Array{Tuple,1}

Cartesian product of collections given in the unary Array argument. Return an Array of Tuple. The number $n$ of output Tuple is equal to the product of sizes of input args

# Example

julia> cart([[1,2,3],["a","b"],[11,12]])
# output
12-element Array{Tuple{Any,Any,Any},1}:
 (1, "a", 11)
 (1, "a", 12)
 (1, "b", 11)
 (1, "b", 12)
 (2, "a", 11)
 (2, "a", 12)
 (2, "b", 11)
 (2, "b", 12)
 (3, "a", 11)
 (3, "a", 12)
 (3, "b", 11)
 (3, "b", 12)
LinearAlgebraicRepresentation.centroidMethod
centroid(P::Lar.LAR)::Array{Float64,1}

Barycenter or centroid of polyhedron P.

Example # unit 3D tetrahedron

julia> V = [0.0 1.0 0.0 0.0; 0.0 0.0 1.0 0.0; 0.0 0.0 0.0 1.0];

julia> FV = [[1, 2, 4], [1, 3, 2], [4, 3, 1], [2, 3, 4]];

julia> P = V,FV;

julia> Lar.centroid(P)
3-element Array{Float64,1}:
 0.25
 0.25
 0.25
LinearAlgebraicRepresentation.chaincomplexMethod
chaincomplex( W::Points, FW::Cells, EW::Cells )
	::Tuple{ Array{Cells,1}, Array{ChainOp,1} }

Chain 3-complex construction from bases of 2- and 1-cells.

From the minimal input, construct the whole two-dimensional chain complex, i.e. the bases for linear spaces C1 and C2 of 1-chains and 2-chains, and the signed coboundary operators from C0 to C1 and from C1 to C2.

Example

julia> cube_1 = ([0 0 0 0 1 1 1 1; 0 0 1 1 0 0 1 1; 0 1 0 1 0 1 0 1],
[[1,2,3,4],[5,6,7,8],[1,2,5,6],[3,4,7,8],[1,3,5,7],[2,4,6,8]],
[[1,2],[3,4],[5,6],[7,8],[1,3],[2,4],[5,7],[6,8],[1,5],[2,6],[3,7],[4,8]] )

julia> cube_2 = Struct([t(0,0,0.5), r(0,0,pi/3), cube_1])

julia> V,FV,EV = struct2lar(Struct([ cube_1, cube_2 ]))

julia> V,bases,coboundaries = chaincomplex(V,FV,EV)

julia> (EV, FV, CV), (cscEV, cscFE, cscCF) = bases,coboundaries

julia> FV # bases[2]
18-element Array{Array{Int64,1},1}:
[1, 3, 4, 6]
[2, 3, 5, 6]
[7, 8, 9, 10]
[1, 2, 3, 7, 8]
[4, 6, 9, 10, 11, 12]
[5, 6, 11, 12]
[1, 4, 7, 9]
[2, 5, 11, 13]
[2, 8, 10, 11, 13]
[2, 3, 14, 15, 16]
[11, 12, 13, 17]
[11, 12, 13, 18, 19, 20]
[2, 3, 13, 17]
[2, 13, 14, 18]
[15, 16, 19, 20]
[3, 6, 12, 15, 19]
[3, 6, 12, 17]
[14, 16, 18, 20]

julia> CV # bases[3]
3-element Array{Array{Int64,1},1}:
[2, 3, 5, 6, 11, 12, 13, 14, 15, 16, 18, 19, 20]
[2, 3, 5, 6, 11, 12, 13, 17]
[1, 2, 3, 4, 6, 7, 8, 9, 10, 11, 12, 13, 17]

julia> cscEV # coboundaries[1]
34×20 SparseMatrixCSC{Int8,Int64} with 68 stored entries: ...

julia> cscFE # coboundaries[2]
18×34 SparseMatrixCSC{Int8,Int64} with 80 stored entries: ...

julia> cscCF # coboundaries[3]
4×18 SparseMatrixCSC{Int8,Int64} with 36 stored entries: ...
LinearAlgebraicRepresentation.chaincomplexMethod
chaincomplex( W::Points, EW::Cells )::Tuple{Array{Cells,1},Array{ChainOp,1}}

Chain 2-complex construction from basis of 1-cells.

From the minimal input, construct the whole two-dimensional chain complex, i.e. the bases for linear spaces C1 and C2 of 1-chains and 2-chains, and the signed coboundary operators from C0 to C1 and from C1 to C2.

Example

julia> W =
[0.0  0.0  0.0  0.0  1.0  1.0  1.0  1.0  2.0  2.0  2.0  2.0  3.0  3.0  3.0  3.0
0.0  1.0  2.0  3.0  0.0  1.0  2.0  3.0  0.0  1.0  2.0  3.0  0.0  1.0  2.0  3.0]
# output
2×16 Array{Float64,2}: ...

julia> EW =
[[1, 2],[2, 3],[3, 4],[5, 6],[6, 7],[7, 8],[9, 10],[10, 11],[11, 12],[13, 14],
[14, 15],[15, 16],[1, 5],[2, 6],[3, 7],[4, 8],[5, 9],[6, 10],[7, 11],[8, 12],
[9, 13],[10, 14],[11, 15],[12, 16]]
# output
24-element Array{Array{Int64,1},1}: ...

julia> V,bases,coboundaries = chaincomplex(W,EW)

julia> bases[1]	# edges
24-element Array{Array{Int64,1},1}: ...

julia> bases[2] # faces -- previously unknown !!
9-element Array{Array{Int64,1},1}: ...

julia> coboundaries[1] # coboundary_1
24×16 SparseMatrixCSC{Int8,Int64} with 48 stored entries: ...

julia> Matrix(coboundaries[2]) # coboundary_1: faces as oriented 1-cycles of edges
9×24 Array{Int8,2}:
-1  0  0  1  0  0  0  0  0  0  0  0  1 -1  0  0  0  0  0  0  0  0  0  0
0 -1  0  0  1  0  0  0  0  0  0  0  0  1 -1  0  0  0  0  0  0  0  0  0
0  0 -1  0  0  1  0  0  0  0  0  0  0  0  1 -1  0  0  0  0  0  0  0  0
0  0  0 -1  0  0  1  0  0  0  0  0  0  0  0  0  1 -1  0  0  0  0  0  0
0  0  0  0 -1  0  0  1  0  0  0  0  0  0  0  0  0  1 -1  0  0  0  0  0
0  0  0  0  0 -1  0  0  1  0  0  0  0  0  0  0  0  0  1 -1  0  0  0  0
0  0  0  0  0  0  0 -1  0  0  1  0  0  0  0  0  0  0  0  0  0  1 -1  0
0  0  0  0  0  0 -1  0  0  1  0  0  0  0  0  0  0  0  0  0  1 -1  0  0
0  0  0  0  0  0  0  0 -1  0  0  1  0  0  0  0  0  0  0  0  0  0  1 -1
LinearAlgebraicRepresentation.characteristicMatrixMethod
characteristicMatrix( FV::Cells )::ChainOp

Binary matrix representing by rows the p-cells of a cellular complex. The input parameter must be of Cells type. Return a sparse binary matrix, providing the basis of a $Chain$ space of given dimension. Notice that the number of columns is equal to the number of vertices (0-cells).

Example

V,(VV,EV,FV,CV) = cuboid([1.,1.,1.], true);

julia> Matrix(characteristicMatrix(FV))
6×8 Array{Int8,2}:
1  1  1  1  0  0  0  0
0  0  0  0  1  1  1  1
1  1  0  0  1  1  0  0
0  0  1  1  0  0  1  1
1  0  1  0  1  0  1  0
0  1  0  1  0  1  0  1

julia> Matrix(characteristicMatrix(CV))
1×8 Array{Int8,2}:
1  1  1  1  1  1  1  1

julia> Matrix(characteristicMatrix(EV))
12×8 Array{Int8,2}:
1  1  0  0  0  0  0  0
0  0  1  1  0  0  0  0
0  0  0  0  1  1  0  0
0  0  0  0  0  0  1  1
1  0  1  0  0  0  0  0
0  1  0  1  0  0  0  0
0  0  0  0  1  0  1  0
0  0  0  0  0  1  0  1
1  0  0  0  1  0  0  0
0  1  0  0  0  1  0  0
0  0  1  0  0  0  1  0
0  0  0  1  0  0  0  1
LinearAlgebraicRepresentation.circleFunction
circle(radius=1.; angle=2*pi)(shape=36)

Compute an approximation of the circunference curve in 2D, centered on the origin.

With default values, i.e. circle()(), return the whole circonference of unit radius, approximated with a shape=36 number of 1-cells.

Example

julia> W,CW = Lar.circle()();

julia> GL.VIEW([
	GL.GLLines(W, CW, GL.COLORS[12]),
	GL.GLFrame
]);
LinearAlgebraicRepresentation.coboundary_1Function
coboundary_1( FV::Cells, EV::Cells)::ChainOp

Generate the signed sparse matrix of the coboundary_1 operator. For each row, start with the first incidence number positive (i.e. assign the orientation of the first edge to the 1-cycle of the face), then bounce back and forth between vertex columns/rows of EV and FE.

Example

julia> Matrix(cscFE) 5×20 Array{Int8,2}: 1 1 1 1 1 0 0 0 1 0 0 0 0 0 0 0 1 0 0 1 1 0 1 0 1 1 1 1 0 1 1 1 0 1 1 1 0 1 1 0 0 0 0 1 0 0 0 1 1 0 0 0 1 0 1 1 1 1 1 1 0 1 0 0 0 1 1 0 0 1 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 1 1 1 1 0 0 0 0 0 0

LinearAlgebraicRepresentation.compute_FVMethod
compute_FV( copEV::Lar.ChainOp, copFE::Lar.ChainOp )::Lar.Cells

Compute the FV array of type Lar.Cells from two Lar.ChainOp, via sparse array product. To be generalized to open 2-manifolds.

LinearAlgebraicRepresentation.cop2larMethod
cop2lar(cop::Lar.ChainOp)::Lar.Cells

Convert a sparse array of type ChainOp into an array of array of type Cells.

Notice that cop2lar is the inverse function of lar2cop. their composition is the identity function.

Example

julia> V,(VV,EV,FV,CV) = Lar.cuboid([1,1,1],true);

julia> Lar.cop2lar(Lar.lar2cop(EV))
12-element Array{Array{Int64,1},1}:
 [1, 2]
 [3, 4]
   ...
 [2, 6]
 [3, 7]
 [4, 8]

julia> Lar.cop2lar(Lar.lar2cop(FV))
6-element Array{Array{Int64,1},1}:
 [1, 2, 3, 4]
 [5, 6, 7, 8]
 [1, 2, 5, 6]
 [3, 4, 7, 8]
 [1, 3, 5, 7]
 [2, 4, 6, 8]

julia> Lar.cop2lar(Lar.lar2cop(CV))
1-element Array{Array{Int64,1},1}:
 [1, 2, 3, 4, 5, 6, 7, 8]
LinearAlgebraicRepresentation.crossingTestMethod

Half-line crossing test. Utility function for pointInPolygonClassification function. Update the count depending of the actual crossing of the tile half-line.

LinearAlgebraicRepresentation.crownFunction
crown(r=1., R=2., angle=2*pi)(shape=[24, 36])

Compute a cellular 2-complex, approximation of a two-dimensional open surface, embedded in a three-dimensional Euclidean space. This open surface is generated as an "half-torus", providing only the external shell.

Example

julia> GL.VIEW([
	GL.GLGrid( Lar.crown()()..., GL.COLORS[1],0.75 ),
	GL.GLFrame
]);
LinearAlgebraicRepresentation.cubicbezier2DMethod
cubicbezier2D(curvePts::Array{Array{Float,1},1})

Produce the two $coordinate functions$ for a $cubic$Bézier curve in $2D$. The input is the array curvePts of 4 control points with 2 coordinates.

Example

julia> curvePts = [[1232.24, 1340.84],[1259.53, 1119.86],
	[1417.91, 1015.16],[1133.47, 1010.63]]

julia> bx,by = cubicbezier2D(curvePts);

julia> (bx(0), by(0))
(1232.24, 1340.84)

julia> (bx(1), by(1))
(1133.47, 1010.63)
LinearAlgebraicRepresentation.cuboidFunction
cuboid(maxpoint::Array, full=false, minpoint::Array=zeros(length(maxpoint)))

Return a $d$-dimensional cube, where $d$ is the common length of arrays minpoint and maxpoint. If flag=true the cells of all dimensions (between 1 and $d$) are generated.

julia> cuboid([-0.5, -0.5])
([0.0 0.0 -0.5 -0.5; 0.0 -0.5 0.0 -0.5], Array{Int64,1}[[1, 2, 3, 4]])

julia> cuboid([-0.5, -0.5, 0], true)
([0.0 0.0 … -0.5 -0.5; 0.0 0.0 … -0.5 -0.5; 0.0 0.0 … 0.0 0.0],
Array{Array{Int64,1},1}[Array{Int64,1}[[1], [2], [3], [4], [5], [6], [7], [8]],
Array{Int64,1}[[1, 2], [3, 4], [5, 6], [7, 8], [1, 3], [2, 4], [5, 7], [6, 8], [1, 5], [2,
6], [3, 7], [4, 8]], Array{Int64,1}[[1, 2, 3, 4], [5, 6, 7, 8], [1, 2, 5, 6], [3, 4, 7,
8], [1, 3, 5, 7], [2, 4, 6, 8]], Array{Int64,1}[[1, 2, 3, 4, 5, 6, 7, 8]]])

julia> V, (VV, EV, FV, CV) = Lar.cuboid([1,1,1], true);

julia> assembly = Lar.Struct([ (V, CV), Lar.t(1.5,0,0), (V, CV) ])

julia> GL.VIEW([
	GL.GLPol( Lar.struct2lar(assembly)..., GL.COLORS[1],0.75 ),
	GL.GLFrame ]);
LinearAlgebraicRepresentation.cuboidGridFunction
cuboidGrid( shape, filled=false )::Union( Cells, Array{Cells,1} )

Multi-dimensional generator function. Generate either a solid $d$-grid of unit $d$-cuboids in $d$-dimensional space, or the array of $p$-skeletons ($0 <=p<= d$), depending on the Boolean variable filled. $0$-cuboids are points, $1$-cuboids are segments, , $2$-cuboids are squares, $3$-cuboids are cubes, etc. The shape=[a,b,c] value determines the number $a x b x c$ of $d$-cells. Notice that d = length(shape)

LinearAlgebraicRepresentation.cylinderFunction
cylinder(radius=.5, height=2., angle=2*pi)(shape=[36, 1])

Compute a cellular 2-complex, approximation of a right circular cylindrical surface in 3D. The open surface has basis on $z=0$ plane and is centered around the $z$ axis.

Example

julia> GL.VIEW([
	GL.GLGrid( Lar.cylinder()()..., GL.COLORS[1],1 ),
	GL.GLFrame
]);
LinearAlgebraicRepresentation.delete_edgesMethod
delete_edges(todel, V::Points, EV::ChainOp)

Delete edges and remove unused vertices from a 2-skeleton.

Loop over the todel edge index list and remove the marked edges from EV. The vertices in V which remained unconnected after the edge deletion are deleted too.

LinearAlgebraicRepresentation.depth_first_searchFunction
depth_first_search(EV::Lar.Cells [, start::Int=1])::Tuple{Lar.Cells, Lar.Cells}

Depth-First-Search algorithm (Tarjan, 1972) to compute a spanning tree of an undirected graph. The graph is entered as an array of edges, given as an array EV of 1-cells. Complexity linear in the number of nodes and edges: O(V, E).

Return a partition of EV into two arrays spanningtree and fronds which define a $palm tree$. The algorithm may start from any vertex index.

Example

julia> EV = [[1,2],[2,3],[3,4],[4,5],[5,6],[1,6],
       [6,7],[4,7],[2,7],[5,8],[1,8],[3,8]];

julia> depth_first_search(EV,8)
(Array{Int64,1}[[8, 1], [1, 2], [2, 3], [3, 4], [4, 5], [5, 6], [6, 7]], Array{Int64,1}[[6, 1], [7, 2], [7, 4], [5, 8], [3, 8]])

julia> depth_first_search(EV)
(Array{Int64,1}[[1, 2], [2, 3], [3, 4], [4, 5], [5, 6], [6, 7], [5, 8]], Array{Int64,1}[[6, 1], [7, 2], [7, 4], [8, 1], [8, 3]])
LinearAlgebraicRepresentation.diskFunction
disk(radius=1., angle=2*pi)(shape=[36, 1])

Compute the cellular complex approximating a circular sector of 2D disk centered on the origin. In geometry, a disk is the region in a plane bounded by a circle. The shape array provides the number of approximating 2-cells.

Example

julia> GL.VIEW([
	GL.GLGrid( Lar.disk()()..., GL.COLORS[1],1 ),
	GL.GLFrame
]);
LinearAlgebraicRepresentation.edge_biconnectMethod
edge_biconnect(EV::Lar.Cells)

Compute maximal 2-edge-connected components of an undirected graph.

Examples

V = hcat([[1.,2],[1,3],[0,3],[1,4],[2,3],[2,2],[1,1],[1,0],[0,1]]...)
EV = [[3,4],[2,4],[2,3],[2,5],[1,2],[5,6],[1,6],[1,7],[7,9],[8,9],[7,8]]
using Plasm
Plasm.view( Plasm.numbering(1.3)((V,[[[k] for k=1:size(V,2)], EV])) )

julia> edge_biconnect(EV::Lar.Cells)
4-element Array{Array{Array{Int64,1},1},1}:
 [[4, 2], [3, 4], [2, 3]]        
 [[6, 1], [5, 6], [2, 5], [1, 2]]
 [[9, 7], [8, 9], [7, 8]]        
 [[1, 7]]                        
julia> V = hcat([[0,4],[4,4],[4,0],[2,0],[0,0],[2,2],[1,2],[3,2],[1,3],[3,3]]...)
julia> EV = [[1,2],[9,10],[1,5],[7,9],[8,10],[2,3],[6,7],[6,8],[4,6],[4,5],[3,4]]
julia> Plasm.view( Plasm.numbering(1.3)((V,[[[k] for k=1:size(V,2)], EV])) )

julia> EVs = edge_biconnect(EV::Lar.Cells)
3-element Array{Array{Array{Int64,1},1},1}:
 [[8, 6], [10, 8], [9, 10], [7, 9], [6, 7]]
 [[4, 6]]                                  
 [[5, 1], [4, 5], [3, 4], [2, 3], [1, 2]]  

julia> map(Lar.depth_first_search, EVs)
LinearAlgebraicRepresentation.extrudeSimplicialMethod
extrudeSimplicial(model::LAR, pattern::Array)::LAR

Algorithm for multimensional extrusion of a simplicial complex. Can be applied to 0-, 1-, 2-, ... simplicial models, to get a 1-, 2-, 3-, .... model. The pattern Array is used to specify how to decompose the added dimension.

A model is a LAR model, i.e. a pair (vertices,cells) to be extruded, whereas pattern is an array of Int64, to be used as lateral measures of the extruded model. pattern elements are assumed as either solid or empty measures, according to their (+/-) sign.

Example

julia> V = [[0,0] [1,0] [2,0] [0,1] [1,1] [2,1] [0,2] [1,2] [2,2]];

julia> FV = [[1,2,4],[2,3,5],[3,5,6],[4,5,7],[5,7,8],[6,8,9]];

julia> pattern = repeat([1,2,-3],outer=4);

julia> model = (V,FV);

julia> W,FW = extrudeSimplicial(model, pattern);

julia> Plasm.view(W,FW)
LinearAlgebraicRepresentation.faces2polygonsMethod
faces2polygons(copEV::ChainOp, copFE::ChainOp)::Array{Array{Array{Int64,1},1},1}

Generate an array of faces (array of polygons, given as oriented cycles of vertex indices).

The input is a chain of two $chain operators$ (coboundaries). The output is a list of list of polygons, given as oriented $cycles of vertex$ indices. $Outer$ polygons are counterclockwise-oriented; $inner$ polygons are clockwise-oriented. The first polygon of each face is the outer one.

# Example Both a polygon with $holes$ and polygons $fitting the holes$ are defined in the following. The reader should note the $mutual orientation$ of holes and polygons fitting them. They are correctly derived from a sparse copFE matrix created by $TGW algorithm$ in 2D.

julia> V = [540.313 540.313 2038.65 2038.65 1990.25 1990.25 583.951 583.951 2038.65 2038.65 1551.97 1990.25 1551.97 1990.25 409.035 2265.24 409.035 2265.24 540.313 540.313 583.951 583.951 934.346 1246.02 934.346 1246.02; -2129.59 -1653.96 -2129.59 -1493.96 -2064.15 -1493.96 -2064.15 -1653.96 -1104.9 -1221.14 -1104.9 -1137.88 -1137.88 -1221.14 -1007.31 -1007.31 -2210.65 -2210.65 -1109.81 -1360.09 -1104.1 -1360.09 -1104.9 -1104.9 -1137.88 -1137.88]

julia> copEV = sparse([1, 2, 1, 8, 2, 3, 3, 7, 4, 5, 4, 7, 5, 6, 6, 8, 9, 10, 9, 13, 10, 14, 11, 12, 11, 14, 12, 13, 15, 16, 15, 17, 16, 18, 17, 18, 19, 21, 19, 22, 20, 21, 20, 22, 23, 24, 23, 26, 24, 25, 25, 26], [1, 1, 2, 2, 3, 3, 4, 4, 5, 5, 6, 6, 7, 7, 8, 8, 9, 9, 10, 10, 11, 11, 12, 12, 13, 13, 14, 14, 15, 15, 16, 16, 17, 17, 18, 18, 19, 19, 20, 20, 21, 21, 22, 22, 23, 23, 24, 24, 25, 25, 26, 26], Int8[-1, -1, 1, -1, 1, -1, 1, -1, -1, -1, 1, 1, 1, -1, 1, 1, -1, -1, 1, -1, 1, -1, -1, -1, 1, 1, 1, 1, -1, -1, 1, -1, 1, -1, 1, 1, -1, -1, 1, -1, -1, 1, 1, 1, -1, -1, 1, -1, 1, -1, 1, 1])

julia> copFE = sparse([1, 3, 1, 3, 1, 3, 1, 3, 1, 3, 1, 3, 1, 3, 1, 3, 2, 3, 2, 3, 2, 3, 2, 3, 2, 3, 2, 3, 3, 3, 3, 3, 3, 4, 3, 4, 3, 4, 3, 4, 3, 5, 3, 5, 3, 5, 3, 5], [1, 1, 2, 2, 3, 3, 4, 4, 5, 5, 6, 6, 7, 7, 8, 8, 9, 9, 10, 10, 11, 11, 12, 12, 13, 13, 14, 14, 15, 16, 17, 18, 19, 19, 20, 20, 21, 21, 22, 22, 23, 23, 24, 24, 25, 25, 26, 26], Int8[-1, 1, 1, -1, 1, -1, -1, 1, 1, -1, 1, -1, 1, -1, -1, 1, -1, 1, 1, -1, -1, 1, 1, -1, -1, 1, 1, -1, -1, 1, -1, 1, -1, 1, 1, -1, 1, -1, -1, 1, 1, -1, -1, 1, -1, 1, 1, -1])

julia> faces2polygons(copEV,copFE)
5-element Array{Array{Array{Int64,1},1},1}:
 [[1, 3, 4, 6, 5, 7, 8, 2]]
 [[9, 11, 13, 12, 14, 10]]
 [[2, 8, 7, 5, 6, 4, 3, 1], [10, 14, 12, 13, 11, 9], [15, 17, 18, 16], [19, 21, 22, 20], [24, 26, 25, 23]]
 [[20, 22, 21, 19]]
 [[23, 25, 26, 24]]

julia> model = (V, [ [[v] for v=1:size(V,2)], cop2lar(copEV) ])
([540.313 540.313 … 934.346 1246.02; -2129.59 -1653.96 … -1137.88 -1137.88], Array{Array{Int64,1},1}[[[1], [2], [3], [4], [5], [6], [7], [8], [9], [10]  …  [17], [18], [19], [20], [21], [22], [23], [24], [25], [26]], [[1, 2], [1, 3], [3, 4], [5, 6], [5, 7], [7, 8], [4, 6], [2, 8], [9, 10], [9, 11]  …  [16, 18], [17, 18], [19, 20], [21, 22], [19, 21], [20, 22], [23, 24], [23, 25], [25, 26], [24, 26]]])

julia> GL.VIEW( GL.numbering(200)( model, GL.COLORS[1], 0.1 ) );
LinearAlgebraicRepresentation.filterByOrderMethod
filterByOrder( n::Int )Array{Array{Array{Int8,1},1},1}

Filter the array of codes (BooleanString) of $n$ bits depending on their integer value (order).

Example

julia> filterByOrder(3)
# output
4-element Array{Array{Array{Int8,1},1},1}:
 Array{Int8,1}[Int8[0, 0, 0]]
 Array{Int8,1}[Int8[0, 0, 1], Int8[0, 1, 0], Int8[1, 0, 0]]
 Array{Int8,1}[Int8[0, 1, 1], Int8[1, 0, 1], Int8[1, 1, 0]]
 Array{Int8,1}[Int8[1, 1, 1]]
LinearAlgebraicRepresentation.firstMomentMethod
firstMoment(P::Lar.LAR)::Array{Float64,1}

First moments as terms of the Euler tensor. Remember that the integration algorithm is a boundary integration. Hence the model must be a boundary model. In this case, a 2-complex of triangles.

Example # unit 3D tetrahedron

julia> V = [0.0 1.0 0.0 0.0; 0.0 0.0 1.0 0.0; 0.0 0.0 0.0 1.0];

julia> FV = [[1, 2, 4], [1, 3, 2], [4, 3, 1], [2, 3, 4]];

julia> P = V,FV;

julia> Lar.firstMoment(P)
3-element Array{Float64,1}:
 0.0416667
 0.0416667
 0.0416667
LinearAlgebraicRepresentation.fix_redundancyMethod
fix_redundancy(target_mat, cscFV,cscEV)

Fix the coboundary_1 matrix, generated by sparse matrix product $FV * EV^t$, for complexes with some non-convex cells. This approach can be used when both EV and FV of the cellular complex are known. It is exact when cells are convex. Maybe non-exact, introducing spurious incidence coefficients ($redundancies$), when adjacent faces share an edge combinatorially, but not geometrically. This happen when an edge is on the boundary of face A, but only its vertices are on the boundary of face B. TODO: Similar situations may appear when computing algebraically CF as product of known CV and FV, with non-convex cells.

In order to remove such $redundancies$, the Euler characteristic of 2-sphere is used, where V-E+F=2. Since we have F=2 (inner and outer face), $V=E$ must hold, and d=E-V is the (non-negative) $defect$ number, called nfixs in the code. It equates the number of columns edges whose sum is greater than 2 for the considered row (face). Remember the in a $d$-complex, including the $outer cell$, all $(d-1)$-faces must be shared by exactly 2 $d$-faces. Note that FVmust include the row of outer shell (exterior face).

Example

FV = [[1,2,3,4,5,17,16,12],
[1,2,3,4,6,7,8,9,10,11,12,13,14,15],
[4,5,9,11,12,13,14,15,16,17],
[2,3,6,7], [8,9,10,11]]

FE = [[1,2,3,4,9,20,17,5],
[1,6,10,7,3,8,11,12,14,15,19,18,16,5],
[4,9,20,17,16,18,19,15,13,8],
[2,10,6,7], [11,12,13,14]]

EV = [[1,2],[2,3],[3,4],[4,5],[1,12],[2,6],[3,7],[4,9],[5,17],[6,7],[8,9],
[8,10],[9,11],[10,11],[11,15],[12,13],[12,16],[13,14],[14,15],[16,17]]

V = [0   2   5   7  10   2   5   3   7  3  7  0  3  3  7  0  10;
    16  16  16  16  16  13  13  11  11  8  8  5  5  2  2  0   0]

cscFE = u_coboundary_1( FV::Cells, EV::Cells, false);
Matrix(cscFE)

Notice that there are two columns (2 and 13) with 3 ones, hence (3-2)+(3-2)=2 defects to fix. The fixed complex can be shown graphically as:


VV = [[k] for k in 1:size(V,2)];
using Plasm
Plasm.view( Plasm.numbering(3)((V,[VV, EV, FV])) )
LinearAlgebraicRepresentation.fragmentlinesMethod
fragmentlines(model::Lar.LAR)::Lar.LAR

Pairwise intersection of 2D line segments.

Example 2D

V,EV = model2d
W, EW = Lar.fragmentlines(model2d) # OK
using Plasm
Plasm.viewexploded(W,EW)(1.2,1.2,1.2)
LinearAlgebraicRepresentation.gridMethod
grid(sequence::Array{Number,1})::Lar.LAR

Generate a 1D LAR model.

Geometry is stored as 1D Points, and Topology is stored as 1D Cells. q() and q()() are used as alias function name.

julia> model1D = Lar.grid(1,-1,1,-1,1,-1,1,-1,1,-1)
([0.0 1.0 … 9.0 10.0], Array{Int64,1}[[1, 2], [3, 4], [5, 6], [7, 8], [9, 10]])

julia> model1D[1]
1×11 Array{Float64,2}:
 0.0  1.0  2.0  3.0  4.0  5.0  6.0  7.0  8.0  9.0  10.0

julia> model1D[2]
5-element Array{Array{Int64,1},1}:
 [1, 2]
 [3, 4]
 [5, 6]
 [7, 8]
 [9, 10]

 julia> mod = Lar.grid(repeat([.1,-.1],outer=5)...)
 ([0.0 0.1 … 0.9 1.0], Array{Int64,1}[[1, 2], [3, 4], [5, 6], [7, 8], [9, 10]])
LinearAlgebraicRepresentation.grid_0Method
grid_0(n::Int)::Array{Int64,2}

Generate a uniform 0D cellular complex. The grid_0 function generates a 0-dimensional uniform complex embedding $n+1$ equally-spaced 0-cells (at unit interval boundaries). It returns by columns the cells of this 0-complex as `Array{Int64,2}.

# Example

julia> grid_0(10)
# output
1×11 Array{Int64,2}:
 0  1  2  3  4  5  6  7  8  9  10
LinearAlgebraicRepresentation.grid_1Method
grid_1(n::Int)::Array{Int64,2}

Generate a uniform 1D cellular complex. The grid_1 function generates a 0-dimensional uniform complex embedding $n+1$ equally-spaced 1-cells (unit intervals). It returns by columns the cells of this 1-complex as Array{Int64,2}.

# Example

julia> grid_1(10)
# output
2×10 Array{Int64,2}:
 0  1  2  3  4  5  6  7  8   9
 1  2  3  4  5  6  7  8  9  10
LinearAlgebraicRepresentation.helicoidFunction
helicoid(R=1., r=0.5, pitch=1., nturns=2)(shape=[36*nturns, 2])

Compute an approximation of the helicoid surface in 3D, with basis on $z=0$ plane and centered around the $z$ axis.

Example

julia> GL.VIEW([
	GL.GLGrid( Lar.helicoid()()..., GL.COLORS[1],1 ),
	GL.GLFrame
]);
LinearAlgebraicRepresentation.helixFunction
helix(radius=1., pitch=1., nturns=2)(shape=36*nturns)

Compute the approximate elix curve in three-dimensional space, with basis on $z=0$ plane and centered around the $z$ axis. The pitch of a helix is the height of one complete helix turn, measured parallel to the axis of the helix.

Example

julia> V, CV = Lar.helix(.1, .1, 10)()
([0.1 0.0984808 … 0.0984808 0.1; 0.0 0.0173648 … -0.0173648 0.0; 0.0 0.0027778 … 0.997222 1.0], Array{Int64,1}[[1, 2], [2, 3], [3, 4], [4, 5], [5, 6], [6, 7], [7, 8], [8, 9], [9, 10], [10, 11]  …  [351, 352], [352, 353], [353, 354], [354, 355], [355, 356], [356, 357], [357, 358], [358, 359], [359, 360], [360, 361]])

julia> GL.VIEW([
	GL.GLLines(V, CV, GL.COLORS[12]),
	GL.GLFrame
]);
LinearAlgebraicRepresentation.hollowBallFunction
hollowBall(r=1., R=2., angle1=pi, angle2=2*pi)(shape=[36, 1, 1])

Compute the cellular 3-complex approximating a 3-sphere. The model is meshed with cubical 3-cells, where the mesh has default decomposition size [24, 36, 8].

Example

julia> V, CV = Lar.hollowBall(1, 2, pi/2, pi/2)([6, 12, 4]);

julia> GL.VIEW([
 	GL.GLPol( V, CV, GL.COLORS[1],0.5 ),
 	GL.GLFrame ]);
...
LinearAlgebraicRepresentation.hollowCylFunction
hollowCyl(r=1., R=2., height=6., angle=2*pi)(shape=[36, 1, 1])

Compute the cellular 3-complex approximating a solid cylinder with a internal axial hole. The model is meshed with cubical 3-cells.

Example

julia> GL.VIEW([
 	GL.GLPol( Lar.hollowCyl()()..., GL.COLORS[1],0.5 ),
 	GL.GLFrame ]);
LinearAlgebraicRepresentation.index2addrMethod
index2addr(shape::Array{Int64,2})(multiIndex::Array{Int,1})::Int

Multi-index to address transformation. Partial function allowing for using both horizontal and vertical vectors for shape parameter. Notice that multi-indices are used here as coordinates of grid points, hence they start from tuples of zeros. Accordingly, the translation formula for multi-index to address transformation is 0-based.

LinearAlgebraicRepresentation.index2addrMethod
index2addr(shape::Array{Int64,1})(multiIndex)::Int

Multi-index to address transformation. Multi-index is a generalization of the concept of an integer index to an ordered tuple of indices. The second-order utility function index2addr transforms a shape list for a multidimensional array into a function that, when applied to a multindex array, i.e. to a list of integer Tuple within the shape's bounds, returns the integer addresses of the corresponding array components within the linear storage of the multidimensional array.

# Example Notice that in the example below, there are $3 x 6$ different multi-index values for the variable index, generated by cart([ 0:2, 0:5 ]).

julia> [index2addr([3,6])(collect(index)) for index in cart([ 0:2, 0:5 ])]'
# output
1×18 RowVector{Int64,Array{Int64,1}}:
 1  2  3  4  5  6  7  8  9  10  11  12  13  14  15  16  17  18

julia> index2addr([3,6])([0,0])
# output
1

julia> index2addr([3,6])([2,5])
# output
18
LinearAlgebraicRepresentation.inertiaMomentMethod
inertiaMoment(P::Lar.LAR)::Array{Float64,1}

Inertia moments of polyhedron P.

Example # unit 3D tetrahedron

julia> V = [0.0 1.0 0.0 0.0; 0.0 0.0 1.0 0.0; 0.0 0.0 0.0 1.0];

julia> FV = [[1, 2, 4], [1, 3, 2], [4, 3, 1], [2, 3, 4]];

julia> P = V,FV;

julia> Lar.inertiaMoment(P)
3-element Array{Float64,1}:
 0.0333333
 0.0333333
 0.0333333
LinearAlgebraicRepresentation.inertiaProductMethod
inertiaProduct(P::Lar.LAR)::Array{Float64,1}

Inertia products as terms of the Euler tensor.

Example # unit 3D tetrahedron

julia> V = [0.0 1.0 0.0 0.0; 0.0 0.0 1.0 0.0; 0.0 0.0 0.0 1.0];

julia> FV = [[1, 2, 4], [1, 3, 2], [4, 3, 1], [2, 3, 4]];

julia> P = V,FV;

julia> Lar.inertiaProduct(P)
3-element Array{Float64,1}:
 0.00833333
 0.00833333
 0.00833333
LinearAlgebraicRepresentation.input_collectionMethod
input_collection(data::Array)::Tuple

Facet selection. Construction of a $(d-1)$-dimensional collection from a $(d-1)$- or $d$-dimensional one. $0-chain$ of LAR type are used as input.

Output is $admissible input$ for algorithms of the 2D/3D arrangement pipeline.

Example 2D

An assembly of geometric objects is generated, and their assembly, including rotated and translated chains, is built producing a collection of input LAR models.

V,(_,EV,FV) = Lar.cuboidGrid([4,4],true);
W,(_,EW,FW) = Lar.cuboidGrid([3,5],true);
mycircle(r,n) = Lar.circle(r)(n)

data2d1 = (V,EV)
data2d2 = Lar.Struct([ Lar.t(2,2), Lar.r(pi/3), Lar.t(-1.5,-2.5), (W,EW) ])
data2d3 = Lar.Struct([ Lar.t(2,2), mycircle(2.5,16) ])
data2d4 = Lar.Struct([ Lar.t(3.5,3.5), mycircle(.25,16) ])
data2d5 = Lar.Struct([ Lar.t(5,3.5), mycircle(.5,16) ])
data2d6 = Lar.Struct([ Lar.t(5,3.5), mycircle(.25,16) ])

model2d = input_collection( [ data2d1, data2d2, data2d3, data2d4, data2d5, data2d6 ] )
V,EV = model2d
VV = [[k] for k in 1:size(V,2)];
using Plasm
Plasm.view( Plasm.numbering(.5)((V,[VV,EV])) )

Note that V,EV is not a cellular complex, since 1-cells intersect out of 0-cells.

Example 3D

V,FV = Lar.sphere(2)([3,4])
EV = Lar.simplexFacets(FV)
mysphere = V,FV,EV

data3d1 = mysphere
data3d2 = Lar.Struct([ Lar.t(0,1,0), mysphere ])
data3d3 = Lar.Struct([ Lar.t(0,0.5,0), Lar.s(0.4,0.4,0.4), mysphere ])
data3d4 = Lar.Struct([ Lar.t(4,0,0), Lar.s(0.8,0.8,0.8), mysphere ])
data3d5 = Lar.Struct([ Lar.t(4,0,0), Lar.s(0.4,0.4,0.4), mysphere ])

model3d = input_collection([ data3d1, data3d2, data3d3, data3d4, data3d5 ])
V,FV,EV = model3d
VV = [[k] for k in 1:size(V,2)];
using Plasm
Plasm.view( Plasm.numbering(1)((V,[VV, EV])) )

Note that V,FV,EV is not a cellular complex, since 1-cells and 2-cells intersect out of 0-cells.

LinearAlgebraicRepresentation.intersectionMethod
intersection(line1,line2)

Intersect two line segments in 2D, by computing the two line parameters of the intersection point.

The line segments intersect if both return parameters α,β are contained in the interval [0,1].

Example

julia> line1 = [[0.,0], [1,2]]
2-element Array{Array{Float64,1},1}:
 [0.0, 0.0]
 [1.0, 2.0]

julia> line2 = [[2.,0], [0,3]]
2-element Array{Array{Float64,1},1}:
 [2.0, 0.0]
 [0.0, 3.0]

julia> Lar.intersection(line1,line2)
(0.8571428571428571, 0.5714285714285714)
LinearAlgebraicRepresentation.lar2copMethod
lar2cop(CV::Lar.Cells)::Lar.ChainOp

Convert an array of array of integer indices to vertices into a sparse matrix.

Examples

For a single 3D unit cube we get:

julia> V,(VV,EV,FV,CV) = Lar.cuboid([1,1,1],true);

julia> Matrix(Lar.lar2cop(EV))
12×8 Array{Int8,2}:
 1  1  0  0  0  0  0  0
 0  0  1  1  0  0  0  0
 0  0  0  0  1  1  0  0
 0  0  0  0  0  0  1  1
 1  0  1  0  0  0  0  0
 0  1  0  1  0  0  0  0
 0  0  0  0  1  0  1  0
 0  0  0  0  0  1  0  1
 1  0  0  0  1  0  0  0
 0  1  0  0  0  1  0  0
 0  0  1  0  0  0  1  0
 0  0  0  1  0  0  0  1

julia> Matrix(Lar.lar2cop(FV))
6×8 Array{Int8,2}:
 1  1  1  1  0  0  0  0
 0  0  0  0  1  1  1  1
 1  1  0  0  1  1  0  0
 0  0  1  1  0  0  1  1
 1  0  1  0  1  0  1  0
 0  1  0  1  0  1  0  1

julia> Matrix(Lar.lar2cop(CV))
1×8 Array{Int8,2}:
 1  1  1  1  1  1  1  1
LinearAlgebraicRepresentation.lar2objMethod
lar2obj(V::Points, cc::ChainComplex)

Triangulated OBJ string representation of the model passed as input.

Use this function to export LAR models into OBJ

Example

	julia> cube_1 = ([0 0 0 0 1 1 1 1; 0 0 1 1 0 0 1 1; 0 1 0 1 0 1 0 1],
	[[1,2,3,4],[5,6,7,8],[1,2,5,6],[3,4,7,8],[1,3,5,7],[2,4,6,8]],
	[[1,2],[3,4],[5,6],[7,8],[1,3],[2,4],[5,7],[6,8],[1,5],[2,6],[3,7],[4,8]] )

	julia> cube_2 = Lar.Struct([Lar.t(0,0,0.5), Lar.r(0,0,pi/3), cube_1])

	julia> V, FV, EV = Lar.struct2lar(Lar.Struct([ cube_1, cube_2 ]))

	julia> V, bases, coboundaries = Lar.chaincomplex(V,FV,EV)

	julia> (EV, FV, CV), (copEV, copFE, copCF) = bases, coboundaries

	julia> FV # bases[2]
	18-element Array{Array{Int64,1},1}:
	 [1, 3, 4, 6]
	 [2, 3, 5, 6]
	 [7, 8, 9, 10]
	 [1, 2, 3, 7, 8]
	 [4, 6, 9, 10, 11, 12]
	 [5, 6, 11, 12]
	 [1, 4, 7, 9]
	 [2, 5, 11, 13]
	 [2, 8, 10, 11, 13]
	 [2, 3, 14, 15, 16]
	 [11, 12, 13, 17]
	 [11, 12, 13, 18, 19, 20]
	 [2, 3, 13, 17]
	 [2, 13, 14, 18]
	 [15, 16, 19, 20]
	 [3, 6, 12, 15, 19]
	 [3, 6, 12, 17]
	 [14, 16, 18, 20]

	julia> CV # bases[3]
	3-element Array{Array{Int64,1},1}:
	 [2, 3, 5, 6, 11, 12, 13, 14, 15, 16, 18, 19, 20]
	 [2, 3, 5, 6, 11, 12, 13, 17]
	 [1, 2, 3, 4, 6, 7, 8, 9, 10, 11, 12, 13, 17]

	julia> copEV # coboundaries[1]
	34×20 SparseMatrixCSC{Int8,Int64} with 68 stored entries: ...

	julia> copFE # coboundaries[2]
	18×34 SparseMatrixCSC{Int8,Int64} with 80 stored entries: ...

	julia> copCF # coboundaries[3]
	4×18 SparseMatrixCSC{Int8,Int64} with 36 stored entries: ...

	objs = Lar.lar2obj(V'::Lar.Points, [coboundaries...])

	open("./two_cubes.obj", "w") do f
    	write(f, objs)
	end
LinearAlgebraicRepresentation.larCellProdMethod
larCellProd(cellLists::Array{Cells,1})::Cells

Generation of grid cells by Cartesian product of 0/1-complexes. The output complex is generated by the product of any number of either 0- or 1-dimensional cell complexes. The product of $d$ 1-complexes generates solid $d$-cells, while the product of $n$ 0-complexes and $d-n$ 1-complexes ($n < d$) generates non-solid $(d-n)$-cells, properly embedded in $d$-space, i.e. with vertices having $d$ coordinates.

Examples

To understand the generation of cuboidal grids from products of 0- or 1-dimensional complexes, below we show a simple example of 2D grids embedded in $R^3$. In particular, v1 = [0. 1. 2. 3.] and v0 = [0. 1. 2.] are two 2-arrays of 1D vertices, c1 = [[0,1],[1,2],[2,3]] and c0 = [[0],[1],[2]] are the LAR representation of one $1$-complex and one $0$-complex, respectively. The solid 2-complex named grid2D is generated in 2D as follows:

julia> v1 = [0. 1. 2. 3.]
1×4 Array{Float64,2}:
 0.0  1.0  2.0  3.0

julia> c1 = [[0,1],[1,2],[2,3]]
3-element Array{Array{Int64,1},1}:
 [0, 1]
 [1, 2]
 [2, 3]

julia> grid2D = larVertProd([v1,v1]),larCellProd([c1,c1])
([0.0 0.0 … 3.0 3.0; 0.0 1.0 … 2.0 3.0], Array{Int64,1}[[1, 2, 5, 6], [2, 3, 6, 7], [3, 4, 7, 8], [5, 6, 9, 10], [6, 7, 10, 11], [7, 8, 11, 12], [9, 10, 13, 14], [10, 11, 14, 15], [11, 12, 15, 16]])

whereas a non-solid$2$-complex in $3D$ is generated as:

julia> v1, c1 = [0. 1. 2. 3.],[[0,1],[1,2],[2,3]]
([0.0 1.0 2.0 3.0], Array{Int64,1}[[0, 1], [1, 2], [2, 3]])

julia> v0, c0 = [0. 1. 2.], [[0],[1],[2]]
([0.0 1.0 2.0], Array{Int64,1}[[0], [1], [2]])

julia> vertGrid = larVertProd([v1, v1, v0])
3×48 Array{Float64,2}:
 0.0  0.0  0.0  0.0  0.0  0.0  …  3.0  3.0  3.0  3.0  3.0  3.0  3.0  3.0  3.0
 0.0  0.0  0.0  1.0  1.0  1.0  …  1.0  1.0  1.0  2.0  2.0  2.0  3.0  3.0  3.0
 0.0  1.0  2.0  0.0  1.0  2.0  …  0.0  1.0  2.0  0.0  1.0  2.0  0.0  1.0  2.0

julia> cellGrid = larCellProd([c1, c1, c0])
27-element Array{Array{Int64,1},1}:
 [1, 4, 13, 16]
 [2, 5, 14, 17]
 ...  ... ...
 [32, 35, 44, 47]
 [33, 36, 45, 48]

julia> grid3D = vertGrid,cellGrid
([0.0 0.0 … 3.0 3.0; 0.0 0.0 … 3.0 3.0; 0.0 1.0 … 1.0 2.0], Array{Int64,1}[[1, 4, 13, 16], [2, 5, 14, 17], … [32, 35, 44, 47], [33, 36, 45, 48]])

julia> using Plasm

julia> Plasm.view(grid3D)
LinearAlgebraicRepresentation.larGridMethod
larGrid(n::Int)(d::Int)::Array{Int64,2}

Generate either a uniform 0D cellular complex or a uniform 1D cellular complex. A larGrid function is given to generate the LAR representation of the cells of either a 0- or a 1-dimensional complex, depending on the value of the d parameter, to take values in the set ${0,1}$, and providing the order of the output complex.

# Example

julia> larGrid(10)(0)
# output
1×11 Array{Int64,2}:
 0  1  2  3  4  5  6  7  8  9  10

julia> larGrid(10)(1)
# output
2×10 Array{Int64,2}:
 0  1  2  3  4  5  6  7  8   9
 1  2  3  4  5  6  7  8  9  10
LinearAlgebraicRepresentation.larGridSkeletonMethod
larGridSkeleton( shape::Array{Int,1} )( d::Int )::Cells

Produce the d-dimensional skeleton (set of d-cells) of a cuboidal grid of given shape.

Example

A shape=[1,1,1] parameter refers to a grid with a single step on the three axes, i.e. to a single 3D unit cube. Below all skeletons of such simplest grid are generated.

julia> Lar.larGridSkeleton([1,1,1])(0)
# output
8-element Array{Array{Int64,1},1}:
[[1], [2], [3], [4], [5], [6], [7], [8]]

julia> Lar.larGridSkeleton([1,1,1])(1)
# output
12-element Array{Array{Int64,1},1}:
[[1,2],[3,4],[5,6],[7,8],[1,3],[2,4],[5,7],[6,8],[1,5],[2,6],[3,7],[4,8]]

julia> Lar.larGridSkeleton([1,1,1])(2)
# output
6-element Array{Array{Int64,1},1}:
[[1,2,3,4], [5,6,7,8], [1,2,5,6], [3,4,7,8], [1,3,5,7], [2,4,6,8]]

julia> Lar.larGridSkeleton([1,1,1])(3)
# output
1-element Array{Array{Int64,1},1}:
 [1, 2, 3, 4, 5, 6, 7, 8]
LinearAlgebraicRepresentation.larImageVertsMethod
larImageVerts(shape::Array{Int,1})::Array{Int64,2}

Linearize the grid of integer vertices, given the shape of a cuboidal grid (typically an image).

Examples

julia> larImageVerts([1024,1024])
# output
2×1050625 Array{Int64,2}:
 0  0  0  0  0  0  0  0  0  0   0   0 … 1024  1024  1024  1024  1024  1024  1024  1024
 0  1  2  3  4  5  6  7  8  9  10  11 … 1017  1018  1019  1020  1021  1022  1023  1024

julia> larImageVerts([1,1,1])
# output
3×8 Array{Int64,2}:
 0  0  0  0  1  1  1  1
 0  0  1  1  0  0  1  1
 0  1  0  1  0  1  0  1
LinearAlgebraicRepresentation.larModelProductMethod
larModelProduct

The larModelProduct function takes as input a pair of LAR models and returns the model of their Cartesian product. Since LAR type is a pair $(geometry,topology)$, the second element of output is the topological product of the input topologies.

Example

Data preparation follows.

julia> geom_0,topol_0 = [0. 1. 2. 3. 4.],[[1],[2],[3],[4],[5]]
([0.0 1.0 … 3.0 4.0], Array{Int64,1}[[1, 2], [2, 3], [3, 4], [4, 5]])

julia> geom_1,topol_1 = [0. 1. 2.], [[1,2],[2,3]]
([0.0 1.0 2.0], Array{Int64,1}[[1, 2], [2, 3]])

julia> mod_0 = (geom_0,topol_0)
([0.0 1.0 … 3.0 4.0], Array{Int64,1}[[1, 2], [2, 3], [3, 4], [4, 5]])

julia> mod_1 = (geom_1,topol_1)
([0.0 1.0 2.0], Array{Int64,1}[[1, 2], [2, 3]])

Generation of a 2D squares model, with 8 two-dimensional cells.

julia> squares = larModelProduct(mod_1,mod_1)
([0.0 0.0 … 4.0 4.0; 0.0 1.0 … 1.0 2.0], Array{Int64,1}[[1, 2, 4, 5], [2, 3, 5, 6], [4, 5, 7, 8], [5, 6, 8, 9], [7, 8, 10, 11], [8, 9, 11, 12], [10, 11, 13, 14], [11, 12, 14, 15]])

julia> squares[1]
2×15 Array{Float64,2}:
 0.0  0.0  0.0  1.0  1.0  1.0  2.0  2.0  2.0  3.0  3.0  3.0  4.0  4.0  4.0
 0.0  1.0  2.0  0.0  1.0  2.0  0.0  1.0  2.0  0.0  1.0  2.0  0.0  1.0  2.0

julia> squares[2]
8-element Array{Array{Int64,1},1}:
[[1,2,4,5], [2,3,5,6], [4,5,7,8], [5,6,8,9], [7,8,10,11], [8,9,11,12], [10,11,13,14], [11,12,14,15]]

Generation of a 3D cubes model, with 32 three-dimensional cells.

julia> cubes = larModelProduct(squares,mod_0)
([0.0 0.0 … 4.0 4.0; 0.0 0.0 … 2.0 2.0; 0.0 1.0 … 3.0 4.0], Array{Int64,1}[[1, 2, 6, 7, 16, 17, 21, 22], [2, 3, 7, 8, 17, 18, 22, 23], [3, 4, 8, 9, 18, 19, 23, 24], [4, 5, 9, 10, 19, 20, 24, 25], … [53, 54, 58, 59, 68, 69, 73, 74], [54, 55, 59, 60, 69, 70, 74, 75]])
LinearAlgebraicRepresentation.larVertProdMethod
larVertProd(vertLists::Array{Points,1})::Points

Generate the integer coordinates of vertices (0-cells) of a multidimensional grid. Grid n-vertices are produced by the larVertProd function, via Cartesian product of vertices of $n$ 0-dimensional arguments (vertex arrays in vertLists), orderly corresponding to $x_1, x_2, ..., x_n$ coordinates in the output points $(x_1, x_2,...,x_n)$ in $R^n$.

# Example

julia> larVertProd([ larGrid(3)(0), larGrid(4)(0) ])
# output
2×20 Array{Int64,2}:
 0  0  0  0  0  1  1  1  1  1  2  2  2  2  2  3  3  3  3  3
 0  1  2  3  4  0  1  2  3  4  0  1  2  3  4  0  1  2  3  4
LinearAlgebraicRepresentation.linefragmentsMethod
linefragments(V,EV,Sigma)

Compute the sequences of ordered parameters fragmenting each input lines.

Extreme parameter values (0.0 and 1.0) are included in each output line. Sigma is the spatial index providing the subset of lines whose containment boxes intersect the box of each input line (given by EV).

julia> V = hcat([[0.,0],[1,0],[1,1],[0,1],[2,1]]...);

julia> EV = [[1,2],[2,3],[3,4],[4,1],[1,5]];

julia> Sigma = Lar.spaceindex((V,EV))
5-element Array{Array{Int64,1},1}:
 [4, 5, 2]
 [1, 3, 5]
 [4, 5, 2]
 [1, 3, 5]
 [4, 1, 3, 2]

julia> Lar.linefragments(V,EV,Sigma)
5-element Array{Any,1}:
 [0.0, 1.0]
 [0.0, 0.5, 1.0]
 [0.0, 1.0]
 [0.0, 1.0]
 [0.0, 0.5, 1.0]
LinearAlgebraicRepresentation.normalizeMethod
normalize(V::Lar.Points; flag=true::Bool)::Lar.Points

2D normalization transformation (isomorphic by defaults) of model vertices to normalized coordinates $[0,1]^2$. Used with SVG importing.

LinearAlgebraicRepresentation.obj2larMethod
obj2lar(path)

Read OBJ file at path and create a 2-skeleton as Tuple{Points, ChainComplex} from it.

This function does care about eventual internal grouping inside the OBJ file.

LinearAlgebraicRepresentation.outputCompMethod
outputComp(stack::Array, u::Int, v::Int)::Array

Internal output of biconnected components. Utility function for DFV_visit graph algorithm (Hopcrof, Tarjan (1973)).

LinearAlgebraicRepresentation.permutationOrbitsMethod
permutationOrbits(perm::OrderedDict)::Array{Array{Int64,1},1}

Compute the $cycles$ of a perm$permutation$ of the first integers (starting from 1).

The perm parameter is an $ordered dictionary$. The output is an array of arrays of integers ($orbits$).

Examples

julia> dict(List) = OrderedDict((i,x) for (i,x) in enumerate(List))
dict (generic function with 1 method)

julia> perm = dict([2, 3, 4, 5, 6, 7, 8, 1]);

julia> permutationOrbits(perm)
1-element Array{Array{Int64,1},1}:
 [2, 3, 4, 5, 6, 7, 8, 1, 2]

julia> perm = dict([3,9,8,12,10,7,2,11,6,4,1,5]);

julia> permutationOrbits(perm)
3-element Array{Array{Int64,1},1}:
 [3, 8, 11, 1, 3]
 [9, 6, 7, 2, 9]
 [12, 5, 10, 4, 12]

julia> permutationOrbits(Dict())
0-element Array{Array{Int64,1},1}

julia> permutationOrbits(Dict(1=>1))
1-element Array{Array{Int64,1},1}:
 [1, 1]
LinearAlgebraicRepresentation.pols2triaMethod
pols2tria(W::Points, copEV::ChainOp,
		copFE::ChainOp, copCF::ChainOp)

Take a chain 3-complex and return arrays of boundary simplicial complexes. Input and ouput Points (embedding geometry) are by columns. Arrays of simplicial complexes are by body, by face, and by boundary of face.

LinearAlgebraicRepresentation.qnMethod
qn(n::Int)(sequence::Array{T,1})::Lar.LAR  where T <: Real

Alias of grid function, with repetition parameter n.

julia> Lar.qn(3)([1.5,-2,0.5])
([0.0 1.5 … 11.5 12.0], Array{Int64,1}[[1, 2], [3, 4], [4, 5], [6, 7], [7, 8], [9, 10]])
LinearAlgebraicRepresentation.quads2trianglesMethod
quads2triangles(quads::Cells)::Cells

Convert an array of quads with type ::Lar.Cells into an array of triangles with the same type.

Examples

The transformation from quads to triangles works for any 2-complex, embedded in any dimensional space

2D example

V,FV = Lar.cuboidGrid([4,5])
triangles = Lar.quads2triangles(FV::Lar.Cells)::Lar.Cells
using Plasm
Plasm.view((V,[triangles]))

3D example

V,(VV,EV,FV,CV) = Lar.cuboidGrid([4,5,3],true)
triangles = Lar.quads2triangles(FV::Lar.Cells)::Lar.Cells
using Plasm
Plasm.view((V,[triangles]))
LinearAlgebraicRepresentation.rMethod
r(args...)

Return an affine transformation Matrix in homogeneous coordinates. Such Rotation Matrix has dimension either equal to 3 or to 4, for 2D and 3D rotation, respectively. The {Number,1} of args either contain a single angle parameter in radiants, or a vector with three elements, whose norm is the rotation angle in 3D and whose normalized value gives the direction of the rotation axis in 3D.

Examples

julia> Lar.r(pi/6)				# 2D rotation of ``π/6`` angle
# return
3×3 Array{Float64,2}:
 0.866025  -0.5       0.0
 0.5        0.866025  0.0
 0.0        0.0       1.0

julia> Lar.r(0,0,pi/4)
# return
4×4 Array{Float64,2}:		# 3D rotation about the ``z`` axis, with ``π/6`` angle
 0.707107  -0.707107  0.0  0.0
 0.707107   0.707107  0.0  0.0
 0.0        0.0       1.0  0.0
 0.0        0.0       0.0  1.0

julia> Lar.r(1,1,1)		# 3D rotation about the ``x=y=z`` axis, with angle ``1.7320508`` angle
# return
4×4 Array{Float64,2}:
  0.226296  -0.183008   0.956712  0.0
  0.956712   0.226296  -0.183008  0.0
 -0.183008   0.956712   1.21332   0.0
  0.0        0.0        0.0       1.0
LinearAlgebraicRepresentation.rayintersectionMethod
rayintersection(point3d::Array{Float64})(V,FV,face::Int)

Compute the intersection point of the vertical line through point3d w face. If the face is parallel to z axis return false.

Example

julia> V,(VV,EV,FV,CV) = Lar.simplex(3,true);

julia> V
3×4 Array{Float64,2}:
 0.0  1.0  0.0  0.0
 0.0  0.0  1.0  0.0
 0.0  0.0  0.0  1.0

julia> FV
4-element Array{Array{Int64,1},1}:
 [1, 2, 3]
 [1, 2, 4]
 [1, 3, 4]
 [2, 3, 4]

 julia> Lar.rayintersection([.333,.333,0])(V,FV,4)
 3-element Array{Float64,1}:
  0.333
  0.333
  0.3340000000000001
LinearAlgebraicRepresentation.removeDupsMethod
removeDups(CW::Cells)::Cells

Remove dublicate cells from Cells object. Then put Cells in canonical form, i.e. with sorted indices of vertices in each (unique) Cells Array element.

LinearAlgebraicRepresentation.ringFunction
ring(r=1., R=2., angle=2*pi)(shape=[36, 1])

Compute the cellular 2-complex approximating a (possibly full) sector of a non-contractible disk. R and r are the external and the internal radiuses, respectively.

Example

julia> GL.VIEW([
	GL.GLGrid( Lar.ring()()..., GL.COLORS[1],1 ),
	GL.GLFrame
]);
LinearAlgebraicRepresentation.rodFunction
rod(radius=1, height=3, angle=2*pi)(shape=[36, 1])

Compute a cellular 3-complex with a single 3-cell starting from a cyclindrical surface generated with the same parameters.

Example

julia> rod()()[1]
# output
3×74 Array{Float64, 2}:
 -0.34202   0.984808   1.0          1.0  …   0.984808  -0.866025  -1.0           0.766044
  0.939693  0.173648  -2.44929e-16  0.0     -0.173648  -0.5        1.22465e-16  -0.642788
  3.0       3.0        0.0          0.0      0.0        3.0        3.0           0.0

julia> rod()()[2]
# output
1-element Array{Array{Int64, 1}, 1}:
 [0, 1, 2, 3, 4, 5, 6, 7, 8, 9  …  64, 65, 66, 67, 68, 69, 70, 71, 72, 73]

julia> GL.VIEW([
 	GL.GLPol( Lar.rod()()..., GL.COLORS[1],0.5 ),
 	GL.GLFrame ]);
LinearAlgebraicRepresentation.sMethod
s(args::Array{Number,1}...)::Matrix

Return an affine transformation Matrix in homogeneous coordinates. Such scaling Matrix has $d+1$ rows and $d+1$ columns, where $d$ is the number of scaling parameters in the args array.

Examples

julia> Lar.s(2,3)			# 2D scaling
# return
3×3 Array{Float64,2}:
 2.0  0.0  0.0
 0.0  3.0  0.0
 0.0  0.0  1.0

julia> Lar.s(2.,3.,4.)		# 3D scaling
# return
4×4 Array{Float64,2}:
 2.0  0.0  0.0  0.0
 0.0  3.0  0.0  0.0
 0.0  0.0  4.0  0.0
 0.0  0.0  0.0  1.0
LinearAlgebraicRepresentation.secondMomentMethod
secondMoment(P::Lar.LAR)::Array{Float64,1}

Second moments as terms of the Euler tensor.

Example # unit 3D tetrahedron

julia> V = [0.0 1.0 0.0 0.0; 0.0 0.0 1.0 0.0; 0.0 0.0 0.0 1.0];

julia> FV = [[1, 2, 4], [1, 3, 2], [4, 3, 1], [2, 3, 4]];

julia> P = V,FV;

julia> Lar.secondMoment(P)
3-element Array{Float64,1}:
 0.0166667
 0.0166667
 0.0166667
LinearAlgebraicRepresentation.setTileMethod
setTile(box)(point)

Set the tileCode of the 2D bbox [b1,b2,b3,b4]:=[ymax,ymin,xmax,xmin]:= x,x,y,y including the 2D point of x,y coordinates. Depending on point position, tileCode ranges in 0:15, and uses bit operators. Used to set the plane tiling depending on position of the query point, in order to subsequently test the tile codes of edges of a 2D polygon, and determine if the query point is either internal, external, or on the boundary of the polygon. Function to be parallelized ...

c1,c2 = tilecode(p1),tilecode(p2)
c_edge, c_un, c_int = c1 ⊻ c2, c1 | c2, c1 & c2
LinearAlgebraicRepresentation.sigmamodelMethod
sigmamodel(V::Lar.Points, copEV::Lar.ChainOP, FV::Lar.Cells, copFE::Lar.ChainOP, 
	sigma::Int, sp_idx::Array)

Transform sigma and sp_idx[sigma] faces to $z=0$ space.

LinearAlgebraicRepresentation.simplexFunction
simplex(n::Int, fullmodel=false::Bool)::Union{Lar.LAR, Lar.LARmodel}

Return a LAR model of the n-dimensional simplex in n-space.

When fullmodel==true return a LARmodel, including the faces, from dimension 1 to n.

Example

using LinearAlgebraicRepresentation, Plasm, LinearAlgebra
Lar = LinearAlgebraicRepresentation

model = Lar.simplex(2)
Plasm.view( Lar.simplex(2) )

V, cells = Lar.simplex(3, true)
Plasm.view(Plasm.numbering(0.5)( (V,cells[1:end-1]) ))
LinearAlgebraicRepresentation.simplexFacetsMethod
simplexFacets(simplices::Cells)::Cells

Compute the (d-1)-skeleton (unoriented set of facets) of a simplicial d-complex.

Example

julia> V,FV = Lar.simplexGrid([1,1]) # 2-dimensional complex
# output
([0 1 0 1; 0 0 1 1], Array{Int64,1}[[1, 2, 3], [2, 3, 4]])

julia> Plasm.view(V,FV)

julia> W,CW = Lar.extrudeSimplicial((V,FV), [1])
([0.0 1.0 … 0.0 1.0; 0.0 0.0 … 1.0 1.0; 0.0 0.0 … 1.0 1.0],
Array{Int64,1}[[1,2,3,5],[2,3,5,6],[3,5,6,7],[2,3,4,6],[3,4,6,7],[4,6,7,8]])

julia> FW = Lar.simplexFacets(CW)
18-element Array{Any,1}:
[[1,3,5],[5,6,7],[3,5,7],[3,6,7],[4,6,7],[4,7,8],[4,6,8],
[6,7,8],[3,5,6],[2,3,5],[2,3,4],[3,4,7],[1,2,3],[2,4,6],[2,5,6],
[1,2,5],[2,3,6],[3,4,6]]

julia> Plasm.view(W,FW)

Example

julia> V,(VV,EV,FV,CV) = Lar.cuboidGrid([3,3,3],true)

julia> TV = Lar.simplexFacets(CV)

julia> Plasm.view(V,TV)
LinearAlgebraicRepresentation.simplexGridMethod
simplexGrid(shape::Array)::LAR

Generate a simplicial complex decomposition of a cubical grid of $d$-cuboids, where $d$ is the length of shape array. Vertices (0-cells) of the grid have Int64 coordinates.

Examples

julia> simplexGrid([0]) # 0-dimensional complex
# output
([0], Array{Int64,1}[])

julia> V,EV = simplexGrid([1]) # 1-dimensional complex
# output
([0 1], Array{Int64,1}[[1, 2]])

julia> V,FV = simplexGrid([1,1]) # 2-dimensional complex
# output
([0 1 0 1; 0 0 1 1], Array{Int64,1}[[1, 2, 3], [2, 3, 4]])

julia> V,CV = simplexGrid([10,10,1]) # 3-dimensional complex
# output
([0 1 … 9 10; 0 0 … 10 10; 0 0 … 1 1], Array{Int64,1}[[1, 2, 12, 122], [2, 12, 122, 123], [12, 122, 123, 133], [2, 12, 13, 123], [12, 13, 123, 133], [13, 123, 133, 134], [2, 3, 13, 123], [3, 13, 123, 124], [13, 123, 124, 134], [3, 13, 14, 124]  …  [119, 229, 230, 240], [109, 119, 120, 230], [119, 120, 230, 240], [120, 230, 240, 241], [109, 110, 120, 230], [110, 120, 230, 231], [120, 230, 231, 241], [110, 120, 121, 231], [120, 121, 231, 241], [121, 231, 241, 242]])

julia> V
# output
3×242 Array{Int64,2}:
 0  1  2  3  4  5  6  7  8  9  10  0  1  2  3  …   1   2   3   4   5   6   7   8   9  10
 0  0  0  0  0  0  0  0  0  0   0  1  1  1  1     10  10  10  10  10  10  10  10  10  10
 0  0  0  0  0  0  0  0  0  0   0  0  0  0  0      1   1   1   1   1   1   1   1   1   1

julia> using Plasm

julia> hpc = Plasm.lar2exploded_hpc(V,CV) # exploded visualization of the grid

julia> Plasm.view(hpc)

julia> V,HV = simplexGrid([1,1,1,1]) # 4-dim cellular complex from the 4D simplex
# output
([0 1 … 0 1; 0 0 … 1 1; 0 0 … 1 1; 0 0 … 1 1], Array{Int64,1}[[1, 2, 3, 5, 9], [2, 3, 5, 9, 10], [3, 5, 9, 10, 11], [5, 9, 10, 11, 13], [2, 3, 5, 6, 10], [3, 5, 6, 10, 11], [5, 6, 10, 11, 13], [6, 10, 11, 13, 14], [3, 5, 6, 7, 11], [5, 6, 7, 11, 13]  …  [4, 6, 10, 11, 12], [6, 10, 11, 12, 14], [3, 4, 6, 7, 11], [4, 6, 7, 11, 12], [6, 7, 11, 12, 14], [7, 11, 12, 14, 15], [4, 6, 7, 8, 12], [6, 7, 8, 12, 14], [7, 8, 12, 14, 15], [8, 12, 14, 15, 16]])
LinearAlgebraicRepresentation.simplifyCellsMethod
W,CW = simplifyCells(V,CV)

Find and remove the duplicated vertices and the incorrect cells. Some vertices may appear two or more times, due to numerical errors on mapped coordinates. Close vertices are identified, according to the PRECISION number of significant digits.

LinearAlgebraicRepresentation.spaceindexMethod
spaceindex(point3d)(model)

Compute the set of face boxes of possible intersection with a point-ray. Work in 3D, where the ray direction is parallel to the z-axis. Return an array of indices of face.

# Example

julia> V,(VV,EV,FV,CV) = Lar.cuboidGrid([1,1,1],true)

julia> spaceindex([.5,.5,.5])((V,FV))
3-element Array{Int64,1}:
 5
 6
LinearAlgebraicRepresentation.spaceindexMethod
spaceindex(model::Lar.LAR)::Array{Array{Int,1},1}

Generation of space indexes for all $(d-1)$-dim cell members of model.

Spatial index made by $d$interval-trees on bounding boxes of $sigma in S_{d−1}$. Spatial queries solved by intersection of $d$ queries on IntervalTrees generated by bounding-boxes of geometric objects (LAR cells).

The return value is an array of arrays of ints, indexing cells whose containment boxes are intersecting the containment box of the first cell. According to Hoffmann, Hopcroft, and Karasick (1989) the worst-case complexity of Boolean ops on such complexes equates the total sum of such numbers.

Examples 2D

julia> V = hcat([[0.,0],[1,0],[1,1],[0,1],[2,1]]...);

julia> EV = [[1,2],[2,3],[3,4],[4,1],[1,5]];

julia> Sigma = Lar.spaceindex((V,EV))
5-element Array{Array{Int64,1},1}:
 [4, 5, 2]
 [1, 3, 5]
 [4, 5, 2]
 [1, 3, 5]
 [4, 1, 3, 2]

From model2d value, available in ?input_collection docstring:

julia> Sigma =  spaceindex(model2d);

Example 3D

model = model3d
Sigma =  spaceindex(model3d);
Sigma
LinearAlgebraicRepresentation.svg2larMethod
svg2lar(filename::String; flag=true)::Lar.LAR

Parse a SVG file to a LAR model (V,EV). Only <line > and <rect > and <path > SVG primitives are currently translated. TODO: interpretation of transformations.

LinearAlgebraicRepresentation.tMethod
t(args::Array{Number,1}...)::Matrix

Return an affine transformation Matrix in homogeneous coordinates. Such translation Matrix has $d+1$ rows and $d+1$ columns, where $d$ is the number of translation parameters in the args array.

Examples

julia> t(1,2)			# 2D translation
# return
3×3 Array{Float64,2}:
 1.0  0.0  1.0
 0.0  1.0  2.0
 0.0  0.0  1.0

julia> Lar.t(1.,2,3)		# 3D translation
# return
4×4 Array{Float64,2}:
 1.0  0.0  0.0  1.0
 0.0  1.0  0.0  2.0
 0.0  0.0  1.0  3.0
 0.0  0.0  0.0  1.0
LinearAlgebraicRepresentation.toroidalFunction
toroidal(r=1., R=2., angle1=2*pi, angle2=2*pi)(shape=[24, 36])

Compute a cellular 2-complex, approximation of the two-dimensional surface, embedded in a three-dimensional Euclidean space. Toroidal is a closed surface having genus one, and therefore possessing a single "hole". It can be constructed from a rectangle by gluing both pairs of opposite edges together with no twists.

Example

julia> GL.VIEW([
	GL.GLGrid( Lar.toroidal()()..., GL.COLORS[1],0.75 ),
	GL.GLFrame
]);
LinearAlgebraicRepresentation.torusFunction
torus(r=1., R=2., h=.5, angle1=2*pi, angle2=2*pi)(shape=[24, 36, 4])

Compute the cellular 3-complex approximating the solid torus in 3D. The model is meshed with cubical 3-cells, where the mesh has default decomposition size [24, 36, 4]. See also: toroidal. h is radius of the circular hole inside the solid.

Example

julia> GL.VIEW([
 	GL.GLPol( Lar.torus(1., 2., .5, pi, pi)()..., GL.COLORS[1],0.5 ),
 	GL.GLFrame ]);
LinearAlgebraicRepresentation.triangulateMethod
triangulate(V::Points, cc::ChainComplex)

Full constrained Delaunnay triangulation of the given 3-dimensional model (given with topology as a ChainComplex)

LinearAlgebraicRepresentation.u_coboundary_1Function
u_coboundary_1( FV::Cells, EV::Cells, convex=true)::ChainOp

Compute the sparse unsigned coboundary1 operator ``C1 -> C_2`. Notice that the output matrix ism x n, wheremis the number of faces, andn` is the number of edges.

Examples

Cellular complex with convex-cells, and without outer cell

julia> V,(VV,EV,FV,CV) = cuboid([1.,1.,1.], true);

julia> u_coboundary_1(FV,EV)
6×12 SparseMatrixCSC{Int8,Int64} with 24 stored entries:
[1 ,  1]  =  1
[3 ,  1]  =  1
[1 ,  2]  =  1
[4 ,  2]  =  1
...		...
[4 , 11]  =  1
[5 , 11]  =  1
[4 , 12]  =  1
[6 , 12]  =  1

julia> Matrix(u_coboundary_1(FV,EV))
6×12 Array{Int8,2}:
1  1  0  0  1  1  0  0  0  0  0  0
0  0  1  1  0  0  1  1  0  0  0  0
1  0  1  0  0  0  0  0  1  1  0  0
0  1  0  1  0  0  0  0  0  0  1  1
0  0  0  0  1  0  1  0  1  0  1  0
0  0  0  0  0  1  0  1  0  1  0  1

julia> unsigned_boundary_2 = u_coboundary_1(FV,EV)';

Compute the Unsignedcoboundary_1 operator matrix as product of two sparse characteristic matrices.

Cellular complex with non-convex cells, and with outer cell

FV = [[1,2,3,4,5,17,16,12], # outer cell
[1,2,3,4,6,7,8,9,10,11,12,13,14,15],
[4,5,9,11,12,13,14,15,16,17],
[2,3,6,7], [8,9,10,11]]

EV = [[1,2],[2,3],[3,4],[4,5],[1,12],[2,6],[3,7],[4,9],[5,17],[6,7],[8,9],
[8,10],[9,11],[10,11],[11,15],[12,13],[12,16],[13,14],[14,15],[16,17]]

out = u_coboundary_1( FV::Cells, EV::Cells, false)

In case of expected 2-chains with non-convex cells, instance the method with convex = false, in order to fix a possible redundancy of incidence values, induced by computation through multiplication of characteristic matrices. (Look at columns 2 and 13 before, generated by default).

LinearAlgebraicRepresentation.u_coboundary_2Function
u_coboundary_2( CV::Cells, FV::Cells[, convex=true::Bool] )::ChainOp

Unsigned 2-coboundary matrix ∂_2 : C_2 -> C_3 from 2-chain to 3-chain space. Compute algebraically the unsigned coboundary matrix ∂_2 from characteristic matrices of CV and FV. Currently usable only with complexes of convex cells.

# Examples

First example

(1) Compute the boundary matrix for a block of 3-cells of size $[32,32,16]$;

(2) compute and show the boundary 2-cell array boundary_2D_cells by decodifying the (mod 2) result of multiplication of the boundary_3 matrix∂_2', transpose of unsigned coboundary_2 matrix times the coordinate vector of the $total$ 3-chain.

julia> using SparseArrays, Plasm
julia> V,(_,_,FV,CV) = cuboidGrid([32,32,16], true)
julia> ∂_2 = u_coboundary_2( CV, FV)
julia> coord_vect_of_all_3D_cells  = ones(size(∂_2,1),1)
julia> coord_vect_of_boundary_2D_cells = ∂_2' * coord_vect_of_all_3D_cells .% 2
julia> out = coord_vect_of_boundary_2D_cells
julia> boundary_2D_cells = [ FV[f] for f in findnz(sparse(out))[1] ]
julia> hpc = Plasm.lar2exploded_hpc(V, boundary_2D_cells)(1.,1.,1.)
julia> Plasm.view(hpc)

Second example example

Using the boundary matrix of the 32 x 32 x 16 "image block" (better if stored on disk) compute the boundary 2-complex of a random sub-image inside the block.

julia> coord_vect_of_segment = [x>0.25 ? 1 : 0  for x in rand(size(∂_2,1)) ]
julia> out = ∂_2' * coord_vect_of_segment .% 2
julia> boundary_2D_cells = [ FV[f] for f in findnz(sparse(out))[1] ]
julia> hpc = Plasm.lar2exploded_hpc(V, boundary_2D_cells)(1.1,1.1,1.1)
julia> Plasm.view(hpc)
LinearAlgebraicRepresentation.verts2vertsMethod
verts2verts(EV::Lar.Cells)::Lar.Cells

Adjacency lists of vertices of a cellular 1-complex.

Example

julia> V,(VV,EV,FV) = Lar.cuboidGrid([3,3],true);

julia> verts2verts(EV::Lar.Cells)::Lar.Cells
16-element Array{Array{Int64,1},1}:
 [2, 5]         
 [1, 3, 6]      
 [2, 4, 7]      
 [3, 8]         
 [1, 6, 9]      
 [2, 5, 7, 10]  
 [3, 6, 8, 11]  
 [4, 7, 12]     
 [5, 10, 13]    
 [6, 9, 11, 14] 
 [7, 10, 12, 15]
 [8, 11, 16]    
 [9, 14]        
 [10, 13, 15]   
 [11, 14, 16]   
 [12, 15]       
LinearAlgebraicRepresentation.volumeMethod
volume(P::Lar.LAR)::Float64

volume integral on polyhedron P.

Example # unit 3D tetrahedron

julia> V = [0.0 1.0 0.0 0.0; 0.0 0.0 1.0 0.0; 0.0 0.0 0.0 1.0]
3×4 Array{Float64,2}:
 0.0  1.0  0.0  0.0
 0.0  0.0  1.0  0.0
 0.0  0.0  0.0  1.0

julia> FV = [[1, 2, 4], [1, 3, 2], [4, 3, 1], [2, 3, 4]]
4-element Array{Array{Int64,1},1}:
 [1, 2, 4]
 [1, 3, 2]
 [4, 3, 1]
 [2, 3, 4]

julia> P = V,FV
([0.0 1.0 0.0 0.0; 0.0 0.0 1.0 0.0; 0.0 0.0 0.0 1.0], 
Array{Int64,1}[[1, 2, 4], [1, 3, 2], [4, 3, 1], [2, 3, 4]])

julia> Lar.volume(P)
0.16666666666666674