DGtal  1.0.beta
Data Structures | Public Types | Static Public Attributes | Private Types | Private Member Functions | Private Attributes | Friends
DGtal::KhalimskySpaceND< dim, TInteger > Class Template Reference

#include <DGtal/topology/KhalimskySpaceND.h>

Inheritance diagram for DGtal::KhalimskySpaceND< dim, TInteger >:
[legend]
Collaboration diagram for DGtal::KhalimskySpaceND< dim, TInteger >:
[legend]

Data Structures

struct  CellMap
 
struct  SCellMap
 
struct  SurfelMap
 

Public Types

enum  Closure { CLOSED, OPEN, PERIODIC }
 
typedef TInteger Integer
 
typedef NumberTraits< Integer >::UnsignedVersion Size
 
typedef SpaceND< dim, IntegerSpace
 
typedef KhalimskySpaceND< dim, IntegerCellularGridSpace
 
typedef KhalimskyPreSpaceND< dim, IntegerPreCellularGridSpace
 
typedef KhalimskyCell< dim, IntegerCell
 
typedef KhalimskyPreCell< dim, IntegerPreCell
 
typedef SignedKhalimskyCell< dim, IntegerSCell
 
typedef SignedKhalimskyPreCell< dim, IntegerSPreCell
 
typedef SCell Surfel
 
typedef bool Sign
 
using DirIterator = typename PreCellularGridSpace::DirIterator
 
typedef PointVector< dim, IntegerPoint
 
typedef PointVector< dim, IntegerVector
 
template<typename CellType >
using AnyCellCollection = typename PreCellularGridSpace::template AnyCellCollection< CellType >
 
typedef AnyCellCollection< CellCells
 
typedef AnyCellCollection< SCellSCells
 
typedef std::set< CellCellSet
 
typedef std::set< SCellSCellSet
 
typedef std::set< SCellSurfelSet
 

Public Member Functions

Standard services
 ~KhalimskySpaceND ()
 
 KhalimskySpaceND ()
 
 KhalimskySpaceND (const KhalimskySpaceND &other)=default
 
KhalimskySpaceNDoperator= (const KhalimskySpaceND &other)=default
 
 KhalimskySpaceND (KhalimskySpaceND &&other)=default
 
KhalimskySpaceNDoperator= (KhalimskySpaceND &&other)=default
 
bool init (const Point &lower, const Point &upper, bool isClosed)
 
bool init (const Point &lower, const Point &upper, Closure closure)
 
bool init (const Point &lower, const Point &upper, const std::array< Closure, dim > &closure)
 
Basic services
Size size (Dimension k) const
 
Integer min (Dimension k) const
 
Integer max (Dimension k) const
 
const PointlowerBound () const
 
const PointupperBound () const
 
const CelllowerCell () const
 
const CellupperCell () const
 
bool uIsValid (const PreCell &c, Dimension k) const
 
bool uIsValid (const PreCell &c) const
 
bool sIsValid (const SPreCell &c, Dimension k) const
 
bool sIsValid (const SPreCell &c) const
 
bool cIsValid (const Point &p, Dimension k) const
 
bool cIsValid (const Point &p) const
 
Closure type query
bool isSpaceClosed () const
 
bool isSpaceClosed (Dimension k) const
 
bool isSpacePeriodic () const
 
bool isSpacePeriodic (Dimension k) const
 
bool isAnyDimensionPeriodic () const
 
Closure getClosure (Dimension k) const
 
Cell creation services
Cell uCell (const PreCell &c) const
 
Cell uCell (const Point &kp) const
 
Cell uCell (Point p, const PreCell &c) const
 
SCell sCell (const SPreCell &c) const
 
SCell sCell (const Point &kp, Sign sign=POS) const
 
SCell sCell (Point p, const SPreCell &c) const
 
Cell uSpel (Point p) const
 
SCell sSpel (Point p, Sign sign=POS) const
 
Cell uPointel (Point p) const
 
SCell sPointel (Point p, Sign sign=POS) const
 
Read accessors to cells
Integer uKCoord (const Cell &c, Dimension k) const
 
Integer uCoord (const Cell &c, Dimension k) const
 
const PointuKCoords (const Cell &c) const
 
Point uCoords (const Cell &c) const
 
Integer sKCoord (const SCell &c, Dimension k) const
 
Integer sCoord (const SCell &c, Dimension k) const
 
const PointsKCoords (const SCell &c) const
 
Point sCoords (const SCell &c) const
 
Sign sSign (const SCell &c) const
 
Write accessors to cells
void uSetKCoord (Cell &c, Dimension k, Integer i) const
 
void sSetKCoord (SCell &c, Dimension k, Integer i) const
 
void uSetCoord (Cell &c, Dimension k, Integer i) const
 
void sSetCoord (SCell &c, Dimension k, Integer i) const
 
void uSetKCoords (Cell &c, const Point &kp) const
 
void sSetKCoords (SCell &c, const Point &kp) const
 
void uSetCoords (Cell &c, const Point &kp) const
 
void sSetCoords (SCell &c, const Point &kp) const
 
void sSetSign (SCell &c, Sign s) const
 
Conversion signed/unsigned
SCell signs (const Cell &p, Sign s) const
 
Cell unsigns (const SCell &p) const
 
SCell sOpp (const SCell &p) const
 
Cell topology services
Integer uTopology (const Cell &p) const
 
Integer sTopology (const SCell &p) const
 
Dimension uDim (const Cell &p) const
 
Dimension sDim (const SCell &p) const
 
bool uIsSurfel (const Cell &b) const
 
bool sIsSurfel (const SCell &b) const
 
bool uIsOpen (const Cell &p, Dimension k) const
 
bool sIsOpen (const SCell &p, Dimension k) const
 
Iterator services for cells
DirIterator uDirs (const Cell &p) const
 
DirIterator sDirs (const SCell &p) const
 
DirIterator uOrthDirs (const Cell &p) const
 
DirIterator sOrthDirs (const SCell &p) const
 
Dimension uOrthDir (const Cell &s) const
 
Dimension sOrthDir (const SCell &s) const
 
Unsigned cell geometry services
Integer uFirst (const PreCell &p, Dimension k) const
 
Cell uFirst (const PreCell &p) const
 
Integer uLast (const PreCell &p, Dimension k) const
 
Cell uLast (const PreCell &p) const
 
Cell uGetIncr (const Cell &p, Dimension k) const
 
bool uIsMax (const Cell &p, Dimension k) const
 
bool uIsInside (const PreCell &p, Dimension k) const
 
bool uIsInside (const PreCell &p) const
 
bool cIsInside (const Point &p, Dimension k) const
 
bool cIsInside (const Point &p) const
 
Cell uGetMax (Cell p, Dimension k) const
 
Cell uGetDecr (const Cell &p, Dimension k) const
 
bool uIsMin (const Cell &p, Dimension k) const
 
Cell uGetMin (Cell p, Dimension k) const
 
Cell uGetAdd (const Cell &p, Dimension k, Integer x) const
 
Cell uGetSub (const Cell &p, Dimension k, Integer x) const
 
Integer uDistanceToMax (const Cell &p, Dimension k) const
 
Integer uDistanceToMin (const Cell &p, Dimension k) const
 
Cell uTranslation (const Cell &p, const Vector &vec) const
 
Cell uProjection (const Cell &p, const Cell &bound, Dimension k) const
 
void uProject (Cell &p, const Cell &bound, Dimension k) const
 
bool uNext (Cell &p, const Cell &lower, const Cell &upper) const
 
Signed cell geometry services
Integer sFirst (const SPreCell &p, Dimension k) const
 
SCell sFirst (const SPreCell &p) const
 
Integer sLast (const SPreCell &p, Dimension k) const
 
SCell sLast (const SPreCell &p) const
 
SCell sGetIncr (const SCell &p, Dimension k) const
 
bool sIsMax (const SCell &p, Dimension k) const
 
bool sIsInside (const SPreCell &p, Dimension k) const
 
bool sIsInside (const SPreCell &p) const
 
SCell sGetMax (SCell p, Dimension k) const
 
SCell sGetDecr (const SCell &p, Dimension k) const
 
bool sIsMin (const SCell &p, Dimension k) const
 
SCell sGetMin (SCell p, Dimension k) const
 
SCell sGetAdd (const SCell &p, Dimension k, Integer x) const
 
SCell sGetSub (const SCell &p, Dimension k, Integer x) const
 
Integer sDistanceToMax (const SCell &p, Dimension k) const
 
Integer sDistanceToMin (const SCell &p, Dimension k) const
 
SCell sTranslation (const SCell &p, const Vector &vec) const
 
SCell sProjection (const SCell &p, const SCell &bound, Dimension k) const
 
void sProject (SCell &p, const SCell &bound, Dimension k) const
 
bool sNext (SCell &p, const SCell &lower, const SCell &upper) const
 
Neighborhood services
Cells uNeighborhood (const Cell &cell) const
 
SCells sNeighborhood (const SCell &cell) const
 
Cells uProperNeighborhood (const Cell &cell) const
 
SCells sProperNeighborhood (const SCell &cell) const
 
Cell uAdjacent (const Cell &p, Dimension k, bool up) const
 
SCell sAdjacent (const SCell &p, Dimension k, bool up) const
 
Incidence services
Cell uIncident (const Cell &c, Dimension k, bool up) const
 
SCell sIncident (const SCell &c, Dimension k, bool up) const
 
