vtf-logo

geom::IndSimpSet< _N, _M, _A, T, V, IS > Class Template Reference

Class for a mesh that stores vertices and indexed simplices. More...

#include <IndSimpSet.h>

Inheritance diagram for geom::IndSimpSet< _N, _M, _A, T, V, IS >:

geom::IndSimpSetIncAdj< _N, _M, _A, T, V, IS > List of all members.

Public Types

typedef T Number
 The number type.
typedef V Vertex
 A vertex.
typedef Simplex< M, VertexSimplex
 A simplex of vertices.
typedef Simplex::Face SimplexFace
 The face of a simplex of vertices.
typedef IS IndexedSimplex
 An indexed simplex. (A simplex of indices.).
typedef IndexedSimplex::Face IndexedSimplexFace
 The face of an indexed simplex.
typedef ads::Array< 1, Vertex,
A > 
VertexContainer
 The vertex container.
typedef VertexContainer::const_iterator VertexConstIterator
 A vertex const iterator.
typedef VertexContainer::iterator VertexIterator
 A vertex iterator.
typedef ads::Array< 1, IndexedSimplex,
A > 
IndexedSimplexContainer
 The indexed simplex container.
typedef IndexedSimplexContainer::const_iterator IndexedSimplexConstIterator
 An indexed simplex const iterator.
typedef IndexedSimplexContainer::iterator IndexedSimplexIterator
 An indexed simplex iterator.
typedef SimplexIterator< IndSimpSetSimplexConstIterator
 A simplex const iterator.
typedef VertexContainer::size_type SizeType
 The size type.
 N = _N
 M = _M
 A = _A
enum  { N = _N, M = _M, A = _A }
 The space dimension, simplex dimension and allocation.

Public Member Functions

Constructors etc.
Suppose that we are dealing with a tetrahedron mesh in 3-D. Below we instantiate a mesh that allocates its own memory for the vertices and indexed simplices. We can construct the mesh from vertices and indexed simplices stored in ADS arrays:
    typedef geom::IndSimpSet<3,3> ISS;
    typedef typename ISS:Vertex Vertex;
    typedef typename ISS:IndexedSimplex IndexedSimplex;
    ads::Array<1,Vertex> vertices(numberOfVertices);
    ads::Array<1,IndexedSimplex> indexedSimplices(numberOfSimplices);
    ...
    geom::IndSimpSet<3,3> mesh(vertices, indexedSimplices);
or use C arrays:
    double* vertices = new[3 * numberOfVertices]
    int* indexedSimplices = new[4 * numberOfSimplices];
    ...
    geom::IndSimpSet<3,3> mesh(numberOfVertices, vertices, numberOfSimplices, simplices);

We can also make meshes that borrow externally allocated memory. One can use ADS arrays:

    geom::IndSimpSet<3,3,false> mesh(vertices, indexedSimplices);
or C arrays:
    geom::IndSimpSet<3,3,false> mesh(numberOfVertices, vertices, numberOfSimplices, indexedSimplices);


 IndSimpSet ()
 Default constructor. Empty simplex set.
template<bool A1, bool A2>
 IndSimpSet (const ads::Array< 1, Vertex, A1 > &vertices, const ads::Array< 1, IndexedSimplex, A2 > &indexedSimplices)
 Construct from arrays of vertices and indexed simplices.
template<bool A1, bool A2>
void build (const ads::Array< 1, Vertex, A1 > &vertices, const ads::Array< 1, IndexedSimplex, A2 > &indexedSimplices)
 Build from arrays of vertices and indexed simplices.
 IndSimpSet (const SizeType numVertices, void *vertices, const SizeType numSimplices, void *indexedSimplices)
 Construct from pointers to the vertices and indexed simplices.
void build (const SizeType numVertices, void *vertices, const SizeType numSimplices, void *indexedSimplices)
 Build from pointers to the vertices and indexed simplices.
 IndSimpSet (const SizeType numVertices, const void *vertices, const SizeType numSimplices, const void *indexedSimplices)
 Construct from pointers to the vertices and indexed simplices.
void build (const SizeType numVertices, const void *vertices, const SizeType numSimplices, const void *indexedSimplices)
 Build from pointers to the vertices and indexed simplices.
 IndSimpSet (const SizeType numVertices, const SizeType numSimplices)
 Construct from the number of vertices and simplices.
