Fish - FiberLib for VISH 0.3
Fish - The Fiber Bundle API for the Vish Visualization Shell

A multidimensional index that is automatically a lower-dimensional index via recursion. More...

#include <MultiIndex.hpp>

Inheritance diagram for Fiber::MultiIndex< Dims >:
Fiber::FragmentLocation< Dims, TheProperty > Fiber::MultiIndexIterator< Dims >

Public Types

enum  { dims = Dims , SIZE = Dims }
using value_type = index_t

Public Member Functions

constexpr MultiIndex (const MultiIndex &M, const MultiIndex &D, const Add &) noexcept
 Computational constructor for adding two multidimensional indices.
constexpr MultiIndex (const MultiIndex &M, const MultiIndex &D, const Sub &) noexcept
 Computational constructor for subtracting two multidimensional indices.
constexpr MultiIndex (const MultiIndex &M, const MultiIndex &D, const Mult &) noexcept
 Computational constructor for component-wise multiplication of two multidimensional indices.
constexpr MultiIndex (const MultiIndex &M, const MultiIndex &D, const Div &)
 Computational constructor for component-wise division of two multidimensional indices.
constexpr MultiIndex (unsigned int bits, const ::Eagle::BinaryAnd &) noexcept
constexpr MultiIndex (const MultiIndex &M, const MultiIndex &D, const ::Eagle::BinaryAnd &) noexcept
constexpr MultiIndex (const MultiIndex &M, const Power2Alignment &) noexcept
constexpr MultiIndex (const std::array< index_t, Dims > &A)
 Construct from std::array of same size.
constexpr MultiIndex (const std::initializer_list< index_t > &A)
 Construct from set of indices.
MultiIndexIterator< Dims > begin () const
 Begin a ranged loop.
MultiIndexIterator< Dims > end () const
 End a ranged loop.
index_t BitIndex () const
 From this given multiindex which is supposed to have index values of either 0 or 1 in each direction, return the corresponding binary pattern that creates this multiindex.
index_t Orientation () const noexcept
 Given a MultiIndex that only contains one element of value 1, all others being zero (such as one created by the Axis() function), return the orientation that creates this.
constexpr int getSignedOrientation (const MultiIndex &B) const
 Returns the dimension, starting with 1, in which this multidimensional index is non-zero.
constexpr MultiIndex () noexcept
 Default constructor (initialization by zero).
constexpr MultiIndex (const index_t &I) noexcept
 Initialize all indices with same value.
constexpr MultiIndex (const MultiIndex< Dims-1 > &Midx, const index_t &I) noexcept
 Multidimensional index as tensor product of subdimension and index.
constexpr MultiIndex (const MultiIndex< Dims-1 > &m, const index_t &Slice, int SliceDirection)
 Multidimensional index as product of subdimension and index for a given slice dimension (create global index from hyperslab index).
constexpr MultiIndex (const index_t &lin_idx, const MultiIndex &Dimens, const CreateFromLinear &)
 Create multidimensional index from one-dimensional index, given the multidimensional dimension on where this index shall count through.
constexpr const MultiIndex< Dims-1 > & subidx () const noexcept
 Return associated constant dimensionator of one dimension less.
constexpr MultiIndex< Dims-1 > & subidx () noexcept
 Return associated dimensionator of one dimension less.
constexpr const index_toperator[] (std::size_t i) const
 Indexing operator to access values for each dimension.
constexpr index_toperator[] (std::size_t i) noexcept
 Indexing operator to access values for each dimension.
Eagle::Assignment< MultiIndex, 0 > operator= (const index_t &i)
 Assignment via comma operator.
constexpr const index_tmaxidx () const noexcept
 Extent of the current (highest) dimension.
constexpr index_tmaxidx () noexcept
 Extent of the current (highest) dimension.
index_t size () const noexcept
 Recursive function to compute the entire number of elements.
constexpr bool operator!= (const MultiIndex &D) const noexcept
 inequality comparision operator
constexpr bool operator== (const MultiIndex &M) const noexcept
 comparison operator
MultiIndexoperator+= (const MultiIndex &D) noexcept
 component-wise self-addition