Cells uLowerIncident (const Cell &c) const
 
Cells uUpperIncident (const Cell &c) const
 
SCells sLowerIncident (const SCell &c) const
 
SCells sUpperIncident (const SCell &c) const
 
Cells uFaces (const Cell &c) const
 
Cells uCoFaces (const Cell &c) const
 
bool sDirect (const SCell &p, Dimension k) const
 
SCell sDirectIncident (const SCell &p, Dimension k) const
 
SCell sIndirectIncident (const SCell &p, Dimension k) const
 
DGtal interface
void selfDisplay (std::ostream &out) const
 
bool isValid () const
 

Static Public Attributes

static const constexpr Dimension dimension = dim
 
static const constexpr Dimension DIM = dim
 
static const constexpr Sign POS = true
 
static const constexpr Sign NEG = false
 

Private Types

typedef KhalimskySpaceNDHelper< KhalimskySpaceND< dim, TInteger > > Helper
 

Private Member Functions

 BOOST_CONCEPT_ASSERT ((concepts::CInteger< TInteger >))
 
Internals
void uAddFaces (Cells &faces, const Cell &c, Dimension axis) const
 
void uAddCoFaces (Cells &cofaces, const Cell &c, Dimension axis) const
 

Private Attributes

Point myLower
 
Point myUpper
 
Cell myCellLower
 
Cell myCellUpper
 
std::array< Closure, dimensionmyClosure
 

Friends

class KhalimskySpaceNDHelper< KhalimskySpaceND< dim, TInteger > >
 

Detailed Description

template<Dimension dim, typename TInteger>
class DGtal::KhalimskySpaceND< dim, TInteger >

Aim: This class is a model of CCellularGridSpaceND. It represents the cubical grid as a cell complex, whose cells are defined as an array of integers. The topology of the cells is defined by the parity of the coordinates (even: closed, odd: open).

Description of template class 'KhalimskySpaceND'

When initializing the space using init(), the user should choose, for each dimension spanned by the space, between a closed and non-periodic (default) cell dimension, an open cell dimension or a periodic cell dimension. The space is generally finite, except for arbitrary size integers and when the space has a periodic dimension.

Supposing that the space has been initialized with digital bounds lower and upper, the methods lowerBound() and upperBound() will always return, respectively, lower and upper. It as also true for periodic dimension, in order to span over the unique digital points of the space.

In the same way, lowerCell() and upperCell() respect the following rules:

The special behavior for periodic dimensions guarantees that each cell has unique Khalimsky coordinates in this range. It is useful to span the space and also for cell-based containers (see e.g. CubicalComplex). Uniqueness also gives meaning to equality tests between cells.

Following this concept, the related methods size(), min(), max(), uFirst(), uLast(), uGetMin(), uGetMax(), uDistanceToMin(), uDistanceToMax(), sFirst(), sLast(), sGetMin(), sGetMax(), sDistanceToMin() and sDistanceToMax() behave for periodic dimensions like for finite dimensions, using the bounds described above.

Thus, if a cell needs to be compared to the bounds, prefer using dedicated tests like uIsMin(), uIsMax(), sIsMin() and sIsMax() that return always false for a periodic dimension, and uIsInside() and sIsInside() that return always true for a periodic dimension.

To be consistent with those choices, each cell returned or modified by a KhalimskySpaceND method will have his Khalimsky coordinates along periodic dimensions between the corresponding coordinates of lowerCell() and upperCell(). But, in order to keep low computational cost, each cell passed by parameter to a KhalimskySpaceND method must follow the same conditions. This validity can be tested with the dedicated methods uIsValid() and sIsValid().

Exceptions exist for uCell(const PreCell &) const and sCell(const SPreCell &) const that are specially featured to correct Khalimsky coordinates of a given cell. In addition, when a method accepts a coordinate as parameter, it is always corrected along periodic dimensions.

Template Parameters
dimthe dimension of the digital space.
TIntegerthe Integer class used to specify the arithmetic computations (default type = int32).
Note
Essentially a backport from ImaGene.
Warning
Periodic Khalimsky space and per-dimension closure specification are new features. Therefore, there is no guarantee that it is compatible with the whole DGtal library.
Examples:
dec/exampleDECSurface.cpp, dec/exampleHeatLaplace.cpp, examples/tutorial-examples/polyhedralizer.cpp, geometry/curves/estimation/exampleCurvature.cpp, geometry/curves/exampleGridCurve3d-2.cpp, geometry/surfaces/dvcm-2d-curvature.cpp, geometry/surfaces/dvcm-3d.cpp, geometry/surfaces/greedy-plane-segmentation-ex2.cpp, geometry/surfaces/greedy-plane-segmentation.cpp, geometry/tools/exampleAlphaShape.cpp, geometry/volumes/distance/exampleFMM3D.cpp, graph/volDistanceTraversal.cpp, io/boards/dgtalBoard3D-2-ks.cpp, io/boards/dgtalBoard3DTo2D-KSCell.cpp, io/viewDualSurface.cpp, io/viewers/viewer3D-10-interaction.cpp, io/viewers/viewer3D-11-extension.cpp, io/viewers/viewer3D-4bis-illustrationMode.cpp, shapes/viewMarchingCubes.cpp, shapes/viewPolygonalMarchingCubes.cpp, topology/3dKSSurfaceExtraction.cpp, topology/area-estimation-with-digital-surface.cpp, topology/area-estimation-with-indexed-digital-surface.cpp, topology/ctopo-1-3d.cpp, topology/ctopo-1.cpp, topology/ctopo-1s-3d.cpp, topology/ctopo-2-3d.cpp, topology/ctopo-2.cpp, topology/ctopo-fillContours.cpp, topology/cubical-complex-collapse.cpp, topology/cubical-complex-illustrations.cpp, topology/digitalSetToCubicalComplexes2D.cpp, topology/digitalSurfaceSlice.cpp, topology/frontierAndBoundary.cpp, topology/khalimskySpaceScanner.cpp, topology/trackImplicitPolynomialSurfaceToOFF.cpp, topology/volBreadthFirstTraversal.cpp, topology/volMarchingCubes.cpp, topology/volScanBoundary.cpp, topology/volToOFF.cpp, topology/volTrackBoundary.cpp, tutorial-examples/AreaSurfaceEstimation-final.cpp, tutorial-examples/FMMErosion.cpp, tutorial-examples/freemanChainFromImage.cpp, and tutorial-examples/polyhedralizer.cpp.

Definition at line 64 of file KhalimskySpaceND.h.

Member Typedef Documentation

template<Dimension dim, typename TInteger>
template<typename CellType >
using DGtal::KhalimskySpaceND< dim, TInteger >::AnyCellCollection = typename PreCellularGridSpace::template AnyCellCollection< CellType >

Definition at line 435 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef KhalimskyCell< dim, Integer > DGtal::KhalimskySpaceND< dim, TInteger >::Cell

Definition at line 414 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef AnyCellCollection<Cell> DGtal::KhalimskySpaceND< dim, TInteger >::Cells

Definition at line 438 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef std::set<Cell> DGtal::KhalimskySpaceND< dim, TInteger >::CellSet

Preferred type for defining a set of Cell(s).

Definition at line 443 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef KhalimskySpaceND<dim, Integer> DGtal::KhalimskySpaceND< dim, TInteger >::CellularGridSpace

Definition at line 410 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
using DGtal::KhalimskySpaceND< dim, TInteger >::DirIterator = typename PreCellularGridSpace::DirIterator

Definition at line 421 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef KhalimskySpaceNDHelper< KhalimskySpaceND< dim, TInteger > > DGtal::KhalimskySpaceND< dim, TInteger >::Helper
private

Features basic operations on coordinates, especially for periodic dimensions.

Definition at line 395 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef TInteger DGtal::KhalimskySpaceND< dim, TInteger >::Integer

Arithmetic ring induced by (+,-,*) and Integer numbers.

Definition at line 403 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef PointVector< dim, Integer > DGtal::KhalimskySpaceND< dim, TInteger >::Point

Definition at line 424 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef KhalimskyPreCell< dim, Integer > DGtal::KhalimskySpaceND< dim, TInteger >::PreCell

Definition at line 415 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef KhalimskyPreSpaceND<dim, Integer> DGtal::KhalimskySpaceND< dim, TInteger >::PreCellularGridSpace

Definition at line 411 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef SignedKhalimskyCell< dim, Integer > DGtal::KhalimskySpaceND< dim, TInteger >::SCell

Definition at line 416 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef AnyCellCollection<SCell> DGtal::KhalimskySpaceND< dim, TInteger >::SCells

Definition at line 439 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef std::set<SCell> DGtal::KhalimskySpaceND< dim, TInteger >::SCellSet

Preferred type for defining a set of SCell(s).

Definition at line 446 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef bool DGtal::KhalimskySpaceND< dim, TInteger >::Sign

Definition at line 420 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef NumberTraits<Integer>::UnsignedVersion DGtal::KhalimskySpaceND< dim, TInteger >::Size

Type used to represent sizes in the digital space.

Definition at line 406 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef SpaceND<dim, Integer> DGtal::KhalimskySpaceND< dim, TInteger >::Space

Definition at line 409 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef SignedKhalimskyPreCell< dim, Integer > DGtal::KhalimskySpaceND< dim, TInteger >::SPreCell

Definition at line 417 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef SCell DGtal::KhalimskySpaceND< dim, TInteger >::Surfel