void build (const SizeType numVertices, const SizeType numSimplices)
 Build from the number of vertices and simplices.
void swap (IndSimpSet &x)
 Swap data with another mesh.
 IndSimpSet (const IndSimpSet &other)
 Copy constructor.
IndSimpSetoperator= (const IndSimpSet &other)
 Assignment operator.
virtual ~IndSimpSet ()
 Destructor. Deletes memory only if it was allocated internally.
Vertex Accessors
int getSpaceDimension () const
 Return the dimension of the space.
SizeType getVerticesSize () const
 Return the number of vertices.
SizeType areVerticesEmpty () const
 Return true if there are no vertices.
VertexConstIterator getVerticesBeginning () const
 Return a const iterator to the beginning of the vertices.
VertexConstIterator getVerticesEnd () const
 Return a const iterator to the end of the vertices.
const VertexgetVertex (const int n) const
 Return a const reference to the n_th vertex.
const VertexContainergetVertices () const
 Return a const reference to the vertex container.
Simplex Accessors
int getSimplexDimension () const
 Return the dimension of the simplices.
SizeType getSimplicesSize () const
 Return the number of simplices.
SizeType areSimplicesEmpty () const
 Return true if there are no simplices.
SizeType isEmpty () const
 Return true if there are no vertices or simplices.
IndexedSimplexConstIterator getIndexedSimplicesBeginning () const
 Return a const iterator to the beginning of the indexed simplices.
IndexedSimplexConstIterator getIndexedSimplicesEnd () const
 Return a const iterator to the end of the indexed simplices.
const IndexedSimplexgetIndexedSimplex (const int n) const
 Return a const reference to the n_th indexed simplex.
const IndexedSimplexContainergetIndexedSimplices () const
 Return a const reference to the indexed simplex container.
SimplexConstIterator getSimplicesBeginning () const
 Return a const iterator to the beginning of the simplices.
SimplexConstIterator getSimplicesEnd () const
 Return a const iterator to the end of the simplices.
const VertexgetSimplexVertex (const int n, const int m) const
 Return a const reference to the m_th vertex of the n_th simplex.
void getSimplex (const int n, Simplex *s) const
 Get the n_th simplex.
void getSimplex (IndexedSimplexConstIterator i, Simplex *s) const
 Get the simplex given an iterator to the indexed simplex.
Vertex Manipulators
VertexIterator getVerticesBeginning ()
 Return an iterator to the beginning of the vertices.
VertexIterator getVerticesEnd ()
 Return an iterator to the end of the vertices.
void setVertex (const int n, const Vertex &vertex)
 Set the specified vertex.
VertexContainergetVertices ()
 Return a reference to the vertex container.
Simplex Manipulators
IndexedSimplexIterator getIndexedSimplicesBeginning ()
 Return an iterator to the beginning of the indexed simplices.
IndexedSimplexIterator getIndexedSimplicesEnd ()
 Return an iterator to the end of the indexed simplices.
IndexedSimplexContainergetIndexedSimplices ()
 Return a reference to the indexed simplex container.
Update the topology.
virtual void updateTopology ()
 Update the data structure following a change in the topology.

Related Functions

(Note that these are not member functions.)

void getFace (const IndSimpSet< N, M, A, T, V, IS > &mesh, int simplexIndex, int vertexIndex, typename IndSimpSet< N, M, A, T, V, IS >::SimplexFace *face)
 Get the face in the simplex that is opposite to the given vertex.
void getIndexedFace (const IndSimpSet< N, M, A, T, V, IS > &mesh, int simplexIndex, int vertexIndex, typename IndSimpSet< N, M, A, T, V, IS >::IndexedSimplexFace *face)
 Get the indexed face in the simplex that is opposite to the given vertex.
void getCentroid (const IndSimpSet< N, M, A, T, V, IS > &mesh, int n, typename IndSimpSet< N, M, A, T, V, IS >::Vertex *x)
 Get the centroid of the specified simplex.
void buildFromQuadMesh (const QuadMesh< N, A, T, V, IF > &quadMesh, IndSimpSet< N, 2, true, T, V, IS > *mesh)
 Build from a quadrilateral mesh.