MultiIndexoperator-= (const MultiIndex &D) noexcept
 component-wise self-subtraction
bool inc (const MultiIndex &Dimens) noexcept
 Increment the current index, if it is larger than the extent a given in the Dimensions parameter, the certain index is reset to zero in this dimension and the next dimension is increased.
bool inc (const MultiIndex &Dimens, const MultiIndex &Increment) noexcept
 Increment the current index by a certain increment.
index_t linear (const MultiIndex &Dimens) const
 Compute the linear index from the given dimensionator.
bool isWithin (const MultiIndex &Range) const noexcept
 Check whether this multidimensional index is within the specified domain.
bool operator< (const MultiIndex &Range) const noexcept
 Check if the current multiindex is smaller than the range.
bool operator> (const MultiIndex &Range) const noexcept
 Check if the current multiindex is larger than the range.
MultiIndex operator+ (const MultiIndex &D) const noexcept
 Add two multidimensional indices.
MultiIndex operator- (const MultiIndex &D) const noexcept
 Subtract multidimensional indices.
MultiIndex operator* (const MultiIndex &D) const noexcept
 Multiply multidimensional indices component-wise.
MultiIndex operator/ (const MultiIndex &D) const noexcept
 Divide multidimensional indices component-wise.
MultiIndex operator& (int i) const
 Shortcut operator to compute the MultiIndex that is just one step aside in the direction as specified by the integer.

Static Public Member Functions

static MultiIndex BitIndex (unsigned int bits) noexcept
 Multidimensional bit indices, that are zero or one in each of the directions as specified by the bits.
static MultiIndex Axis (unsigned int orientation) noexcept
 Return a MultiIndex that points just in the given orientation.
static int log2 (index_t N) noexcept
 Logarithm of basis 2 where we expect a maximum value of 1<<(Dims-1) here.
template<class Functor>
static bool ForEach (Functor &F, const MultiIndex &Start, const MultiIndex &End, const MultiIndex &Increment=MultiIndex(1))

Static Protected Member Functions

template<class Functor, Dims_t SuperDims>
static bool ForEachRecursion (Functor &F, const MultiIndex< SuperDims > &SuperIndex, const MultiIndex &Start, const MultiIndex &End, const MultiIndex &Increment)

Friends

MultiIndex distance (const MultiIndex &M0, const MultiIndex &M1)
 Compute the distance between two multiindices, which is always positive in each index (in contrast to the subtraction operation).
MultiIndex bitalign (const MultiIndex &M) noexcept
MultiIndex operator| (const MultiIndex &M, int i)
 Axis index subtraction operator, similar to adding an integer to a multidimensional index, but opposite direction.

(Note that these are not member symbols.)

template<Dims_t Dims>
constexpr FixedArray< double, Dims > div (const MultiIndex< Dims > &Divisor, const MultiIndex< Dims > &Dividend)
 Compute the rational division of two MultiIndex es.
template<class VectorType, class Domain, class scalar_t>
Eagle::DomainVector< VectorType, Domain, scalar_t > CellSize (const Eagle::DomainVector< VectorType, Domain, scalar_t > &V, const MultiIndex< VectorType::SIZE > &Dividend, const MultiIndex< VectorType::SIZE > &Stride=MultiIndex< VectorType::SIZE >(1))
 Compute the size of a "cell" given a distance between a set of points.
template<class VectorType, class Domain, class scalar_t>
Eagle::DomainVector< VectorType, Domain, scalar_t > CellSize0 (const Eagle::DomainVector< VectorType, Domain, scalar_t > &V, const MultiIndex< VectorType::SIZE > &Dividend, const MultiIndex< VectorType::SIZE > &Stride=MultiIndex< VectorType::SIZE >(1))
 Compute the size of a "cell" given a distance between a set of points.
template<Dims_t Dims>
MultiIndex< Dims > operator+ (index_t C, const MultiIndex< Dims > &I)
 Add constant plus multi-index.
template<Dims_t Dims>
MultiIndex< Dims > operator+ (const MultiIndex< Dims > &I, index_t C)
 Add multi-index plus constant.
template<Dims_t Dims>
FixedArray< double, Dims > operator* (const FixedArray< double, Dims > &a, const MultiIndex< Dims > &I)
 Component-wise multiplication.