Definition at line 419 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef std::set<SCell> DGtal::KhalimskySpaceND< dim, TInteger >::SurfelSet

Preferred type for defining a set of surfels (always signed cells).

Definition at line 449 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
typedef PointVector< dim, Integer > DGtal::KhalimskySpaceND< dim, TInteger >::Vector

Definition at line 425 of file KhalimskySpaceND.h.

Member Enumeration Documentation

template<Dimension dim, typename TInteger>
enum DGtal::KhalimskySpaceND::Closure

Boundaries closure type.

Enumerator
CLOSED 

The dimension is closed and non-periodic.

OPEN 

The dimension is open.

PERIODIC 

The dimension is periodic.

Definition at line 470 of file KhalimskySpaceND.h.

471  {
472  CLOSED,
473  OPEN,
474  PERIODIC
475  };
The dimension is open.
The dimension is periodic.
The dimension is closed and non-periodic.

Constructor & Destructor Documentation

template<Dimension dim, typename TInteger>
DGtal::KhalimskySpaceND< dim, TInteger >::~KhalimskySpaceND ( )

Destructor.

template<Dimension dim, typename TInteger>
DGtal::KhalimskySpaceND< dim, TInteger >::KhalimskySpaceND ( )

Default constructor.

template<Dimension dim, typename TInteger>
DGtal::KhalimskySpaceND< dim, TInteger >::KhalimskySpaceND ( const KhalimskySpaceND< dim, TInteger > &  other)
default

Copy constructor.

Parameters
otherthe object to clone.
template<Dimension dim, typename TInteger>
DGtal::KhalimskySpaceND< dim, TInteger >::KhalimskySpaceND ( KhalimskySpaceND< dim, TInteger > &&  other)
default

Move constructor.

Parameters
otherthe object to clone.

Member Function Documentation

template<Dimension dim, typename TInteger>
DGtal::KhalimskySpaceND< dim, TInteger >::BOOST_CONCEPT_ASSERT ( (concepts::CInteger< TInteger >)  )
private
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::cIsInside ( const Point p,
Dimension  k 
) const

Useful to check if you are going out of the space.

Parameters
pany integer point (Khalimsky coordinates).
kthe tested coordinate.
Returns
true if [p] has its [k]-coordinate within the allowed bounds.
Note
It returns always true for periodic dimension.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::cIsInside ( const Point p) const

Useful to check if you are going out of the space.

Parameters
pany integer point (Khalimsky coordinates).
Returns
true if [p] has its coordinates within the allowed bounds.
Note
Only the non-periodic dimensions are checked.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::cIsValid ( const Point p,
Dimension  k 
) const
Parameters
pan integer point (Khalimsky coordinates of cell).
ka dimension.
Returns
true if the given cell has his k-th Khalimsky coordinate between those of the cells returned by lowerCell and upperCell.
Note
For periodic dimension, even if there is no bounds, it guarantees that each cell has unique coordinates.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::cIsValid ( const Point p) const
Parameters
pan integer point (Khalimsky coordinates of cell).
Returns
true if the given cell has Khalimsky coordinates between those of the cells returned by lowerCell and upperCell.
Note
For periodic dimension, even if there is no bounds, it guarantees that each cell has unique coordinates.
template<Dimension dim, typename TInteger>
Closure DGtal::KhalimskySpaceND< dim, TInteger >::getClosure ( Dimension  k) const

Gets closure type.

Parameters
kthe dimension.
Returns
closure type along the specified dimension.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::init ( const Point lower,
const Point upper,
bool  isClosed 
)

Specifies the upper and lower bounds for the maximal cells in this space.

Parameters
lowerthe lowest point in this space (digital coords)
upperthe upper point in this space (digital coords)
isClosed'true' if this space is closed and non-periodic in every dimension, 'false' if open.
Returns
true if the initialization was valid (ie, such bounds are representable with these integers).
Examples:
dec/exampleDECSurface.cpp, dec/exampleHeatLaplace.cpp, examples/tutorial-examples/polyhedralizer.cpp, geometry/curves/estimation/exampleCurvature.cpp, geometry/curves/exampleGridCurve2d.cpp, geometry/curves/exampleGridCurve3d.cpp, geometry/surfaces/dvcm-2d-curvature.cpp, geometry/surfaces/dvcm-3d.cpp, geometry/surfaces/greedy-plane-segmentation-ex2.cpp, geometry/surfaces/greedy-plane-segmentation.cpp, geometry/tools/exampleAlphaShape.cpp, geometry/volumes/distance/exampleFMM3D.cpp, graph/volDistanceTraversal.cpp, io/boards/dgtalBoard3D-2-ks.cpp, io/boards/dgtalBoard3DTo2D-KSCell.cpp, io/viewDualSurface.cpp, io/viewers/viewer3D-10-interaction.cpp, io/viewers/viewer3D-11-extension.cpp, io/viewers/viewer3D-4bis-illustrationMode.cpp, shapes/viewMarchingCubes.cpp, shapes/viewPolygonalMarchingCubes.cpp, topology/3dKSSurfaceExtraction.cpp, topology/area-estimation-with-digital-surface.cpp, topology/area-estimation-with-indexed-digital-surface.cpp, topology/ctopo-1-3d.cpp, topology/ctopo-1.cpp, topology/ctopo-1s-3d.cpp, topology/ctopo-2-3d.cpp, topology/ctopo-2.cpp, topology/cubical-complex-collapse.cpp, topology/cubical-complex-illustrations.cpp, topology/digitalSetToCubicalComplexes2D.cpp, topology/digitalSurfaceSlice.cpp, topology/frontierAndBoundary.cpp, topology/khalimskySpaceScanner.cpp, topology/trackImplicitPolynomialSurfaceToOFF.cpp, topology/volBreadthFirstTraversal.cpp, topology/volMarchingCubes.cpp, topology/volScanBoundary.cpp, topology/volToOFF.cpp, topology/volTrackBoundary.cpp, tutorial-examples/AreaSurfaceEstimation-final.cpp, tutorial-examples/FMMErosion.cpp, tutorial-examples/freemanChainFromImage.cpp, and tutorial-examples/polyhedralizer.cpp.

Referenced by DGtal::functions::generateVoxelComplexTable().

template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::init ( const Point lower,
const Point upper,
Closure  closure 
)

Specifies the upper and lower bounds for the maximal cells in this space.

Parameters
lowerthe lowest point in this space (digital coords)
upperthe upper point in this space (digital coords)
closureCLOSED, OPEN or PERIODIC if this space is resp. closed (and non-periodic), open or periodic in every dimension.
Returns
true if the initialization was valid (ie, such bounds are representable with these integers).
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::init ( const Point lower,
const Point upper,
const std::array< Closure, dim > &  closure 
)

Specifies the upper and lower bounds for the maximal cells in this space.

Parameters
lowerthe lowest point in this space (digital coords)
upperthe upper point in this space (digital coords)
closurean array of CLOSED, OPEN or PERIODIC if this space is resp. closed (and non-periodic), open or periodic in the corresponding dimension.
Returns
true if the initialization was valid (ie, such bounds are representable with these integers).
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::isAnyDimensionPeriodic ( ) const
Returns
'true' iff the space is periodic along at least one dimension.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::isSpaceClosed ( ) const
Returns
'true' iff the space is closed or periodic along every dimension.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::isSpaceClosed ( Dimension  k) const
Parameters
kthe dimension.
Returns
'true' iff the space is closed or periodic along the specified dimension.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::isSpacePeriodic ( ) const
Returns
'true' iff the space is periodic along every dimension.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::isSpacePeriodic ( Dimension  k) const
Parameters
kthe dimension.
Returns
'true' iff the space is periodic along the specified dimension.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::isValid ( ) const

Checks the validity/consistency of the object.

Returns
'true' if the object is valid, 'false' otherwise.
template<Dimension dim, typename TInteger>
const Point& DGtal::KhalimskySpaceND< dim, TInteger >::lowerBound ( ) const
Returns
the lower bound for digital points in this space.
Note
For periodic dimension, it returns the lower digital point with unique coordinates. See also the class documentation.
Examples:
shapes/viewMarchingCubes.cpp, shapes/viewPolygonalMarchingCubes.cpp, topology/ctopo-2.cpp, topology/digitalSurfaceSlice.cpp, topology/volMarchingCubes.cpp, topology/volScanBoundary.cpp, and topology/volToOFF.cpp.
template<Dimension dim, typename TInteger>
const Cell& DGtal::KhalimskySpaceND< dim, TInteger >::lowerCell ( ) const
Returns
the lower bound for cells in this space.
Note
For periodic dimension, it returns the upper cell with unique coordinates. See also the class documentation.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::max ( Dimension  k) const
Parameters
ka coordinate.
Returns
the maximal digital coordinate in the k-dimension.
Note
For periodic dimension, it returns the maximal unique digital coordinate along that dimension.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::min ( Dimension  k) const
Parameters
ka dimension.
Returns
the minimal digital coordinate in the k-dimension.
Note
For periodic dimension, it returns the minimal unique digital coordinate along that dimension.
template<Dimension dim, typename TInteger>
KhalimskySpaceND& DGtal::KhalimskySpaceND< dim, TInteger >::operator= ( const KhalimskySpaceND< dim, TInteger > &  other)
default

Copy operator.