void buildFromSubsetVertices (const IndSimpSet< N, M, A, T, V, IS > &in, IntForIter verticesBeginning, IntForIter verticesEnd, IndSimpSet< N, M, true, T, V, IS > *out)
 Make a mesh from a subset of vertices of a mesh.
void buildFromSubsetSimplices (const IndSimpSet< N, M, A, T, V, IS > &in, IntForIter simplicesBeginning, IntForIter simplicesEnd, IndSimpSet< N, M, true, T, V, IS > *out)
 Make a new mesh from the subset of simplices.
void buildFromVerticesInside (const IndSimpSet< N, M, A, T, V, IS > &in, const LSF &f, IndSimpSet< N, M, true, T, V, IS > *out)
 Make a mesh by selecting vertices from the input mesh that are inside the object.
void buildFromSimplicesInside (const IndSimpSet< N, M, A, T, V, IS > &in, const LSF &f, IndSimpSet< N, M, true, T, V, IS > *out)
 Make a mesh by selecting simplices from the input mesh that are inside the object.
void buildBoundary (const IndSimpSetIncAdj< N, M, A, T, V, ISimp > &in, IndSimpSet< N, M-1, true, T, V, IFace > *out, IntOutputIterator usedVertexIndices)
 Make a mesh that is the boundary of the input mesh.
void buildBoundary (const IndSimpSetIncAdj< N, M, A, T, V, ISimp > &in, IndSimpSet< N, M-1, true, T, V, IFace > *out)
 Make a mesh that is the boundary of the input mesh.
void buildBoundaryWithoutPacking (const IndSimpSetIncAdj< N, M, A, T, V, ISimp > &in, IndSimpSet< N, M-1, true, T, V, IFace > *out, IntOutputIterator incidentSimplices)
 Make a mesh that is the boundary of the input mesh.
void buildBoundaryWithoutPacking (const IndSimpSetIncAdj< N, M, A, T, V, ISimp > &in, IndSimpSet< N, M-1, true, T, V, IFace > *out)
 Make a mesh that is the boundary of the input mesh.
void buildBoundaryOfComponentsWithoutPacking (const IndSimpSetIncAdj< N, M, A, T, V, ISimp > &in, IndSimpSet< N, M-1, true, T, V, IFace > *out, IntOutputIterator1 delimiterIterator, IntOutputIterator2 incidentSimplices)
 Make a mesh (separated into connected components) that is the boundary of the input mesh.
void buildBoundaryOfComponentsWithoutPacking (const IndSimpSetIncAdj< N, M, A, T, V, ISimp > &in, IndSimpSet< N, M-1, true, T, V, IFace > *out, IntOutputIterator delimiterIterator)
 Make a mesh (separated into connected components) that is the boundary of the input mesh.
void buildBoundaryOfComponents (const IndSimpSetIncAdj< N, M, A, T, V, ISimp > &in, IndSimpSet< N, M-1, true, T, V, IFace > *out, IntOutputIterator1 delimiterIterator, IntOutputIterator2 incidentSimplices)
 Make a mesh (separated into connected components) that is the boundary of the input mesh.
void buildBoundaryOfComponents (const IndSimpSetIncAdj< N, M, A, T, V, ISimp > &in, IndSimpSet< N, M-1, true, T, V, IFace > *out, IntOutputIterator delimiterIterator)
 Make a mesh (separated into connected components) that is the boundary of the input mesh.
void centerPointMesh (const IndSimpSet< N, M-1, A, T, V, IFace > &boundary, IndSimpSet< N, M, true, T, V, ISimp > *mesh)
 Make a mesh by connecting the boundary nodes to a new center point.
void merge (MeshInputIterator beginning, MeshInputIterator end, IndSimpSet< N, M, true, T, V, IS > *out)
 Merge a range of meshes to make a single mesh.
void merge2 (const IndSimpSet< N, M, A, T, V, IS > &a, const IndSimpSet< N, M, A, T, V, IS > &b, IndSimpSet< N, M, true, T, V, IS > *out)
 Merge two meshes to make a single mesh.
