Private functions

Base.epsMethod
Base.eps(T::Type{<:AbstractFloat},grid::Grid)

Return the machine roundoff error of a Grid. It returns eps() of the maximum length of the bounding box of the Grid.

source
Base.splitFunction
split(p::Polyhedron,plane;kwargs...)

It splits a polyhedron by a plane into two polyhedra.

It returns a tuple of Union{Polyhedron,Nothing}. If one side is empty, it returns nothing for that side.

Optional keyword arguments

  • side::Symbol=:both: It returns :both sides, the :left side, or the :right side
  • invert::Bool=false: It inverts the plane
source
GridapEmbedded.Interfaces.cutMethod
cut([cutter::STLCutter,]model::DiscreteModel,geo::STLGeometry)

Main interface of GridapEmbedded. It intersect each cell in model with the geometry geo. It returns an embedded discretization.

Usage

The basic usage is to call cut with a model and a geometry.

cutgeo = cut(model,geo)

We can also set custom keyword arguments, see subtriangulate.

cutter = STLCutter(;tolfactor=1e3)
cutgeo = cut(cutter,model,geo)

Then we can extract the domain and the boundary of the embedded discretization.

Ω = Triangulation(cutgeo)
Γ = EmbeddedBoundary(cutgeo)
source
STLCutters._group_verticesMethod
_group_vertices(grid::Grid;atol)

Group a list of vertices which are closer than atol. It returns a vector of vectors with the indices of grouped vertices.

To perform the classification, this function utilizes an auxiliar CartesianGrid of size h>atol.

source
STLCutters._map_equal_vertices_from_cloudMethod
_map_equal_vertices_from_cloud(stl::DiscreteModel;atol])

Find vertices closer than a tolerance atol. It returns a vector from each vector to the index of a vertex which is closer than atol.

In order to find the vertices without quadratic complexity, it requires to group the vertices in a grid with an arbitrary tolerance > atol.

source
STLCutters.bisector_planeMethod
bisector_plane(e::Face,Π1::AbstractPlane,Π2::AbstractPlane)

It returns the bisector plane between two planes. The bisector contains the edge e.

source
STLCutters.bisector_planeMethod
bisector_plane(stl::STL,d::Integer,dface::Integer,Πf::AbstractArray)

Compute the bisector plane of two facets in a STL model. The two facets must be connected by a STL edge (dface).

Arguments

  • stl::STL: The surface model.
  • d::Integer: The dimension of the edge (0 for 2D, 1 for 3D).
  • dface::Integer: The index of the edge.
  • Πf::AbstractArray: The array of planes of the facets of the STL
source
STLCutters.clipMethod
clip(p::Polyhedron,planes;kwargs...) -> Polyhedron

It clips a polyhedron by the halfspace of a plane or a set of planes.

Optional keyword arguments

  • inside::Bool=true: It clips the polyhedron by the inside ourside of the union of halfspace.
  • inout::Vector{Bool}=trues(length(planes)): In reverses the halfspaces with inside[i]=false.
  • boundary::nothing:
    • If boundary=true, it preserves the vertices on the planes (zero distance)
    • If boundary=false, it removes the vertices on the planes (zero distance)
    • If boundary=nothing, it preserves the vertices on the planes if the normal of the plane points inwards.
source
STLCutters.closest_point!Function
closest_point!(cache,point::Point,stl::DiscreteModel,[faces])

It returns the closest point to point in the stl.

Optional arguments

  • faces::Vector{Int}: List of faces to search the closest point.
source
STLCutters.closest_pointMethod
closest_point(points,stl::DiscreteModel)

It returns an array of Point which are the closest_point! in the stl for each element in points .

Optional arguments

  • point_to_faces: Vector of vectors of the stl faces to search for each point.
source
STLCutters.closest_pointMethod
closest_point(p::Point,f::Face)

It returns the closest point in a face to a point. It returns the projection if it is inside the face, otherwise it returns the projection to the boundary.

source
STLCutters.complete_in_or_out!Method
complete_in_or_out!(a::AbstractEmbeddedDiscretization,b::AbstractEmbeddedDiscretization,acells,bcells)