Parameters
otherthe object to copy.
Returns
a reference on 'this'.
template<Dimension dim, typename TInteger>
KhalimskySpaceND& DGtal::KhalimskySpaceND< dim, TInteger >::operator= ( KhalimskySpaceND< dim, TInteger > &&  other)
default

Move operator.

Parameters
otherthe object to copy.
Returns
a reference on 'this'.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sAdjacent ( const SCell p,
Dimension  k,
bool  up 
) const
Parameters
pany cell.
kthe coordinate that is changed.
upif 'true' the orientation is forward along axis [k], otherwise backward.
Returns
the adjacent element to [p] along axis [k] in the given direction and orientation.
Note
It is an alias to 'up ? sGetIncr( p, k ) : sGetDecr( p, k )'.
Precondition
sIsValid(p) and not sIs(Min|Max)(p, k) depending on up value.
Postcondition
`sIsValid(sAdjacent(p, k, up)) is true.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sCell ( const SPreCell c) const

From a signed cell, returns a signed cell lying into this Khalismky space.

Along a non-periodic dimension, if the given Khalimsky coordinate lies outside the space, it replaces it by the nearest valid coordinate.

Along a periodic dimension, the Khalimsky coordinate is corrected (by periodicity) to lie between the coordinates of lowerCell() and upperCell().

Parameters
ca signed cell.
Returns
the same signed cell with appropriate Khalimsky coordinates along periodic dimensions.
Postcondition
sIsValid(sCell(c)) is true.
sCell(c) == c if sIsValid(c).
Examples:
io/boards/dgtalBoard3DTo2D-KSCell.cpp, io/viewers/viewer3D-10-interaction.cpp, and topology/ctopo-1s-3d.cpp.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sCell ( const Point kp,
Sign  sign = POS 
) const

From the Khalimsky coordinates of a cell and a sign, builds the corresponding signed cell lying into this Khalismky space.

Along a non-periodic dimension, if the given Khalimsky coordinate lies outside the space, it is replaced by the nearest valid coordinate.

Along a periodic dimension, the Khalimsky coordinate is corrected (by periodicity) to lie between the coordinates of lowerCell() and upperCell().

Parameters
kpan integer point (Khalimsky coordinates of cell).
signthe sign of the cell (either POS or NEG).
Returns
the signed cell.
Postcondition
sIsValid(sCell(kp, sign)) is true.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sCell ( Point  p,
const SPreCell c 
) const

From the digital coordinates of a point in Zn and a signed cell type, builds the corresponding signed cell lying into this Khalismky space.

Along a non-periodic dimension, if the given digital coordinate lies outside the space, it is replaced by the nearest valid coordinate.

Along a periodic dimension, the digital coordinate is corrected (by periodicity) to lie between the coordinates of lowerCell() and upperCell().

Parameters
pan integer point (digital coordinates of cell).
canother cell defining the topology and sign.
Returns
the cell having the topology and sign of [c] and the given digital coordinates [p].
Postcondition
sIsValid(sCell(p, c)) is true.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::sCoord ( const SCell c,
Dimension  k 
) const
Parameters
cany signed cell.
kany valid dimension.
Returns
its digital coordinate along [k].
Precondition
uIsValid(c) is true.
template<Dimension dim, typename TInteger>
Point DGtal::KhalimskySpaceND< dim, TInteger >::sCoords ( const SCell c) const
template<Dimension dim, typename TInteger>
Dimension DGtal::KhalimskySpaceND< dim, TInteger >::sDim ( const SCell p) const
Parameters
pany signed cell.
Returns
the dimension of the cell [p].

Referenced by DGtal::functors::SCellToPoint< KSpace >::operator()().

template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::sDirect ( const SCell p,
Dimension  k 
) const

Return 'true' if the direct orientation of [p] along [k] is in the positive coordinate direction. The direct orientation in a direction allows to go from positive incident cells to positive incident cells. This means that

K.sSign( K.sIncident( p, k, K.sDirect( p, k ) ) ) == K.POS

is always true.

Parameters
pany signed cell.
kany coordinate.
Returns
the direct orientation of [p] along [k] (true is upward, false is backward).
Examples:
topology/area-estimation-with-digital-surface.cpp, and topology/area-estimation-with-indexed-digital-surface.cpp.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sDirectIncident ( const SCell p,
Dimension  k 
) const
Parameters
pany signed cell.
kany coordinate.
Returns
the direct incident cell of [p] along [k] (the incident cell along [k])
Precondition
sIsValid(p) and the cell should have a direct incident cell in this direction.
Postcondition
sIsValid(c) is true for every returned cell c. ose sign is positive).
Examples:
examples/tutorial-examples/polyhedralizer.cpp, geometry/surfaces/greedy-plane-segmentation-ex2.cpp, geometry/surfaces/greedy-plane-segmentation.cpp, tutorial-examples/AreaSurfaceEstimation-final.cpp, and tutorial-examples/polyhedralizer.cpp.

Referenced by DGtal::functors::SCellToInnerPoint< KSpace >::operator()(), and DGtal::functors::SCellToIncidentPoints< KSpace >::operator()().

template<Dimension dim, typename TInteger>
DirIterator DGtal::KhalimskySpaceND< dim, TInteger >::sDirs ( const SCell p) const

Given a signed cell [p], returns an iterator to iterate over each coordinate the cell spans. (A spel spans all coordinates; a surfel all but one, etc). Example:

...
for ( KSpace::DirIterator q = ks.uDirs( p ); q != 0; ++q )
{
Dimension dir = *q;
...
}
Parameters
pany signed cell.
Returns
an iterator that points on the first coordinate spanned by the cell.
Examples:
geometry/curves/exampleGridCurve3d-2.cpp, and topology/digitalSurfaceSlice.cpp.

Referenced by DGtal::functors::SCellToPoint< KSpace >::operator()(), DGtal::functors::SCellToArrow< KSpace >::operator()(), and DGtal::functors::SCellToCode< KSpace >::operator()().

template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::sDistanceToMax ( const SCell p,
Dimension  k 
) const

Useful to check if you are going out of the space (for non-periodic dimensions).

Parameters
pany cell.
kthe coordinate that is tested.
Returns
the number of increment to do to reach the maximum value.
Precondition
sIsValid(p) is true.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::sDistanceToMin ( const SCell p,
Dimension  k 
) const

Useful to check if you are going out of the space (for non-periodic dimensions).

Parameters
pany cell.
kthe coordinate that is tested.
Returns
the number of decrement to do to reach the minimum value.
Precondition
sIsValid(p) is true.
template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::selfDisplay ( std::ostream &  out) const

Writes/Displays the object on an output stream.

Parameters
outthe output stream where the object is written.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::sFirst ( const SPreCell p,
Dimension  k 
) const
Returns
the k-th coordinate of the first cell of the space with the same type as [p].
Note
For periodic dimension, it returns the first unique coordinate of a cell of same type as p.
Postcondition
The returned coordinate is between lowerCell()[k] and upperCell()[k].
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sFirst ( const SPreCell p) const
Returns
the first cell of the space with the same type as [p].
Note
Along periodic dimensions, it returns the first unique coordinate of a cell of same type as p.
Postcondition
sIsValid(sFirst(p)) is true.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sGetAdd ( const SCell p,
Dimension  k,
Integer  x 
) const
Parameters
pany cell.
kthe coordinate that is changed.
xthe increment.
Returns
the same element as [p] except for a coordinate [k] incremented with x.
Precondition
sIsValid(p) and ( x <= sDistanceToMax(p, k) or isSpacePeriodic(k) ).
Postcondition
sIsValid(sGetAdd(p, k, x)) is true.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sGetDecr ( const SCell p,
Dimension  k 
) const
Parameters
pany cell.
kthe coordinate that is changed.
Returns
the same element as [p] except for an decremented coordinate [k].
Precondition
sIsValid(p) and not sIsMin(p).
Postcondition
sIsValid(sGetDecr(p, k)) is true.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sGetIncr ( const SCell p,
Dimension  k 
) const
Parameters
pany cell.
kthe coordinate that is changed.
Returns
the same element as [p] except for the incremented coordinate [k].
Precondition
sIsValid(p) and not sIsMax(p).
Postcondition
sIsValid(sGetIncr(p, k)) is true.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sGetMax ( SCell  p,
Dimension  k 
) const

Useful to check if you are going out of the space.

Parameters
pany cell.
kthe concerned coordinate.
Returns
the cell similar to [p] but with the maximum allowed [k]-coordinate.
Precondition
sIsValid(p, d) is true for each dimension d different than k.
Postcondition
sIsValid(sGetMax(p, k)) is true.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sGetMin ( SCell  p,
Dimension  k 
) const

Useful to check if you are going out of the space.

Parameters
pany cell.
kthe concerned coordinate.
Returns
the cell similar to [p] but with the minimum allowed [k]-coordinate.
Precondition
sIsValid(p, d) is true for each dimension d different than k.
Postcondition
sIsValid(sGetMin(p, k)) is true.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sGetSub ( const SCell p,
Dimension  k,
Integer  x 
) const
Parameters
pany cell.
kthe coordinate that is changed.
xthe decrement.
Returns
the same element as [p] except for a coordinate [k] decremented with x.
Precondition
suIsValid(p)and (x <= sDistanceToMin(p, k)orisSpacePeriodic(k)). @postsIsValid(sGetSub(p, k, x))` is true.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::signs ( const Cell p,
Sign  s 
) const

Creates a signed cell from an unsigned one and a given sign.