void buildFromSimplices (VertexForIter verticesBeginning, VertexForIter verticesEnd, IndSimpSet< N, M, true, T, V, IS > *mesh)
 Build from the simplices.
void computeMeanRatio (const IndSimpSet< N, M, A, T, V, IS > &iss, OutputIterator output)
 Calculate the mean ratio function for each simplex in the mesh.
void computeModifiedMeanRatio (const IndSimpSet< N, M, A, T, V, IS > &iss, OutputIterator output)
 Calculate the modified mean ratio function for each simplex in the mesh.
void computeConditionNumber (const IndSimpSet< N, M, A, T, V, IS > &iss, OutputIterator output)
 Calculate the condition number function for each simplex in the mesh.
void computeModifiedConditionNumber (const IndSimpSet< N, M, A, T, V, IS > &iss, OutputIterator output)
 Calculate the modified condition number function for each simplex in the mesh.
void computeContent (const IndSimpSet< N, M, A, T, V, IS > &iss, OutputIterator output)
 Calculate the content for each simplex in the mesh.
int removeContact (const IndSimpSet< N, N-1, A, T, V, IS > &surface, VertexForwardIterator verticesBeginning, VertexForwardIterator verticesEnd)
 Move the vertices to remove contact.
bool operator== (const IndSimpSet< N, M, A1, T, V1, IS1 > &x, const IndSimpSet< N, M, A2, T, V2, IS2 > &y)
 Return true if the vertices and indexed simplices are equal.
bool operator!= (const IndSimpSet< N, M, A1, T, V1, IS1 > &x, const IndSimpSet< N, M, A2, T, V2, IS2 > &y)
 Return true if the vertices and indexed simplices are not equal.
void writeAscii (std::ostream &out, const IndSimpSet< N, M, A, T, V, IS > &x)
 Write an indexed simplex set in ascii format.
void writeBinary (std::ostream &out, const IndSimpSet< N, M, A, T, V, IS > &x)
 Write an indexed simplex set in binary format.
void readAscii (std::istream &in, IndSimpSet< N, M, true, T, V, IS > *x)
 Read an indexed simplex set in ascii format.
void readBinary (std::istream &in, IndSimpSet< N, M, true, T, V, IS > *x)
 Read an indexed simplex set in binary format.
void writeVtkXml (std::ostream &out, const IndSimpSet< N, M, A, T, V, IS > &x)
 Write in VTK XML unstructured grid format.
void writeVtkXml (std::ostream &out, const IndSimpSet< N, M, A, T, V, IS > &x, const ads::Array< 1, F, A2 > &cellData, std::string dataName)
 Write the mesh and the cell data in VTK XML unstructured grid format.
void writeVtkXml (std::ostream &out, const IndSimpSet< N, M, A, T, V, IS > &x, ContainerIter cellDataContainersBeginning, ContainerIter cellDataContainersEnd, StringIter dataNamesBeginning, StringIter dataNamesEnd)
 Write the mesh and the cell data in VTK XML unstructured grid format.
void writeVtkLegacy (std::ostream &out, const IndSimpSet< N, M, A, T, V, IS > &x, std::string title="")
 Write in legacy VTK unstructured grid format.
void fit (IndSimpSetIncAdj< 2, 2, A, T, V, IS > *mesh, const ISS_SignedDistance< ISS, 2 > &signedDistance, T deviationTangent, int numSweeps)
 Fit the boundary of a mesh to a level-set description.
void determineVerticesOnManifold (const IndSimpSet< N, M, A, T, V, IS > &mesh, const IndSimpSet< N, MM, MA, MT, MV, MIS > &manifold, IntOutIter indexIterator, T epsilon=std::sqrt(std::numeric_limits< T >::epsilon()))
 Get the vertices on the manifold.
void determineVerticesOnManifold (const IndSimpSet< N, M, A, T, V, IS > &mesh, IntInIter indicesBeginning, IntInIter indicesEnd, const IndSimpSet< N, MM, MA, MT, MV, MIS > &manifold, IntOutIter indexIterator, T epsilon=std::sqrt(std::numeric_limits< T >::epsilon()))
 Get the vertices (from the set of vertices) on the manifold.