template<Dims_t Dims>
FixedArray< double, Dims > operator* (const MultiIndex< Dims > &I, const FixedArray< double, Dims > &a)
 Component-wise multiplication.
template<Dims_t Dims>
FixedArray< double, Dims > operator/ (const FixedArray< double, Dims > &a, const MultiIndex< Dims > &I)
 Component-wise division.
template<Dims_t Dims>
FixedArray< double, Dims > & operator/= (FixedArray< double, Dims > &a, const MultiIndex< Dims > &I)
 Component-wise division.
template<Dims_t Dims>
FixedArray< double, Dims > & operator*= (FixedArray< double, Dims > &a, const MultiIndex< Dims > &I)
 Component-wise multiplication.
template<Dims_t Dims>
FixedArray< int, Dims > MultiIndexToFixedArray (const MultiIndex< Dims > &Value)
template<Dims_t Dims>
MultiIndex< Dims > FixedArrayToMultiIndex (const FixedArray< int, Dims > &Value)
template<Dims_t Dims>
MultiIndex< Dims > minMultiIndex (const MultiIndex< Dims > &A, const MultiIndex< Dims > &B)
template<Dims_t Dims>
MultiIndex< Dims > maxMultiIndex (const MultiIndex< Dims > &A, const MultiIndex< Dims > &B)

Detailed Description

template<Dims_t Dims>
class Fiber::MultiIndex< Dims >

A multidimensional index that is automatically a lower-dimensional index via recursion.

Recursively defined multidimensional index type.

Template Parameters
DimsNumber of dimensions (Dims >= 1).

Overview

MultiIndex<Dims> represents a D-dimensional index where each dimension stores a single integer coordinate. The class is defined recursively:

  • MultiIndex<1> stores a single coordinate `idx`.
  • MultiIndex<D> (D > 1) inherits from MultiIndex<D-1> and adds one additional coordinate `idx` representing the highest dimension.

Thus, a MultiIndex<D> contains D coordinates, accessible via operator[]. The lowest dimension is stored in the base class, the highest in the derived class.

A recursively defined multidimensional index. The base class corresponds to the lowest dimension, each new child class adds a new dimension.

A MultiIndex ("multidimensional Index") can be: A Type Product between two Indices An Index with Dimensionator

Semantics

  • `idx` (and therefore `maxidx()`) always refers to the coordinate of the *current* (highest) dimension in this recursion layer.
  • `subidx()` returns the (D-1)-dimensional prefix.
  • `size()` returns the total number of elements in the D-dimensional space: size() = maxidx() * subidx().size() For D=1, size() == maxidx().
  • Component-wise arithmetic (Add, Sub, Mult, Div) is implemented by applying the operation to the highest dimension and recursively to all lower dimensions.
  • Linearization and de-linearization follow standard row-major order: linear = idx * subidx().size() + subidx().linear(...) de-linearization splits the linear index into high and low parts.

Construction

MultiIndex<D> can be constructed in several ways:

Access

Coordinates are accessed via:

idx = MI[i]; // 0 <= i < Dims

where i = 0 refers to the lowest dimension and i = Dims-1 to the highest.

Arithmetic

MultiIndex supports component-wise:

MI + MJ MI - MJ MI * MJ MI / MJ

and their in-place variants. Equality and inequality compare all dimensions recursively.

Iteration

MultiIndex provides:

begin(), end()

returning MultiIndexIterator<Dims> for range-based iteration over the highest dimension. Full multidimensional iteration is provided by:

MultiIndex<Dims>::ForEach(Functor, Start, End, Increment)

which recursively iterates all dimensions in row-major order.

Operations

MultiIndex supports bit-index operations for hypercube traversal:

  • BitIndex(bits): construct a D-dimensional 0/1 vector from bitmask.
  • Axis(i): return a unit vector along dimension i.
  • BitIndex(): convert a 0/1 MultiIndex back into a bitmask.
  • Orientation(): return the highest dimension whose bit is set.

Linearization