This function considers two discretizations of the same background model on the cells acells and bcells. These sets of cells have null intersection. If only one of the discretizations is not cut it sets the other one as in or out.

source
STLCutters.compute_cartesian_descriptorMethod
compute_cartesian_descriptor(pmin,pmax;nmin,nmax)

Compute CartesianDescriptor in a bounding box with the same cell size (h) in all directions:

h = min( max((pmin-pmax)/nmax), min((pmax-pmin)/nmin) )

source
STLCutters.compute_cell_to_facetsFunction
compute_cell_to_facets(a::UnstructuredGrid,b::Grid[,a_mask])

compute_cell_to_facets computes the cells in b colliding each cell in a:UnstructuredGrid. The output is a vector of vectors.

Note

This function uses a CartesianGrid internally. It is not optimized for higly irregular grids.

Note

This function allows false positives.

source
STLCutters.compute_cell_to_facetsFunction
compute_cell_to_facets(a::CartesianGrid,b::Grid[,a_mask,b_mask])

compute_cell_to_facets computes the cells in b colliding each cell in a using optimizations for a CartesianGrid. The output is a vector of vectors.

!!! note This function allows false positives.

source
STLCutters.compute_cell_to_facetsMethod
compute_cell_to_facets(a::DiscreteModel,b::Grid,args...)

Compute a map of cells in a to the cells in b which potentially intersect. It is designed to filter the STL faces (in b) colliding each background cell in a. The output is a vector of vectors.

Note

This function allows false positives.

source
STLCutters.compute_cell_to_isboundaryFunction
compute_cell_to_isboundary(model::DiscreteModel,indices[,d=0])

It returns a mask whether a cell touches the subdomain interface.

Arguments

  • model::DiscreteModel: Model of the subdomain.
  • indices::AbstractLocalIndices: Partition indices
  • d::Integer=0: Dimemension of the d-faces touching the subdomain interface.
source
STLCutters.compute_distances!Method
compute_distances!(p::Polyhedron,planes,plane_ids[;atol=0])

Compute the distances from the vertices of the polyhedron to each plane of the list planes. If the distance is below atol, it sets the distance to zero.

source
STLCutters.compute_face_neighbor_to_inoutcutFunction
compute_face_neighbor_to_inoutcut(model::DiscreteModel,indices,d,face_to_inoutcut[,face_neighbors])

It returns a whether interfaces of dimension $d$ are in, out, cut or undefined.

The number of interfaces coincides with the number of neigbors given by compute_face_neighbors

Note

If the subdomain is not cut, the neighbors are considered undefined.

source
STLCutters.compute_face_neighborsMethod
compute_face_neighbors(model::DiscreteModel,indices,d)

It returns a neighboring graph of the subdomain's neighbors through the interfaces of dimension d.

source
STLCutters.compute_gridMethod
compute_grid(coordinates,connectivity,polytope)

Compute a Grid from coordinates, connectivity and polytope. It assumes that all the cells of the grid have the same polytope.

source
STLCutters.compute_modelMethod
compute_model(coordinates,connectivity,polytope)

Compute a DiscreteModel from coordinates, connectivity and polytope. It assumes that all the cells of the model have the same polytope.

source
STLCutters.compute_stl_modelMethod
compute_stl_model(coordinates,connectivity)

Compute a DiscreteModel from coordinates and connectivity. It assumes that the model is surface composed of simplex faces.

source
STLCutters.displaceFunction
displace(plane::AbstractPlane,dist[,oriented=true])

Move a plane a distance dist along its normal vector. If oriented is false, move a distance -dist.

source
STLCutters.distanceMethod
distance(a,b)

It returns the minimum absolute distance between two entities (Face or Point).

source
STLCutters.expand_planesMethod
explan_planes(planes::AbstractVector{<:AbstractPlane},mask,dist)

Move a list of planes a distance dist along their normal vectors. If the mask entries are false, the planes are moved a distance -dist.

source
STLCutters.find_root_partMethod
find_root_part(cuts::AbstractArray{<:AbstractEmbeddedDiscretization},cells)