void countAdjacencies (const IndSimpSetIncAdj< N, M, A, T, V, IS > &iss, ads::FixedArray< M+2, int > *counts)
 Calculate the adjacency counts for the simplices in the mesh.
computeMinimumEdgeLength (const IndSimpSet< N, M, A, T, V, IS > &mesh)
 Return the minimum edge length.
computeMaximumEdgeLength (const IndSimpSet< N, M, A, T, V, IS > &mesh)
 Return the maximum edge length.
computeContent (const IndSimpSet< N, M, A, T, V, IS > &iss)
 Return the total content of the simplices in the mesh.
void computeContentStatistics (const IndSimpSet< N, M, A, T, V, IS > &iss, T *minimumContent, T *maximumContent, T *meanContent)
 Calculate content (hypervolume) statistics for the simplices in the mesh.
void computeDeterminantStatistics (const IndSimpSet< N, M, A, T, V, IS > &iss, T *minimumDeterminant, T *maximumDeterminant, T *meanDeterminant)
 Calculate determinant statistics for the simplices in the mesh.
void computeModifiedMeanRatioStatistics (const IndSimpSet< N, M, A, T, V, IS > &iss, T *minimumModMeanRatio, T *maximumModMeanRatio, T *meanModMeanRatio)
 Calculate modified mean ratio function statistics for the simplices in the mesh.
void computeModifiedConditionNumberStatistics (const IndSimpSet< N, M, A, T, V, IS > &iss, T *minimumModCondNum, T *maximumModCondNum, T *meanModCondNum)
 Calculate modified condition number function statistics for the simplices in the mesh.
void computeQualityStatistics (const IndSimpSet< N, M, A, T, V, IS > &iss, T *minimumContent, T *maximumContent, T *meanContent, T *minimumDeterminant, T *maximumDeterminant, T *meanDeterminant, T *minimumModMeanRatio, T *maximumModMeanRatio, T *meanModMeanRatio, T *minimumModCondNum, T *maximumModCondNum, T *meanModCondNum)
 Calculate quality statistics for the simplices in the mesh.
void printQualityStatistics (std::ostream &out, const IndSimpSet< N, M, A, T, V, IS > &mesh)
 Print quality statistics for the simplices in the mesh.
void printInformation (std::ostream &out, const IndSimpSet< N, M, A, T, V, IS > &mesh)
 Print information about the mesh.
void determineVerticesInside (const IndSimpSet< N, M, A, T, V, IS > &mesh, const LSF &f, IntOutIter indexIterator)
 Determine the vertices that are inside the object.
void determineSimplicesInside (const IndSimpSet< N, M, A, T, V, IS > &mesh, const LSF &f, IntOutIter indexIterator)
 Determine the simplices that are inside the object.
void determineSimplicesThatSatisfyCondition (const IndSimpSet< N, M, A, T, V, IS > &mesh, const UnaryFunction &f, IntOutIter indexIterator)
 Determine the simplices which satisfy the specified condition.
void determineOverlappingSimplices (const IndSimpSet< N, M, A, T, V, IS > &mesh, const BBox< N, T > &domain, IntOutIter indexIterator)
 Determine the simplices whose bounding boxes overlap the domain.
void determineIncidentVertices (const IndSimpSet< N, M, A, T, V, IS > &mesh, IntInIter simplexIndicesBeginning, IntInIter simplexIndicesEnd, IntOutIter vertexIndicesIterator)
 Add the vertices which are incident to the simplices.
void subdivide (const IndSimpSet< N, 1, A, T, V, IS > &in, IndSimpSet< N, 1, true, T, V, IS > *out)
 Subdivide by splitting each simplex in half.
void subdivide (const IndSimpSetIncAdj< N, 2, A, T, V, IS > &in, IndSimpSet< N, 2, true, T, V, IS > *out)
 Subdivide by splitting each simplex into four similar simplices.
void tile (const BBox< 2, T > &domain, T length, const LSF &f, IndSimpSet< 2, 2, true, T, V, IS > *mesh)
 Tile the object with equilateral triangles.
void tile (const BBox< 2, T > &domain, const T length, IndSimpSet< 2, 2, true, T, V, IS > *mesh)
 Tile the rectangular region with equilateral triangles.