Given a dimension descriptor Dimens (another MultiIndex<Dims> whose coordinates specify extents), the following hold:

  • `linear(Dimens)` returns the row-major linear index.
  • `MultiIndex(lin, Dimens, CreateFromLinear{})` reconstructs the multidimensional index from a linear index.

Checking

  • `isWithin(Range)` returns true if all coordinates are strictly less than the corresponding coordinates of Range.
  • `operator<` and `operator>` use isWithin() semantics, not lexicographic comparison.

Increment

  • `inc(Dimens)` increments the multidimensional index in row-major order, wrapping lower dimensions first.
  • `inc(Dimens, Increment)` increments by a multidimensional step vector.

Summary

MultiIndex<Dims> is a compact, recursive, constexpr-capable representation of multidimensional indices supporting:

  • Component-wise arithmetic
  • Row-major linearization and reconstruction
  • Bit-index operations for hypercube traversal
  • Recursive iteration across all dimensions
  • Full coordinate access via operator[]

The recursive structure ensures that all multidimensional operations naturally decompose into operations on lower-dimensional indices.

Examples
ColoredLines.cpp, CreatorExample.cpp, InterpolatedVectorDerivative.cpp, Interpolation.cpp, ModifyingPointCloudRange.cpp, XF_LineSetHierarchicalRegularized.cpp, and uniform_scalar.cpp.

Member Enumeration Documentation

◆ anonymous enum

template<Dims_t Dims>
anonymous enum
Enumerator
dims 

The dimension of this multidimensional index.

SIZE 

Export an SIZE enum for treatment like an FixedArray.

Constructor & Destructor Documentation

◆ MultiIndex() [1/3]

template<Dims_t Dims>
Fiber::MultiIndex< Dims >::MultiIndex ( const std::initializer_list< index_t > & A)
inlineconstexpr

Construct from set of indices.

MultiIndex<3> MI = { 33, 44, 66 };
constexpr MultiIndex(const MultiIndex &M, const MultiIndex &D, const Add &) noexcept
Computational constructor for adding two multidimensional indices.
Definition MultiIndex.hpp:458

◆ MultiIndex() [2/3]

template<Dims_t Dims>
Fiber::MultiIndex< Dims >::MultiIndex ( const MultiIndex< Dims-1 > & m,
const index_t & Slice,
int SliceDirection )
inlineconstexpr

Multidimensional index as product of subdimension and index for a given slice dimension (create global index from hyperslab index).

Todo
Optimize this to peform recursive construction!

◆ MultiIndex() [3/3]

template<Dims_t Dims>
Fiber::MultiIndex< Dims >::MultiIndex ( const index_t & lin_idx,
const MultiIndex< Dims > & Dimens,
const CreateFromLinear &  )
inlineconstexpr

Create multidimensional index from one-dimensional index, given the multidimensional dimension on where this index shall count through.

This is the inverse operation to the linear() member function.

Test code:

{
MultiIndex<3> Number = MIndex(11, 17, 23);
for(index_t I = 11; I<331; I+=29)
{
MultiIndex<3> I3(I, Number);
index_t Irec = I3.linear( Number );
cout << "Index " << I << " with I[0]=" << I3[0] << " is " << I3 << " was " << Irec << endl;
}
}
basic_ostream< _CharT, _Traits > & endl(basic_ostream< _CharT, _Traits > &__os)
ostream cout
IndexTypeConfig< sizeof(void *)>::index_t index_t
Define the index type as according to the size of a pointer, i.e.
Definition Index.hpp:22

Member Function Documentation

◆ Axis()

template<Dims_t Dims>
MultiIndex Fiber::MultiIndex< Dims >::Axis ( unsigned int orientation)
inlinestaticnoexcept

Return a MultiIndex that points just in the given orientation.

It is given by

BitIndex( 1<< orientation )
index_t BitIndex() const
From this given multiindex which is supposed to have index values of either 0 or 1 in each direction,...
Definition MultiIndex.hpp:607

with orientation < Dims . The inverse operation is the Orientation() member function.

Referenced by Fiber::RegularTopology::ComputeEdgeIDfromVertexAndOrientation(), Fiber::RegularTopology::ComputeFirstEdgeVertex(), Fiber::RegularTopology::EdgeOrientation(), Fiber::MultiIndex< 2 >::operator&(), Fiber::MultiIndex< 2 >::operator|, and Fiber::PartialDerivative().