Parameters
pany unsigned cell.
sa sign.
Returns
the signed version of the cell [p] with sign [s].
Precondition
uIsValid(p) is true.
Postcondition
sIsValid(signs(p, s)) is true.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sIncident ( const SCell c,
Dimension  k,
bool  up 
) const
Parameters
cany signed cell.
kany coordinate.
upif 'true' the orientation is forward along axis [k], otherwise backward.
Returns
the forward or backward signed cell incident to [c] along axis [k], depending on [up]. It is worthy to note that the forward and backward cell have opposite sign. Furthermore, the sign of these cells is defined so as to satisfy a boundary operator.
Note
It may be a lower incident cell if [c] is open along axis [k], else an upper incident cell.
Precondition
sIsValid(c) and the cell should have an incident cell in this direction/orientation.
Postcondition
sIsValid(sIncident(c, k, up)) is true.
Examples:
io/boards/dgtalBoard3D-2-ks.cpp, io/viewers/viewer3D-4bis-illustrationMode.cpp, and topology/frontierAndBoundary.cpp.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sIndirectIncident ( const SCell p,
Dimension  k 
) const
Parameters
pany signed cell.
kany coordinate.
Returns
the indirect incident cell of [p] along [k] (the incident cell along [k] whose sign is negative).
Precondition
sIsValid(p) and the cell should have an indirect incident cell in this direction.
Postcondition
sIsValid(c) is true for every returned cell c. ose sign is positive).

Referenced by DGtal::functors::SCellToPoint< KSpace >::operator()(), DGtal::functors::SCellToArrow< KSpace >::operator()(), DGtal::functors::SCellToOuterPoint< KSpace >::operator()(), DGtal::functors::SCellToIncidentPoints< KSpace >::operator()(), and DGtal::functors::SCellToCode< KSpace >::operator()().

template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::sIsInside ( const SPreCell p,
Dimension  k 
) const

Useful to check if you are going out of the space.

Parameters
pany cell.
kthe tested coordinate.
Returns
true if [p] has its [k]-coordinate within the allowed bounds.
Note
It returns always true for periodic dimension.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::sIsInside ( const SPreCell p) const

Useful to check if you are going out of the space.

Parameters
pany cell.
Returns
true if [p] has its coordinates within the allowed bounds.
Note
Only the non-periodic dimensions are checked.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::sIsMax ( const SCell p,
Dimension  k 
) const

Useful to check if you are going out of the space.

Parameters
pany cell.
kthe tested coordinate.
Returns
true if [p] cannot have its [k]-coordinate augmented without leaving the space.
Note
It returns always false for periodic dimension.
Precondition
sIsInside(p) is true.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::sIsMin ( const SCell p,
Dimension  k 
) const

Useful to check if you are going out of the space.

Parameters
pany cell.
kthe tested coordinate.
Returns
true if [p] cannot have its [k]-coordinate decreased without leaving the space.
Note
It returns always false for periodic dimension.
Precondition
sIsInside(p) is true.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::sIsOpen ( const SCell p,
Dimension  k 
) const
Parameters
pany signed cell.
kany direction.
Returns
'true' if [p] is open along the direction [k].
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::sIsSurfel ( const SCell b) const
Parameters
bany signed cell.
Returns
'true' if [b] is a surfel (spans all but one coordinate).
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::sIsValid ( const SPreCell c,
Dimension  k 
) const
Parameters
ca signed cell.
ka dimension.
Returns
true if the given signed cell his k-th Khalimsky coordinate between those of the cells returned by lowerCell and upperCell.
Note
For periodic dimension, even if there is no bounds, it guarantees that each cell has unique coordinates.
See also
sCell(const SCell &) const to correct Khalimsky coordinates along periodic dimensions.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::sIsValid ( const SPreCell c) const
Parameters
ca signed cell.
Returns
true if the given signed cell has Khalimsky coordinates between those of the cells returned by lowerCell and upperCell.
Note
For periodic dimension, even if there is no bounds, it guarantees that each cell has unique coordinates.
See also
sCell(const SCell &) const to correct Khalimsky coordinates along periodic dimensions.
template<Dimension dim, typename TInteger>
Size DGtal::KhalimskySpaceND< dim, TInteger >::size ( Dimension  k) const
Parameters
ka dimension.
Returns
the width of the space in the k-dimension.
Note
For periodic dimension, it returns the number of unique coordinates along that dimension.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::sKCoord ( const SCell c,
Dimension  k 
) const
Parameters
cany signed cell.
kany valid dimension.
Returns
its Khalimsky coordinate along [k].
Precondition
uIsValid(c) is true.
template<Dimension dim, typename TInteger>
const Point& DGtal::KhalimskySpaceND< dim, TInteger >::sKCoords ( const SCell c) const
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::sLast ( const SPreCell p,
Dimension  k 
) const
Returns
the k-th Khalimsky coordinate of the last cell of the space with the same type as [p].
Note
For periodic dimension, it returns the last unique coordinate of a cell of same type as p.
Postcondition
The returned coordinate is between lowerCell()[k] and upperCell()[k].
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sLast ( const SPreCell p) const
Returns
the last cell of the space with the same type as [p].
Note
Along periodic dimensions, it returns the last unique coordinate of a cell of same type as p.
Postcondition
sIsValid(sLast(p)) is true.
template<Dimension dim, typename TInteger>
SCells DGtal::KhalimskySpaceND< dim, TInteger >::sLowerIncident ( const SCell c) const
Parameters
cany signed cell.
Returns
the signed cells directly low incident to c in this space.
Note
it is the lower boundary of c expressed as a list of signed cells.
Precondition
sIsValid(c) is true.
Postcondition
sIsValid(b) is true for every returned cell b.
template<Dimension dim, typename TInteger>
SCells DGtal::KhalimskySpaceND< dim, TInteger >::sNeighborhood ( const SCell cell) const

Computes the 1-neighborhood of the cell [c] and returns it. It is the set of cells with same topology that are adjacent to [c] and which are within the bounds of this space.

Parameters
cellthe signed cell of interest.
Returns
the cells of the 1-neighborhood of [cell].
Precondition
sIsValid(cell) is true.
Postcondition
sIsValid(c) is true for every returned cell c.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::sNext ( SCell p,
const SCell lower,
const SCell upper 
) const

Increment the cell [p] to its next position (as classically done in a scanning). Example:

Cell first, last; // lower and upper bounds
Cell p = first;
do
{ // ... whatever [p] is the current cell
}
while ( K.uNext( p, first, last ) );
Parameters
pany cell.
lowerthe lower bound.
upperthe upper bound.
Returns
true if p is still within the bounds, false if the scanning is finished.
Precondition
sIsValid(p) and sIsValid(lower) and sIsValid(upper) and sTopology(p) == sTopology(lower) == sTopology(upper).
Postcondition
sIsValid(p) is true.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sOpp ( const SCell p) const

Creates the signed cell with the inverse sign of [p].

Parameters
pany signed cell.
Returns
the cell [p] with opposite sign.
Precondition
sIsValid(p) is true.
Postcondition
sIsValid(sOpp(p)) is true.
template<Dimension dim, typename TInteger>
Dimension DGtal::KhalimskySpaceND< dim, TInteger >::sOrthDir ( const SCell s) const
template<Dimension dim, typename TInteger>
DirIterator DGtal::KhalimskySpaceND< dim, TInteger >::sOrthDirs ( const SCell p) const

Given a signed cell [p], returns an iterator to iterate over each coordinate the cell does not span. (A spel spans all coordinates; a surfel all but one, etc). Example:

...
for ( KSpace::DirIterator q = ks.uOrthDirs( p ); q != 0; ++q )
{
Dimension dir = *q;
...
}
Parameters
pany signed cell.
Returns
an iterator that points on the first coordinate spanned by the cell.

Referenced by DGtal::functors::SCellToInnerPoint< KSpace >::operator()(), DGtal::functors::SCellToOuterPoint< KSpace >::operator()(), and DGtal::functors::SCellToIncidentPoints< KSpace >::operator()().

template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sPointel ( Point  p,
Sign  sign = POS 
) const

From the digital coordinates of a point in Zn, builds the corresponding pointel (cell of dimension 0) lying into this Khalismky space.

Along a non-periodic dimension, if the given digital coordinate lies outside the space, it is replaced by the nearest valid coordinate.

Along a periodic dimension, the digital coordinate is corrected (by periodicity) to lie between the coordinates of lowerCell() and upperCell().

Parameters
pan integer point (digital coordinates of cell).
signthe sign of the cell (either POS or NEG).
Returns
the signed pointel having the given digital coordinates [p].
Postcondition
uIsValid(uSpel(p)) is true.
Examples:
topology/ctopo-1s-3d.cpp.
template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::sProject ( SCell p,
const SCell bound,
Dimension  k 
) const

Projects [p] along the [k]th direction toward [bound]. Otherwise said, p[ k ] == bound[ k ] afterwards.

Parameters
pany cell.
boundthe element acting as bound (same topology as p).
kthe concerned coordinate.
Precondition
sIsValid(p) and sIsValid(bound) and sIsOpen(p, k) == sIsOpen(bound, k)
Postcondition
sIsValid(p).
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sProjection ( const SCell p,
const SCell bound,
Dimension  k 
) const

Return the projection of [p] along the [k]th direction toward [bound]. Otherwise said, p[ k ] == bound[ k ] afterwards.