void tile (const BBox< 3, T > &domain, const T length, const LSF &f, IndSimpSet< 3, 3, true, T, V, IS > *mesh)
 Tile the object with a body-centered cubic lattice.
void tile (const BBox< 3, T > &domain, const T length, IndSimpSet< 3, 3, true, T, V, IS > *mesh)
 Tile the rectilinear region with a body-centered cubic lattice.
void transferIndices (const IndSimpSet< N, M, A, T, V, IS > &mesh, const PointArray &points, IndexArray *indices)
 Determine the simplex indices in an ISS for transfering fields.
void transfer (const IndSimpSet< N, M, A, T, V, IS > &mesh, const SourceFieldArray &sourceFields, const PointArray &points, TargetFieldArray *targetFields)
 Transfer fields for indexed simplex sets.
void pack (IndSimpSet< N, M, true, T, V, IS > *mesh, IntOutputIterator usedVertexIndices)
 Pack the ISS to get rid of unused vertices.
void pack (IndSimpSet< N, M, true, T, V, IS > *mesh)
 Pack the ISS to get rid of unused vertices.
void orientPositive (IndSimpSet< N, N, A, T, V, IS > *mesh)
 Orient each simplex so it has non-negative volume.
void reverseOrientation (IndSimpSet< N, M, A, T, V, IS > *mesh)
 Reverse the orientation of each simplex.
void transform (IndSimpSet< N, M, A, T, V, IS > *mesh, IntForIter beginning, IntForIter end, const UnaryFunction &f)
 Transform each vertex in the range with the specified function.
void transform (IndSimpSet< N, M, A, T, V, IS > mesh, const UnaryFunction &f)
 Transform each vertex in the mesh with the specified function.

Detailed Description

template<int _N, int _M = _N, bool _A = true, typename T = double, typename V = ads::FixedArray<_N,T>, typename IS = Simplex<_M,int>>
class geom::IndSimpSet< _N, _M, _A, T, V, IS >

Class for a mesh that stores vertices and indexed simplices.

Parameters:
N is the space dimension.
M is the simplex dimension By default it is N.
A determines whether the mesh will allocate its own memory or use externally allocated memory. By default A is true.
T is the number type. By default it is double.
V is the vertex type, an N-tuple of the number type. It must be subscriptable. By default it is ads::FixedArray<N,T>.
IS is the indexed simplex type, a tuple of M+1 integers. It must be subscriptable. By default it is Simplex<M,int>.
Note that the indices for indexed simplices follow the C convention of starting at 0.

The free functions that operate on this class are grouped into the following categories:


Constructor & Destructor Documentation

template<int _N, int _M = _N, bool _A = true, typename T = double, typename V = ads::FixedArray<_N,T>, typename IS = Simplex<_M,int>>
template<bool A1, bool A2>
geom::IndSimpSet< _N, _M, _A, T, V, IS >::IndSimpSet ( const ads::Array< 1, Vertex, A1 > &  vertices,
const ads::Array< 1, IndexedSimplex, A2 > &  indexedSimplices 
) [inline]

Construct from arrays of vertices and indexed simplices.

If A is true, the arrays will be copied. If A is false, the memory is adopted. It will not be freed when the class destructor is called.

Parameters:
vertices is the array of vertices.
indexedSimplices is the array of indexed simplices.

template<int _N, int _M = _N, bool _A = true, typename T = double, typename V = ads::FixedArray<_N,T>, typename IS = Simplex<_M,int>>
geom::IndSimpSet< _N, _M, _A, T, V, IS >::IndSimpSet ( const SizeType  numVertices,
void *  vertices,
const SizeType  numSimplices,
void *  indexedSimplices 
) [inline]

Construct from pointers to the vertices and indexed simplices.

If A is true, the arrays will be copied. If A is false, the memory is adopted. It will not be freed when the class destructor is called. The objects to which vertices and indexedSimplices point will be cast to Vertex and IndexedSimplex, respectively. Thus they should have the same memory layout as these classes.

Parameters:
numVertices is the number of vertices.
vertices points to the data for the vertices.
numSimplices is the number of simplices.
indexedSimplices points to the data for the indexed simplices.