◆ BitIndex() [1/2]

template<Dims_t Dims>
index_t Fiber::MultiIndex< Dims >::BitIndex ( ) const
inline

From this given multiindex which is supposed to have index values of either 0 or 1 in each direction, return the corresponding binary pattern that creates this multiindex.

cout << "Check bit indices, all corners of a cube: " << endl;
for(unsigned i=0; i<8; i++)
{
cout << i << " --> " << MultiIndex<3>::BitIndex( i )
<< " --> " << MultiIndex<3>::BitIndex( i ).BitIndex()
<< endl;
}

Referenced by Fiber::MultiIndex< 2 >::Axis(), Fiber::RegularTopology::ComputeSecondEdgeVertex(), Fiber::RegularTopology::NumberOfEdges(), Fiber::MultiIndex< 2 >::Orientation(), and Fiber::Polygonize().

◆ BitIndex() [2/2]

template<Dims_t Dims>
MultiIndex Fiber::MultiIndex< Dims >::BitIndex ( unsigned int bits)
inlinestaticnoexcept

Multidimensional bit indices, that are zero or one in each of the directions as specified by the bits.

They can be used to compute the corners of an n-dimensional hypercube or the direction of the n-th axis.

cout << "Check bit indices, all corners of a cube: " << endl;
for(unsigned i=0; i<8; i++)
{
cout << i << " --> " << MultiIndex<3>::BitIndex( i ) << endl;
}
cout << "Check bit indices, all axis of a cube: " << endl;
for(unsigned i=0; i<3; i++)
{
cout << i << " --> " << MultiIndex<3>::BitIndex( 1<<i ) << endl;
}

Referenced by Fiber::renderVolume().

◆ getSignedOrientation()

template<Dims_t Dims>
int Fiber::MultiIndex< Dims >::getSignedOrientation ( const MultiIndex< Dims > & B) const
inlineconstexpr

Returns the dimension, starting with 1, in which this multidimensional index is non-zero.

The result will be negative if the specified index is smaller than the current one, and zero, if both indices are identical. If there is more than one non-zero index, then it will return the highest dimension where an index difference is found.

Referenced by Fiber::ComputeInterpolationEdgeAndWeight(), and Fiber::MultiIndex< 2 >::getSignedOrientation().

◆ inc()

template<Dims_t Dims>
bool Fiber::MultiIndex< Dims >::inc ( const MultiIndex< Dims > & Dimens,
const MultiIndex< Dims > & Increment )
inlinenoexcept

Increment the current index by a certain increment.

If it is larger than the extent a given in the Dimensions parameter, the certain index is reset to zero in this dimension and the next dimension is increased.

◆ linear()

◆ operator&()

template<Dims_t Dims>
MultiIndex Fiber::MultiIndex< Dims >::operator& ( int i) const
inline

Shortcut operator to compute the MultiIndex that is just one step aside in the direction as specified by the integer.

Adding 0 yields the original vector, adding +1 moves one element forward in the first axis, adding -1 moves one element backward along the first axis, etc. . The idea stems from Massimo Pierro's FermiQCD package.

◆ operator[]()

template<Dims_t Dims>
const index_t & Fiber::MultiIndex< Dims >::operator[] ( std::size_t i) const
inlineconstexpr

Indexing operator to access values for each dimension.

The current class member's index corresponds to the highest dimension, the base class to the lowest dimension.

◆ Orientation()

template<Dims_t Dims>
index_t Fiber::MultiIndex< Dims >::Orientation ( ) const
inlinenoexcept

Given a MultiIndex that only contains one element of value 1, all others being zero (such as one created by the Axis() function), return the orientation that creates this.

cout << "Check bit indices, all axis of a cube: " << endl;
for(unsigned i=0; i<3; i++)
{
cout << i << " --> " << MultiIndex<3>::BitIndex( 1<<i )
<< " --> " << MultiIndex<3>::BitIndex( 1<<i ).Orientation()
<< endl;
}
Note
If more than one index is defined, then the orientation of the highest one will be returned. If there are values other than one or zero, the result is undefined.