Parameters
pany cell.
boundthe element acting as bound (same topology as p).
kthe concerned coordinate.
Returns
the projection.
Precondition
sIsValid(p) and sIsValid(bound) and sIsOpen(p, k) == sIsOpen(bound, k)
Postcondition
sIsValid(sProjection(p, bound, k)) and sTopology(p) == sTopology(sProjection(p, bound, k)).
template<Dimension dim, typename TInteger>
SCells DGtal::KhalimskySpaceND< dim, TInteger >::sProperNeighborhood ( const SCell cell) const

Computes the proper 1-neighborhood of the cell [c] and returns it. It is the set of cells with same topology that are adjacent to [c], different from [c] and which are within the bounds of this space.

Parameters
cellthe signed cell of interest.
Returns
the cells of the proper 1-neighborhood of [cell].
Precondition
sIsValid(cell) is true.
Postcondition
sIsValid(c) is true for every returned cell c.
template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::sSetCoord ( SCell c,
Dimension  k,
Integer  i 
) const

Sets the [k]-th digital coordinate of [c] to [i].

Parameters
cany signed cell.
kany valid dimension.
ian integer coordinate within the space.
Note
The coordinate i will be corrected if k is a periodic dimension.
Precondition
sIsValid(c) is true and i is within the space bounds if k is a non-periodic dimension.
Postcondition
sIsValid(c) is true.
template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::sSetCoords ( SCell c,
const Point kp 
) const

Sets the digital coordinates of [c] to [kp].

Parameters
cany signed cell.
kpthe new Khalimsky coordinates for [c].
Note
The coordinates kp will be corrected for periodic dimensions.
Precondition
kp is within the space bounds for non-periodic dimensions.
Postcondition
sIsValid(c) is true.
template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::sSetKCoord ( SCell c,
Dimension  k,
Integer  i 
) const

Sets the [k]-th Khalimsky coordinate of [c] to [i].

Parameters
cany signed cell.
kany valid dimension.
ian integer coordinate within the space.
Note
The coordinate i will be corrected if k is a periodic dimension.
Precondition
sIsValid(c) is true and i is within the space bounds if k is a non-periodic dimension.
Postcondition
sIsValid(c) is true.
template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::sSetKCoords ( SCell c,
const Point kp 
) const

Sets the Khalimsky coordinates of [c] to [kp].

Parameters
cany signed cell.
kpthe new Khalimsky coordinates for [c].
Note
The coordinates kp will be corrected for periodic dimensions.
Precondition
kp is within the space bounds for non-periodic dimensions.
Postcondition
sIsValid(c) is true.
template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::sSetSign ( SCell c,
Sign  s 
) const

Sets the sign of the cell.

Parameters
c(modified) any signed cell.
sany sign.
template<Dimension dim, typename TInteger>
Sign DGtal::KhalimskySpaceND< dim, TInteger >::sSign ( const SCell c) const
Parameters
cany signed cell.
Returns
its sign.
Precondition
uIsValid(c) is true.
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sSpel ( Point  p,
Sign  sign = POS 
) const

From the digital coordinates of a point in Zn, builds the corresponding spel (cell of maximal dimension) lying into this Khalismky space.

Along a non-periodic dimension, if the given digital coordinate lies outside the space, it is replaced by the nearest valid coordinate.

Along a periodic dimension, the digital coordinate is corrected (by periodicity) to lie between the coordinates of lowerCell() and upperCell().

Parameters
pan integer point (digital coordinates of cell).
signthe sign of the cell (either POS or NEG).
Returns
the signed spel having the given digital coordinates [p].
Postcondition
sIsValid(sSpel(p, sign)) is true.
Examples:
io/boards/dgtalBoard3D-2-ks.cpp, io/viewers/viewer3D-4bis-illustrationMode.cpp, and topology/frontierAndBoundary.cpp.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::sTopology ( const SCell p) const
Parameters
pany signed cell.
Returns
the topology word of [p].
template<Dimension dim, typename TInteger>
SCell DGtal::KhalimskySpaceND< dim, TInteger >::sTranslation ( const SCell p,
const Vector vec 
) const

Add the vector [vec] to [p].

Parameters
pany cell.
vecany pointel.
Returns
the signed code of the cell [p] translated by [coord].
Precondition
sIsValid(p, k) and ( sDistanceToMin(p, k) <= vec[k] <= sDistanceToMax(p, k) or isPeriodicSpace(k) ) for each dimension k.
Postcondition
sIsValid(sTranslation(p, vec)) is true.
Examples:
topology/ctopo-1s-3d.cpp.
template<Dimension dim, typename TInteger>
SCells DGtal::KhalimskySpaceND< dim, TInteger >::sUpperIncident ( const SCell c) const
Parameters
cany signed cell.
Returns
the signed cells directly up incident to c in this space.
Note
it is the upper boundary of c expressed as a list of signed cells.
Precondition
sIsValid(c) is true.
Postcondition
sIsValid(b) is true for every returned cell b.
template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::uAddCoFaces ( Cells cofaces,
const Cell c,
Dimension  axis 
) const
private

Used by uCoFaces for computing incident cofaces.

template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::uAddFaces ( Cells faces,
const Cell c,
Dimension  axis 
) const
private

Used by uFaces for computing incident faces.