template<int _N, int _M = _N, bool _A = true, typename T = double, typename V = ads::FixedArray<_N,T>, typename IS = Simplex<_M,int>>
geom::IndSimpSet< _N, _M, _A, T, V, IS >::IndSimpSet ( const SizeType  numVertices,
const void *  vertices,
const SizeType  numSimplices,
const void *  indexedSimplices 
) [inline]

Construct from pointers to the vertices and indexed simplices.

A must be true to use this constructor.

The objects to which vertices and indexedSimplices point will be cast to Vertex and IndexedSimplex, respectively. Thus they should have the same memory layout as these classes.

Parameters:
numVertices is the number of vertices.
vertices points to the data for the vertices.
numSimplices is the number of simplices.
indexedSimplices points to the data for the indexed simplices.

template<int _N, int _M = _N, bool _A = true, typename T = double, typename V = ads::FixedArray<_N,T>, typename IS = Simplex<_M,int>>
geom::IndSimpSet< _N, _M, _A, T, V, IS >::IndSimpSet ( const SizeType  numVertices,
const SizeType  numSimplices 
) [inline]

Construct from the number of vertices and simplices.

A must be true.

The vertices and indexed simplices are left uninitialized.

Parameters:
numVertices is the number of vertices.
numSimplices is the number of simplices.


Member Function Documentation

template<int _N, int _M = _N, bool _A = true, typename T = double, typename V = ads::FixedArray<_N,T>, typename IS = Simplex<_M,int>>
void geom::IndSimpSet< _N, _M, _A, T, V, IS >::build ( const SizeType  numVertices,
const SizeType  numSimplices 
) [inline]

Build from the number of vertices and simplices.

Performs same actions as the constructor.

Parameters:
numVertices is the number of vertices.
numSimplices is the number of simplices.

template<int _N, int _M = _N, bool _A = true, typename T = double, typename V = ads::FixedArray<_N,T>, typename IS = Simplex<_M,int>>
void geom::IndSimpSet< _N, _M, _A, T, V, IS >::build ( const SizeType  numVertices,
const void *  vertices,
const SizeType  numSimplices,
const void *  indexedSimplices 
)

Build from pointers to the vertices and indexed simplices.

Performs same actions as the constructor.

Parameters:
numVertices is the number of vertices.
vertices points to the data for the vertices.
numSimplices is the number of simplices.
indexedSimplices points to the data for the indexed simplices.

template<int _N, int _M = _N, bool _A = true, typename T = double, typename V = ads::FixedArray<_N,T>, typename IS = Simplex<_M,int>>
void geom::IndSimpSet< _N, _M, _A, T, V, IS >::build ( const SizeType  numVertices,
void *  vertices,
const SizeType  numSimplices,
void *  indexedSimplices 
)

Build from pointers to the vertices and indexed simplices.

Performs same actions as the constructor.

Parameters:
numVertices is the number of vertices.
vertices points to the data for the vertices.
numSimplices is the number of simplices.
indexedSimplices points to the data for the indexed simplices.

template<int _N, int _M = _N, bool _A = true, typename T = double, typename V = ads::FixedArray<_N,T>, typename IS = Simplex<_M,int>>
template<bool A1, bool A2>
void geom::IndSimpSet< _N, _M, _A, T, V, IS >::build ( const ads::Array< 1, Vertex, A1 > &  vertices,
const ads::Array< 1, IndexedSimplex, A2 > &  indexedSimplices 
) [inline]

Build from arrays of vertices and indexed simplices.

Performs same actions as the constructor.

Parameters:
vertices is the array of vertices.
indexedSimplices is the array of indexed simplices.

template<int _N, int _M = _N, bool _A = true, typename T = double, typename V = ads::FixedArray<_N,T>, typename IS = Simplex<_M,int>>
virtual void geom::IndSimpSet< _N, _M, _A, T, V, IS >::updateTopology (  )  [inline, virtual]

Update the data structure following a change in the topology.

For this class, this function does nothing. For derived classes, it updates data structures that hold auxillary topological information.

Reimplemented in geom::IndSimpSetIncAdj< _N, _M, _A, T, V, IS >.


The documentation for this class was generated from the following files:
Generated on Fri Aug 24 12:56:01 2007 for Computational Geometry Package by  doxygen 1.4.7