This functions sets the root processor to a potentially idling processor. It returns the first non-cut part.

source
STLCutters.get_cell_planesMethod
get_cell_planes(p::Polytope,coords::AbstractVector{<:Point})

It returns a list of planes bounding a cell. The cell is defind by a Polytope and a list of coordinates. If the cell is a voxel, the bounding planes are CartesianPlane.

source
STLCutters.get_dfaceMethod
get_dface(f::Face,i::Integer,::Val{d})

It returns the ith d-face of a face f. E.g., vertex, edge, facet, cell.

source
STLCutters.intersection_pointMethod
intersection_point(a,b)

It returns the intersection point of two faces, or a face and a plane. The intersection point is only well defined if num_dims(a)+num_dims(b)-D=0, where D is the dimension of the space.

source
STLCutters.istouchedFunction
istouched(cut::AbstractEmbeddedDiscretization[,cells])

Check is an embedded discretization cut (or a set of cells) have information about the geometry. In other words, it checks if any cell is cut.

source
STLCutters.measureMethod
measure(a::Grid,mask)

Compute the sum of the cell measures of a grid. It only considers the cells which are true in the mask.

source
STLCutters.measureMethod
measure(f::Face)

It returns the measure of a face. E.g., the length, the area or the volume.

source
STLCutters.merge_nodesMethod
merge_nodes(stl::DiscreteModel [,atol])

Merge close nodes of the stl which are closer than a tolerance atol.

source
STLCutters.min_heightMethod
min_height(f::Face)

It returns the minimum height of a face, i.e., the minimum distance between two opposite entities of the face.

source
STLCutters.normalMethod
normal(a::Face[,i])

It returns the normal vector of a face (or its i-th facet). Note, that the normal is only defined for D-1 faces.

source
STLCutters.preprocess_small_facetsMethod
preprocess_small_facets(stl::DiscreteModel;atol)

Fix small facet problems of a quasi-degenerated stl.

If an edge and an opposite vertes are closer than atol, the edge and the face are split. The conformity is kept in the surrounding faces.

source
STLCutters.propagate_inout!Method
propagate_inout!(bgmocel,bgcell_to_inoutcut,bgnode_to_inoutcut,in_or_out)

Propagates IN or OUT (in_or_out) label of each backgroun cell (bgcell_to_inoutcut) through the nodes of the background model (bgnode_to_inoutcut). It overwrites bgcell_to_inoutcut`.

It is implemented with a deep-first search algorithm that stops at CUT cells.

source
STLCutters.read_stlMethod
read_stl(filename)

Read STL file with MeshIO.jl and Gridap.jl arrays:

  • vertex_to_coordinates is a vector of Point{3,Float64}
  • facet_to_vertices is a Gridap.Table of Int
  • facet_to_normals is a vector of VectorValue{3,Float64}
source
STLCutters.save_as_stlMethod
save_as_stl(stl::Vector{<:DiscreteModel},filename)

Save a list of DiscreteModel{2,3} into STL files. It save the files as filename_1.stl, filename_2.stl,..., filename_n.stl`.

source
STLCutters.set_in_or_out!Method
set_in_or_out!(cut::DistributedEmbeddedDiscretization,in_or_out[,cells,facets])

Sets all the background cells and facets (or a set of cells and facets) as in or out.

source
STLCutters.simplex_faceMethod
simplex_face(x::Point...)

It returns a simplex face from a list (or a tuple) of points. The dimension of the simplex is N+1, where N is the number of points.

source
STLCutters.simplexify_boundaryMethod
simplexify_boundary(p::Polyhedron,t::GridTopology)

It generates a simplex mesh of the surface of the polyhedron that touches the grid topology.

It returns the coordinates of the vertices, the connectivities and a mapping to the topology facets.

source
STLCutters.voxel_intersectionMethod
voxel_intersection(f::Face,pmin::Point,pmax::Point,p::Polytope)

Predicate that checks if a face intersects by a voxel. Here, a voxel is a hypercube defined by its extrema (pmin and pmax).

Note

This function can return false positives.

source