template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uAdjacent ( const Cell p,
Dimension  k,
bool  up 
) const
Parameters
pany cell.
kthe coordinate that is changed.
upif 'true' the orientation is forward along axis [k], otherwise backward.
Returns
the adjacent element to [p] along axis [k] in the given direction and orientation.
Note
It is an alias to 'up ? uGetIncr( p, k ) : uGetDecr( p, k )'.
Precondition
uIsValid(p) and not uIs(Min|Max)(p, k) depending on up value.
Postcondition
`uIsValid(uAdjacent(p, k, up)) is true.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uCell ( const PreCell c) const

From an unsigned cell, returns an unsigned cell lying into this Khalismky space.

Along a non-periodic dimension, if the given Khalimsky coordinate lies outside the space, it replaces it by the nearest valid coordinate.

Along a periodic dimension, the Khalimsky coordinate is corrected (by periodicity) to lie between the coordinates of lowerCell() and upperCell().

Parameters
ca cell.
Returns
the same cell with appropriate Khalimsky coordinates along periodic dimensions.
Postcondition
uIsValid(uCell(c)) is true.
uCell(c) == c if uIsValid(c).
Examples:
io/boards/dgtalBoard3DTo2D-KSCell.cpp, topology/ctopo-1-3d.cpp, topology/ctopo-1.cpp, topology/cubical-complex-collapse.cpp, topology/cubical-complex-illustrations.cpp, and topology/khalimskySpaceScanner.cpp.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uCell ( const Point kp) const

From the Khalimsky coordinates of a cell, builds the corresponding unsigned cell lying into this Khalismky space.

Along a non-periodic dimension, if the given Khalimsky coordinate lies outside the space, it is replaced by the nearest valid coordinate.

Along a periodic dimension, the Khalimsky coordinate is corrected (by periodicity) to lie between the coordinates of lowerCell() and upperCell().

Parameters
kpan integer point (Khalimsky coordinates of cell).
Returns
the unsigned cell.
Postcondition
uIsValid(uCell(kp)) is true.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uCell ( Point  p,
const PreCell c 
) const

From the digital coordinates of a point in Zn and a cell type, builds the corresponding unsigned cell lying into this Khalismky space.

Along a non-periodic dimension, if the given digital coordinate lies outside the space, it is replaced by the nearest valid coordinate.

Along a periodic dimension, the digital coordinate is corrected (by periodicity) to lie between the coordinates of lowerCell() and upperCell().

Parameters
pan integer point (digital coordinates of cell).
canother cell defining the topology.
Returns
the cell having the topology of [c] and the given digital coordinates [p].
Postcondition
uIsValid(uCell(p, c)) is true.
template<Dimension dim, typename TInteger>
Cells DGtal::KhalimskySpaceND< dim, TInteger >::uCoFaces ( const Cell c) const
Parameters
cany unsigned cell.
Returns
the proper cofaces of [c] (chain of upper incidence) that belong to the space.
Precondition
uIsValid(c) is true.
Postcondition
uIsValid(b) is true for every returned cell b.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::uCoord ( const Cell c,
Dimension  k 
) const
Parameters
cany unsigned cell.
kany valid dimension.
Returns
its digital coordinate along [k].
Precondition
uIsValid(c) is true.
template<Dimension dim, typename TInteger>
Point DGtal::KhalimskySpaceND< dim, TInteger >::uCoords ( const Cell c) const
Parameters
cany unsigned cell.
Returns
its digital coordinates
Precondition
uIsValid(c) is true.
Examples:
topology/khalimskySpaceScanner.cpp.
template<Dimension dim, typename TInteger>
Dimension DGtal::KhalimskySpaceND< dim, TInteger >::uDim ( const Cell p) const
Parameters
pany unsigned cell.
Returns
the dimension of the cell [p].
template<Dimension dim, typename TInteger>
DirIterator DGtal::KhalimskySpaceND< dim, TInteger >::uDirs ( const Cell p) const

Given an unsigned cell [p], returns an iterator to iterate over each coordinate the cell spans. (A spel spans all coordinates; a surfel all but one, etc). Example:

...
for ( KSpace::DirIterator q = ks.uDirs( p ); q != 0; ++q )
{
Dimension dir = *q;
...
}
Parameters
pany unsigned cell.
Returns
an iterator that points on the first coordinate spanned by the cell.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::uDistanceToMax ( const Cell p,
Dimension  k 
) const

Useful to check if you are going out of the space (for non-periodic dimensions).

Parameters
pany cell.
kthe coordinate that is tested.
Returns
the number of increment to do to reach the maximum value.
Precondition
uIsValid(p) is true.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::uDistanceToMin ( const Cell p,
Dimension  k 
) const

Useful to check if you are going out of the space (for non-periodic dimensions).

Parameters
pany cell.
kthe coordinate that is tested.
Returns
the number of decrement to do to reach the minimum value.
Precondition
uIsValid(p) is true.
template<Dimension dim, typename TInteger>
Cells DGtal::KhalimskySpaceND< dim, TInteger >::uFaces ( const Cell c) const
Parameters
cany unsigned cell.
Returns
the proper faces of [c] (chain of lower incidence) that belong to the space.
Precondition
uIsValid(c) is true.
Postcondition
uIsValid(b) is true for every returned cell b.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::uFirst ( const PreCell p,
Dimension  k 
) const
Returns
the k-th Khalimsky coordinate of the first cell of the space with the same type as [p].
Note
For periodic dimension, it returns the first unique coordinate of a cell of same type as p.
Postcondition
The returned coordinate is between lowerCell()[k] and upperCell()[k].
Examples:
topology/khalimskySpaceScanner.cpp.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uFirst ( const PreCell p) const
Returns
the first cell of the space with the same type as [p].
Note
Along periodic dimensions, it returns the first unique coordinate of a cell of same type as p.
Postcondition
uIsValid(uFirst(p)) is true.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uGetAdd ( const Cell p,
Dimension  k,
Integer  x 
) const
Parameters
pany cell.
kthe coordinate that is changed.
xthe increment.
Returns
the same element as [p] except for a coordinate [k] incremented with x.
Precondition
uIsValid(p) and ( x <= uDistanceToMax(p, k) or isSpacePeriodic(k) ).
Postcondition
uIsValid(uGetAdd(p, k, x)) is true.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uGetDecr ( const Cell p,
Dimension  k 
) const
Parameters
pany cell.
kthe coordinate that is changed.
Returns
the same element as [p] except for an decremented coordinate [k].
Precondition
uIsValid(p) and not uIsMin(p).
Postcondition
uIsValid(uGetDecr(p, k)) is true.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uGetIncr ( const Cell p,
Dimension  k 
) const
Parameters
pany cell.
kthe coordinate that is changed.
Returns
the same element as [p] except for the incremented coordinate [k].
Precondition
uIsValid(p) and not uIsMax(p).
Postcondition
uIsValid(uGetIncr(p, k)) is true.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uGetMax ( Cell  p,
Dimension  k 
) const

Useful to check if you are going out of the space.

Parameters
pany cell.
kthe concerned coordinate.
Returns
the cell similar to [p] but with the maximum allowed [k]-coordinate.
Precondition
uIsValid(p, d) is true for each dimension d different than k.
Postcondition
uIsValid(uGetMax(p, k)) is true.
Examples:
topology/khalimskySpaceScanner.cpp.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uGetMin ( Cell  p,
Dimension  k 
) const

Useful to check if you are going out of the space.

Parameters
pany cell.
kthe concerned coordinate.
Returns
the cell similar to [p] but with the minimum allowed [k]-coordinate.
Precondition
uIsValid(p, d) is true for each dimension d different than k.
Postcondition
uIsValid(uGetMin(p, k)) is true.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uGetSub ( const Cell p,
Dimension  k,
Integer  x 
) const
Parameters
pany cell.
kthe coordinate that is changed.
xthe decrement.
Returns
the same element as [p] except for a coordinate [k] decremented with x.
Precondition
uIsValid(p) and ( x <= uDistanceToMin(p, k) or isSpacePeriodic(k) ).
Postcondition
uIsValid(uGetSub(p, k, x)) is true.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uIncident ( const Cell c,
Dimension  k,
bool  up 
) const
Parameters
cany unsigned cell.
kany coordinate.
upif 'true' the orientation is forward along axis [k], otherwise backward.
Returns
the forward or backward unsigned cell incident to [c] along axis [k], depending on [up].
Note
It may be a lower incident cell if [c] is open along axis [k], else an upper incident cell.
Precondition
uIsValid(c) and the cell should have an incident cell in this direction/orientation.
Postcondition
uIsValid(uIncident(c, k, up)) is true.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::uIsInside ( const PreCell p,
Dimension  k 
) const

Useful to check if you are going out of the space.

Parameters
pany cell.
kthe tested coordinate.
Returns
true if [p] has its [k]-coordinate within the allowed bounds.
Note
It returns always true for periodic dimension.
Examples:
topology/khalimskySpaceScanner.cpp.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::uIsInside ( const PreCell p) const

Useful to check if you are going out of the space.

Parameters
pany cell.
Returns
true if [p] has its coordinates within the allowed bounds.
Note
Only the non-periodic dimensions are checked.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::uIsMax ( const Cell p,
Dimension  k 
) const

Useful to check if you are going out of the space.

Parameters
pany cell.
kthe tested coordinate.
Returns
true if [p] cannot have its [k]-coordinate augmented without leaving the space.
Note
It returns always false for periodic dimension.
Precondition
uIsInside(p) is true.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::uIsMin ( const Cell p,
Dimension  k 
) const

Useful to check if you are going out of the space.

Parameters
pany cell.
kthe tested coordinate.
Returns
true if [p] cannot have its [k]-coordinate decreased without leaving the space.
Note
It returns always false for periodic dimension.
Precondition
uIsInside(p) is true.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::uIsOpen ( const Cell p,
Dimension  k 
) const
Parameters
pany cell.
kany direction.
Returns
'true' if [p] is open along the direction [k].
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::uIsSurfel ( const Cell b) const
Parameters
bany unsigned cell.
Returns
'true' if [b] is a surfel (spans all but one coordinate).
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::uIsValid ( const PreCell c,
Dimension  k 
) const
Parameters
can unsigned cell.
ka dimension.
Returns
true if the given unsigned cell has his k-th Khalimsky coordinate between those of the cells returned by lowerCell and upperCell.
Note
For periodic dimension, even if there is no bounds, it guarantees that each cell has unique coordinates.
See also
uCell(const Cell &) const to correct Khalimsky coordinates along periodic dimensions.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::uIsValid ( const PreCell c) const
Parameters
ca unsigned cell.
Returns
true if the given unsigned cell has Khalimsky coordinates between those of the cells returned by lowerCell and upperCell.
Note
For periodic dimension, even if there is no bounds, it guarantees that each cell has unique coordinates.
See also
uCell(const Cell &) const to correct Khalimsky coordinates along periodic dimensions.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::uKCoord ( const Cell c,
Dimension  k 
) const
Parameters
cany unsigned cell.
kany valid dimension.
Returns
its Khalimsky coordinate along [k].
Precondition
uIsValid(c) is true.
template<Dimension dim, typename TInteger>
const Point& DGtal::KhalimskySpaceND< dim, TInteger >::uKCoords ( const Cell c) const
Parameters
cany unsigned cell.
Returns
its Khalimsky coordinates.
Precondition
uIsValid(c) is true.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::uLast ( const PreCell p,
Dimension  k 
) const
Returns
the k-th Khalimsky coordinate of the last cell of the space with the same type as [p].
Note
For periodic dimension, it returns the last unique coordinate of a cell of same type as p.
Postcondition
The returned coordinate is between lowerCell()[k] and upperCell()[k].
Examples:
topology/khalimskySpaceScanner.cpp.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uLast ( const PreCell p) const
Returns
the last cell of the space with the same type as [p].
Note
Along periodic dimensions, it returns the last unique coordinate of a cell of same type as p.
Postcondition
uIsValid(uLast(p)) is true.
template<Dimension dim, typename TInteger>
Cells DGtal::KhalimskySpaceND< dim, TInteger >::uLowerIncident ( const Cell c) const
Parameters
cany unsigned cell.
Returns
the cells directly low incident to c in this space.
Precondition
uIsValid(c) is true.
Postcondition
uIsValid(b) is true for every returned cell b.
template<Dimension dim, typename TInteger>
Cells DGtal::KhalimskySpaceND< dim, TInteger >::uNeighborhood ( const Cell cell) const

Computes the 1-neighborhood of the cell [c] and returns it. It is the set of cells with same topology that are adjacent to [c] and which are within the bounds of this space.

Parameters
cellthe unsigned cell of interest.
Returns
the cells of the 1-neighborhood of [cell].
Precondition
uIsValid(cell) is true.
Postcondition
uIsValid(c) is true for every returned cell c.
template<Dimension dim, typename TInteger>
bool DGtal::KhalimskySpaceND< dim, TInteger >::uNext ( Cell p,
const Cell lower,
const Cell upper 
) const

Increment the cell [p] to its next position (as classically done in a scanning). Example:

Cell first, last; // lower and upper bounds
Cell p = first;
do
{ // ... whatever [p] is the current cell
}
while ( K.uNext( p, first, last ) );
Parameters
pany cell.
lowerthe lower bound.
upperthe upper bound.
Returns
true if p is still within the bounds, false if the scanning is finished.
Precondition
uIsValid(p) and uIsValid(lower) and uIsValid(upper) and uTopology(p) == uTopology(lower) == uTopology(upper).
Postcondition
uIsValid(p) is true.
Examples:
topology/khalimskySpaceScanner.cpp.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::unsigns ( const SCell p) const

Creates an unsigned cell from a signed one.

Parameters
pany signed cell.
Returns
the unsigned version of the cell [p].
Precondition
sIsValid(p) is true.
Postcondition
uIsValid(unsigns(p)) is true.
Examples:
dec/exampleHeatLaplace.cpp, geometry/surfaces/dvcm-3d.cpp, geometry/surfaces/greedy-plane-segmentation-ex2.cpp, geometry/surfaces/greedy-plane-segmentation.cpp, and topology/volBreadthFirstTraversal.cpp.
template<Dimension dim, typename TInteger>
Dimension DGtal::KhalimskySpaceND< dim, TInteger >::uOrthDir ( const Cell s) const

Given an unsigned surfel [s], returns its orthogonal direction (ie, the coordinate where the surfel is closed).

Parameters
san unsigned surfel
Returns
the orthogonal direction of [s]
template<Dimension dim, typename TInteger>
DirIterator DGtal::KhalimskySpaceND< dim, TInteger >::uOrthDirs ( const Cell p) const

Given an unsigned cell [p], returns an iterator to iterate over each coordinate the cell does not span. (A spel spans all coordinates; a surfel all but one, etc). Example:

...
for ( KSpace::DirIterator q = ks.uOrthDirs( p ); q != 0; ++q )
{
Dimension dir = *q;
...
}
Parameters
pany unsigned cell.
Returns
an iterator that points on the first coordinate spanned by the cell.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uPointel ( Point  p) const

From the digital coordinates of a point in Zn, builds the corresponding pointel (cell of dimension 0) lying into this Khalismky space.

Along a non-periodic dimension, if the given digital coordinate lies outside the space, it is replaced by the nearest valid coordinate.

Along a periodic dimension, the digital coordinate is corrected (by periodicity) to lie between the coordinates of lowerCell() and upperCell().

Parameters
pan integer point (digital coordinates of cell).
Returns
the pointel having the given digital coordinates [p].
Postcondition
uIsValid(uSpel(p)) is true.
Examples:
topology/ctopo-1-3d.cpp, and topology/ctopo-1.cpp.
template<Dimension dim, typename TInteger>
const Point& DGtal::KhalimskySpaceND< dim, TInteger >::upperBound ( ) const
Returns
the upper bound for digital points in this space.
Note
For periodic dimension, it returns the upper digital point with unique coordinates. See also the class documentation.
Examples:
topology/ctopo-2.cpp.
template<Dimension dim, typename TInteger>
const Cell& DGtal::KhalimskySpaceND< dim, TInteger >::upperCell ( ) const
Returns
the upper bound for cells in this space.
Note
For periodic dimension, it returns the upper cell with unique coordinates. See also the class documentation.
template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::uProject ( Cell p,
const Cell bound,
Dimension  k 
) const

Projects [p] along the [k]th direction toward [bound]. Otherwise said, p[ k ] == bound[ k ] afterwards

Parameters
[in,out]pany cell.
[in]boundthe element acting as bound (same topology as p).
[in]kthe concerned coordinate.
Precondition
uIsValid(p) and uIsValid(bound) and uIsOpen(p, k) == uIsOpen(bound, k)
Postcondition
uIsValid(p).
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uProjection ( const Cell p,
const Cell bound,
Dimension  k 
) const

Return the projection of [p] along the [k]th direction toward [bound]. Otherwise said, p[ k ] == bound[ k ] afterwards.

Parameters
pany cell.
boundthe element acting as bound (same topology as p).
kthe concerned coordinate.
Returns
the projection.
Precondition
uIsValid(p) and uIsValid(bound) and uIsOpen(p, k) == uIsOpen(bound, k)
Postcondition
uIsValid(uProjection(p, bound, k)) and uTopology(p) == uTopology(uProjection(p, bound, k)).
template<Dimension dim, typename TInteger>
Cells DGtal::KhalimskySpaceND< dim, TInteger >::uProperNeighborhood ( const Cell cell) const

Computes the proper 1-neighborhood of the cell [c] and returns it. It is the set of cells with same topology that are adjacent to [c], different from [c] and which are within the bounds of this space.

Parameters
cellthe unsigned cell of interest.
Returns
the cells of the proper 1-neighborhood of [cell].
Precondition
uIsValid(cell) is true.
Postcondition
uIsValid(c) is true for every returned cell c.
template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::uSetCoord ( Cell c,
Dimension  k,
Integer  i 
) const

Sets the [k]-th digital coordinate of [c] to [i].

Parameters
cany unsigned cell.
kany valid dimension.
ian integer coordinate within the space.
Note
The coordinate i will be corrected if k is a periodic dimension.
Precondition
uIsValid(c) is true and i is within the space bounds if k is a non-periodic dimension.
Postcondition
uIsValid(c) is true.
template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::uSetCoords ( Cell c,
const Point kp 
) const

Sets the digital coordinates of [c] to [kp].

Parameters
cany unsigned cell.
kpthe new Khalimsky coordinates for [c].
Note
The coordinates kp will be corrected for periodic dimensions.
Precondition
kp is within the space bounds for non-periodic dimensions.
Postcondition
uIsValid(c) is true.
template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::uSetKCoord ( Cell c,
Dimension  k,
Integer  i 
) const

Sets the [k]-th Khalimsky coordinate of [c] to [i].

Parameters
cany unsigned cell.
kany valid dimension.
ian integer coordinate within the space.
Note
The coordinate i will be corrected if k is a periodic dimension.
Precondition
uIsValid(c) is true and i is within the space bounds if k is a non-periodic dimension.
Postcondition
uIsValid(c) is true.
template<Dimension dim, typename TInteger>
void DGtal::KhalimskySpaceND< dim, TInteger >::uSetKCoords ( Cell c,
const Point kp 
) const

Sets the Khalimsky coordinates of [c] to [kp].

Parameters
cany unsigned cell.
kpthe new Khalimsky coordinates for [c].
Note
The coordinates kp will be corrected for periodic dimensions.
Precondition
kp is within the space bounds for non-periodic dimensions.
Postcondition
uIsValid(c) is true.
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uSpel ( Point  p) const

From the digital coordinates of a point in Zn, builds the corresponding spel (cell of maximal dimension) lying into this Khalismky space.

Along a non-periodic dimension, if the given digital coordinate lies outside the space, it is replaced by the nearest valid coordinate.

Along a periodic dimension, the digital coordinate is corrected (by periodicity) to lie between the coordinates of lowerCell() and upperCell().

Parameters
pan integer point (digital coordinates of cell).
Returns
the spel having the given digital coordinates [p].
Postcondition
uIsValid(uSpel(p)) is true.
Examples:
topology/ctopo-1.cpp, topology/cubical-complex-illustrations.cpp, and topology/khalimskySpaceScanner.cpp.
template<Dimension dim, typename TInteger>
Integer DGtal::KhalimskySpaceND< dim, TInteger >::uTopology ( const Cell p) const
Parameters
pany unsigned cell.
Returns
the topology word of [p].
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::uTranslation ( const Cell p,
const Vector vec 
) const

Add the vector [vec] to [p].

Parameters
pany cell.
vecany pointel.
Returns
the unsigned copy of the cell [p] translated by [coord].
Precondition
uIsValid(p, k) and ( uDistanceToMin(p, k) <= vec[k] <= uDistanceToMax(p, k) or isPeriodicSpace(k) ) for each dimension k.
Postcondition
uIsValid(uTranslation(p, vec)) is true.
Examples:
io/viewDualSurface.cpp, topology/ctopo-1-3d.cpp, and topology/ctopo-1.cpp.
template<Dimension dim, typename TInteger>
Cells DGtal::KhalimskySpaceND< dim, TInteger >::uUpperIncident ( const Cell c) const
Parameters
cany unsigned cell.
Returns
the cells directly up incident to c in this space.
Precondition
uIsValid(c) is true.
Postcondition
uIsValid(b) is true for every returned cell b.

Friends And Related Function Documentation

template<Dimension dim, typename TInteger>
friend class KhalimskySpaceNDHelper< KhalimskySpaceND< dim, TInteger > >
friend

Definition at line 396 of file KhalimskySpaceND.h.

Field Documentation

template<Dimension dim, typename TInteger>
const constexpr Dimension DGtal::KhalimskySpaceND< dim, TInteger >::DIM = dim
static

Definition at line 430 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
const constexpr Dimension DGtal::KhalimskySpaceND< dim, TInteger >::dimension = dim
static
template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::myCellLower
private

Definition at line 1932 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
Cell DGtal::KhalimskySpaceND< dim, TInteger >::myCellUpper
private

Definition at line 1933 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
std::array<Closure, dimension> DGtal::KhalimskySpaceND< dim, TInteger >::myClosure
private

Definition at line 1934 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
Point DGtal::KhalimskySpaceND< dim, TInteger >::myLower
private

Definition at line 1930 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
Point DGtal::KhalimskySpaceND< dim, TInteger >::myUpper
private

Definition at line 1931 of file KhalimskySpaceND.h.

template<Dimension dim, typename TInteger>
const constexpr Sign DGtal::KhalimskySpaceND< dim, TInteger >::NEG = false
static
template<Dimension dim, typename TInteger>
const constexpr Sign DGtal::KhalimskySpaceND< dim, TInteger >::POS = true
static

The documentation for this class was generated from the following file: