CGAL 5.2 - Bounding Volumes
CGAL::Approximate_min_ellipsoid_d< Traits > Class Template Reference

#include <CGAL/Approximate_min_ellipsoid_d.h>

Definition

An object of class Approximate_min_ellipsoid_d is an approximation to the ellipsoid of smallest volume enclosing a finite multiset of points in \( d\)-dimensional Euclidean space \( \E^d\), \( d\ge 2\).

An ellipsoid in \( \E^d\) is a Cartesian pointset of the form \( \{ x\in\E^d \mid x^T E x + x^T e + \eta\leq 0 \}\), where \( E\) is some positive definite matrix from the set \( \mathbb{R}^{d\times d}\), \( e\) is some real \( d\)-vector, and \( \eta\in\mathbb{R}\). A pointset \( P\subseteq \E^d\) is called full-dimensional if its affine hull has dimension \( d\). For a finite, full-dimensional pointset \( P\) we denote by \( (P)\) the smallest ellipsoid that contains all points of \( P\); this ellipsoid exists and is unique.

For a given finite and full-dimensional pointset \( P\subset \E^d\) and a real number \( \epsilon\ge 0\), we say that an ellipsoid \( {\cal E}\subset\E^d\) is an \( (1+\epsilon)\)-appoximation to \( (P)\) if \( P\subset {\cal E}\) and \( ({\cal E}) \leq (1+\epsilon) ((P))\). In other words, an \( (1+\epsilon)\)-approximation to \( (P)\) is an enclosing ellipsoid whose volume is by at most a factor of \( 1+\epsilon\) larger than the volume of the smallest enclosing ellipsoid of \( P\).

Given this notation, an object of class Approximate_min_ellipsoid_d represents an \( (1+\epsilon)\)-approximation to \( (P)\) for a given finite and full-dimensional multiset of points \( P\subset\E^d\) and a real constant \( \epsilon>0\).A multiset is a set where elements may have multiplicity greater than \( 1\). When an Approximate_min_ellipsoid_d<Traits> object is constructed, an iterator over the points \( P\) and the number \( \epsilon\) have to be specified; the number \( \epsilon\) defines the desired approximation ratio \( 1+\epsilon\). The underlying algorithm will then try to compute an \( (1+\epsilon)\)-approximation to \( (P)\), and one of the following two cases takes place.

  • The algorithm determines that \( P\) is not full-dimensional (see is_full_dimensional() below).

    Important note: due to rounding errors, the algorithm cannot in all cases decide correctly whether \( P\) is full-dimensional or not. If is_full_dimensional() returns false, the points lie in such a "thin" subspace of \( \E^d\) that the algorithm is incapable of computing an approximation to \( (P)\). More precisely, if is_full_dimensional() returns false, there exist two parallel hyperplanes in \( \E^d\) with the points \( P\) in between so that the distance \( \delta\) between the hyperplanes is very small, possible zero. (If \( \delta=0\) then \( P\) is not full-dimensional.)

    If \( P\) is not full-dimensional, linear algebra techniques should be used to determine an affine subspace \( S\) of \( \E^d\) that contains the points \( P\) as a (w.r.t. \( S\)) full-dimensional pointset; once \( S\) is determined, the algorithm can be invoked again to compute an approximation to (the lower-dimensional) \( (P)\) in \( S\). Since is_full_dimensional() might (due to rounding errors, see above) return false even though \( P\) is full-dimensional, the lower-dimensional subspace \( S\) containing \( P\) need not exist. Therefore, it might be more advisable to fit a hyperplane \( H\) through the pointset \( P\), project \( P\) onto this affine subspace \( H\), and compute an approximation to the minimum-volume enclosing ellipsoid of the projected points within \( H\); the fitting can be done for instance using the linear_least_squares_fitting() function from the CGAL package Principal_component_analysis.

  • The algorithm determines that \( P\) is full-dimensional. In this case, it provides an approximation \( {\cal E}\) to \( (P)\), but depending on the input problem (i.e., on the pair \( (P,\epsilon)\)), it may not have achieved the desired approximation ratio but merely some worse approximation ratio \( 1+\epsilon'>1+\epsilon\). The achieved approximation ratio \( 1+\epsilon'\) can be queried using achieved_epsilon(), which returns \( \epsilon'\). The ellipsoid \( {\cal E}\) itself can be queried via the methods defining_matrix(), defining_vector(), and defining_scalar().

The ellipsoid \( {\cal E}\) computed by the algorithm satisfies the inclusions

\[ \frac{1}{(1+\epsilon')d} {\cal E} \subseteq \mathop{\rm conv}\nolimits(P) \subseteq {\cal E} \]

where \( f {\cal E}\) denotes the ellipsoid \( {\cal E}\) scaled by the factor \( f\in\mathbb{R}^+\) with respect to its center, and where \( \mathop{\rm conv}\nolimits(A)\) denotes the convex hull of a pointset \( A\subset \E^d\).

The underlying algorithm can cope with all kinds of inputs (multisets \( P\), \( \epsilon\in[0,\infty)\)) and terminates in all cases. There is, however, no guarantee that any desired approximation ratio is actually achieved; the performance of the algorithm in this respect highly depends on the input pointset. Values of at least \( 0.01\) for \( \epsilon\) are usually handled without problems.

Internally, the algorithm represents the input points' Cartesian coordinates as double's. For this conversion to work, the input point coordinates must be convertible to double. Also, in order to compute the achieved epsilon \( \epsilon'\) mentioned above, the algorithm requires a number type ET that provides exact arithmetic. (Both these aspects are discussed in the documentation of the concept ApproximateMinEllipsoid_d_Traits_d.)

Template Parameters
Traitsmust be a model for ApproximateMinEllipsoid_d_Traits_d.

We provide the model CGAL::Approximate_min_ellipsoid_d_traits_d<K> using the \( d\)-dimensional CGAL kernel; the models CGAL::Approximate_min_ellipsoid_d_traits_2<K> and CGAL::Approximate_min_ellipsoid_d_traits_3<K> are for use with the \( 2\)- and \( 3\)-dimensional CGAL kernel, respectively.

See also
CGAL::Min_ellipse_2<Traits>

Implementation

We implement Khachyian's algorithm for rounding polytopes [10]. Internally, we use double-arithmetic and (initially a single) Cholesky-decomposition. The algorithm's running time is \( {\cal O}(nd^2(\epsilon^{-1}+\ln d + \ln\ln(n)))\), where \( n=|P|\) and \( 1+\epsilon\) is the desired approximation ratio.

Example

To illustrate the usage of Approximate_min_ellipsoid_d we give two examples in 2D. The first program generates a random set \( P\subset\E^2\) and outputs the points and a \( 1.01\)-approximation of \( (P)\) as an EPS-file, which you can view using gv, for instance. (In both examples you can change the variables n and d to experiment with the code.)


File Approximate_min_ellipsoid_d/ellipsoid.cpp

#include <CGAL/Cartesian_d.h>
#include <CGAL/MP_Float.h>
#include <CGAL/point_generators_d.h>
#include <CGAL/Approximate_min_ellipsoid_d.h>
#include <CGAL/Approximate_min_ellipsoid_d_traits_d.h>
#include <vector>
#include <iostream>
typedef CGAL::Cartesian_d<double> Kernel;
typedef CGAL::MP_Float ET;
typedef Traits::Point Point;
typedef std::vector<Point> Point_list;
int main()
{
const int n = 1000; // number of points
const int d = 2; // dimension
const double eps = 0.01; // approximation ratio is (1+eps)
// create a set of random points:
Point_list P;
CGAL::Random_points_in_cube_d<Point> rpg(d,100.0);
for (int i = 0; i < n; ++i) {
P.push_back(*rpg);
++rpg;
}
// compute approximation:
Traits traits;
AME ame(eps, P.begin(), P.end(), traits);
// write EPS file:
if (ame.is_full_dimensional() && d == 2)
ame.write_eps("example.eps");
// output center coordinates:
std::cout << "Cartesian center coordinates: ";
for (AME::Center_coordinate_iterator c_it = ame.center_cartesian_begin();
c_it != ame.center_cartesian_end();
++c_it)
std::cout << *c_it << ' ';
std::cout << ".\n";
if (d == 2 || d == 3) {
// output axes:
AME::Axes_lengths_iterator axes = ame.axes_lengths_begin();
for (int i = 0; i < d; ++i) {
std::cout << "Semiaxis " << i << " has length " << *axes++ << "\n"
<< "and Cartesian coordinates ";
for (AME::Axes_direction_coordinate_iterator
d_it = ame.axis_direction_cartesian_begin(i);
d_it != ame.axis_direction_cartesian_end(i); ++d_it)
std::cout << *d_it << ' ';
std::cout << ".\n";
}
}
}

The second program outputs the approximation in a format suitable for display in Maplesoft's Maple.


File Approximate_min_ellipsoid_d/ellipsoid_for_maple.cpp

// Usage: ./maple_example > maple.text
// Then enter in Maple 'read "maple.text";'.
#include <CGAL/Cartesian_d.h>
#include <CGAL/MP_Float.h>
#include <CGAL/point_generators_d.h>
#include <CGAL/Approximate_min_ellipsoid_d.h>
#include <CGAL/Approximate_min_ellipsoid_d_traits_d.h>
#include <vector>
#include <iostream>
#include <iomanip>
typedef CGAL::Cartesian_d<double> Kernel;
typedef CGAL::MP_Float ET;
typedef Traits::Point Point;
typedef std::vector<Point> Point_list;
int main()
{
const int n = 100; // number of points
const int d = 2; // dimension
const double eps = 0.01; // approximation ratio is (1+eps)
// create a set of random points:
Point_list P;
CGAL::Random_points_in_cube_d<Point> rpg(d,1.0);
for (int i = 0; i < n; ++i) {
P.push_back(*rpg);
++rpg;
}
// compute approximation:
Traits traits;
AME mel(eps, P.begin(), P.end(), traits);
// output for Maple:
if (mel.is_full_dimensional() && d == 2) {
const double alpha = (1+mel.achieved_epsilon())*(d+1);
// output points:
using std::cout;
cout << "restart;\n"
<< "with(LinearAlgebra):\n"
<< "with(plottools):\n"
<< "n:= " << n << ":\n"
<< "P:= Matrix(" << d << "," << n << "):\n";
for (int i=0; i<n; ++i)
for (int j=0; j<d; ++j)
cout << "P[" << j+1 << "," << i+1 << "] := "
<< std::setiosflags(std::ios::scientific)
<< std::setprecision(20) << P[i][j] << ":\n";
cout << "\n";
// output defining equation:
cout << "Mp:= Matrix([\n";
for (int i=0; i<d; ++i) {
cout << " [";
for (int j=0; j<d; ++j) {
cout << mel.defining_matrix(i,j)/alpha;
if (j<d-1)
cout << ",";
}
cout << "]";
if (i<d-1)
cout << ",";
cout << "\n";
}
cout << "]);\n" << "mp:= Vector([";
for (int i=0; i<d; ++i) {
cout << mel.defining_vector(i)/alpha;
if (i<d-1)
cout << ",";
}
cout << "]);\n"
<< "eta:= " << (mel.defining_scalar()/alpha-1.0) << ";\n"
<< "v:= Vector([x,y]):\n"
<< "e:= Transpose(v).Mp.v+Transpose(v).mp+eta;\n"
<< "plots[display]({seq(point([P[1,i],P[2,i]]),i=1..n),\n"
<< " plots[implicitplot](e,x=-5..5,y=-5..5,numpoints=10000)},\n"
<< " scaling=CONSTRAINED);\n";
}
}
Note
This class requires the Eigen library.
Examples:
Approximate_min_ellipsoid_d/ellipsoid.cpp, and Approximate_min_ellipsoid_d/ellipsoid_for_maple.cpp.

Types

typedef unspecified_type FT
 typedef Traits::FT FT (which is always a typedef to double).
 
typedef unspecified_type ET
 typedef Traits::ET ET (which is an exact number type used for exact computation like for example in achieved_epsilon()).
 
typedef unspecified_type Point
 typedef Traits::Point Point
 
typedef unspecified_type Cartesian_const_iterator
 typedef Traits::Cartesian_const_iterator Cartesian_const_iterator
 
typedef unspecified_type Center_coordinate_iterator
 A model of STL concept RandomAccessIterator with value type double that is used to iterate over the Cartesian center coordinates of the computed ellipsoid, see center_cartesian_begin().
 
typedef unspecified_type Axes_lengths_iterator
 A model of STL concept RandomAccessIterator with value type double that is used to iterate over the lengths of the semiaxes of the computed ellipsoid, see axes_lengths_begin().
 
typedef unspecified_type Axes_direction_coordinate_iterator
 A model of STL concept RandomAccessIterator with value type double that is used to iterate over the Cartesian coordinates of the direction of a fixed axis of the computed ellipsoid, see axis_direction_cartesian_begin().
 

Creation

An object of type Approximate_min_ellipsoid_d can be created from an arbitrary point set \( P\) and some nonnegative double value eps.

template<class Iterator >
 Approximate_min_ellipsoid_d (double eps, Iterator first, Iterator last, const Traits &traits=Traits())
 initializes ame to an \( (1+\epsilon)\)-approximation of \( (P)\) with \( P\) being the set of points in the range [first,last). More...
 

Access Functions

The following methods can be used to query the achieved approximation ratio \( 1+\epsilon'\) and the computed ellipsoid \( {\cal E} = \{ x\in\E^d \mid x^T E x + x^T e + \eta\leq 0 \}\).

The methods defining_matrix(), defining_vector(), and defining_scalar() do not return \( E\), \( e\), and \( \eta\) directly but yield multiples of these quantities that are exactly representable using the double type. (This is necessary because the parameters \( E\), \( e\), and \( \eta\) of the computed approximation ellipsoid \( {\cal E}\) might not be exactly representable as double numbers.) In order to access the center and semiaxes of the computed approximation ellipsoid, the functions center_cartesian_begin(), axes_lengths_begin(), and axis_direction_cartesian_begin() can be used. In constrast to the above access functions achieved_epsilon(), defining_matrix(), defining_vector(), and defining_scalar(), which return the described quantities exactly, the routines below return numerical approximations to the real center and real semiaxes of the computed ellipsoid; the comprised relative error may be larger than zero, and there are no guarantees for the returned quantities.

unsigned int number_of_points () const
 returns the number of points of ame, i.e., \( |P|\).
 
double achieved_epsilon () const
 returns a number \( \epsilon'\) such that the computed approximation is (under exact arithmetic) guaranteed to be an \( (1+\epsilon')\)-approximation to \( (P)\). More...
 
double defining_matrix (int i, int j) const
 gives access to the \( (i,j)\)th entry of the matrix \( E\) in the representation \( \{ x\in\E^d \mid x^T E x + x^T e + \eta\leq0 \}\) of the computed approximation ellipsoid \( {\cal E}\). More...
 
double defining_vector (int i) const
 gives access to the \( i\)th entry of the vector \( e\) in the representation \( \{ x\in\E^d \mid x^T E x + x^T e + \eta\leq0 \}\) of the computed approximation ellipsoid \( {\cal E}\). More...
 
double defining_scalar () const
 gives access to the scalar \( \eta\) from the representation \( \{ x\in\E^d \mid x^T E x + x^T e + \eta\leq0 \}\) of the computed approximation ellipsoid \( {\cal E}\). More...
 
const Traits & traits () const
 returns a const reference to the traits class object.
 
int dimension () const
 returns the dimension of the ambient space, i.e., the dimension of the points \( P\).
 
Center_coordinate_iterator center_cartesian_begin ()
 returns an iterator pointing to the first of the \( d\) Cartesian coordinates of the computed ellipsoid's center. More...
 
Center_coordinate_iterator center_cartesian_end ()
 returns the past-the-end iterator corresponding to center_cartesian_begin(). More...
 
Axes_lengths_iterator axes_lengths_begin ()
 returns an iterator pointing to the first of the \( d\) descendantly sorted lengths of the computed ellipsoid's axes. More...
 
Axes_lengths_iterator axes_lengths_end ()
 returns the past-the-end iterator corresponding to axes_lengths_begin(). More...
 
Axes_direction_coordinate_iterator axis_direction_cartesian_begin (int i)
 returns an iterator pointing to the first of the \( d\) Cartesian coordinates of the computed ellipsoid's \( i\)th axis direction (i.e., unit vector in direction of the ellipsoid's \( i\)th axis). More...
 
Axes_direction_coordinate_iterator axis_direction_cartesian_end (int i)
 returns the past-the-end iterator corresponding to axis_direction_cartesian_begin(). More...
 

Predicates

bool is_full_dimensional () const
 returns whether \( P\) is full-dimensional or not, i.e., returns true if and only if \( P\) is full-dimensional. More...
 

Validity Check

An object ame is valid iff

  • ame contains all points of its defining set \( P\),
  • ame is an \( (1+\epsilon')\)-approximation to the smallest ellipsoid \( (P)\) of \( P\),
  • The ellipsoid represented by ame fulfills the inclusion ( eqapproximate_min_ellipsoid_incl ).
bool is_valid (bool verbose=false) const
 returns true iff ame is valid according to the above definition. More...
 

Miscellaneous

void write_eps (const std::string &name) const
 Writes the points \( P\) and the computed approximation to \( (P)\) as an EPS-file under pathname name. More...
 

Constructor & Destructor Documentation

◆ Approximate_min_ellipsoid_d()

template<typename Traits >
template<class Iterator >
CGAL::Approximate_min_ellipsoid_d< Traits >::Approximate_min_ellipsoid_d ( double  eps,
Iterator  first,
Iterator  last,
const Traits &  traits = Traits() 
)

initializes ame to an \( (1+\epsilon)\)-approximation of \( (P)\) with \( P\) being the set of points in the range [first,last).

The number \( \epsilon\) in this will be at most eps, if possible. However, due to the limited precision in the algorithm's underlying arithmetic, it can happen that the computed approximation ellipsoid has a worse approximation ratio (and \( \epsilon\) can thus be larger than eps in general). In any case, the number \( \epsilon\) (and with this, the achived approximation \( 1+\epsilon\)) can be queried by calling the routine achieved_epsilon() discussed below.

Template Parameters
Iteratormust be a model of InputIterator with Point as value type.
Precondition
The dimension \( d\) of the input points must be at least \( 2\), and \( \epsilon>0\).

Member Function Documentation

◆ achieved_epsilon()

template<typename Traits >
double CGAL::Approximate_min_ellipsoid_d< Traits >::achieved_epsilon ( ) const

returns a number \( \epsilon'\) such that the computed approximation is (under exact arithmetic) guaranteed to be an \( (1+\epsilon')\)-approximation to \( (P)\).

Precondition
ame.is_full_dimensional() == true.
Postcondition
\( \epsilon'>0\).

◆ axes_lengths_begin()

template<typename Traits >
Axes_lengths_iterator CGAL::Approximate_min_ellipsoid_d< Traits >::axes_lengths_begin ( )

returns an iterator pointing to the first of the \( d\) descendantly sorted lengths of the computed ellipsoid's axes.

The \( d\) returned numbers are floating-point approximations to the exact axes-lengths of the computed ellipsoid; no guarantee is given w.r.t. the involved relative error. (See also method axes_direction_cartesian_begin().)

Precondition
ame.is_full_dimensional() == true, and \( d\in\{2,3\}\).

◆ axes_lengths_end()

template<typename Traits >
Axes_lengths_iterator CGAL::Approximate_min_ellipsoid_d< Traits >::axes_lengths_end ( )

returns the past-the-end iterator corresponding to axes_lengths_begin().

Precondition
ame.is_full_dimensional() == true, and \( d\in\{2,3\}\).

◆ axis_direction_cartesian_begin()

template<typename Traits >
Axes_direction_coordinate_iterator CGAL::Approximate_min_ellipsoid_d< Traits >::axis_direction_cartesian_begin ( int  i)

returns an iterator pointing to the first of the \( d\) Cartesian coordinates of the computed ellipsoid's \( i\)th axis direction (i.e., unit vector in direction of the ellipsoid's \( i\)th axis).

The direction described by this iterator is a floating-point approximation to the exact axis direction of the computed ellipsoid; no guarantee is given w.r.t. the involved relative error. An approximation to the length of axis \( i\) is given by the \( i\)th entry of axes_lengths_begin().

Precondition
ame.is_full_dimensional() == true, and \( d\in\{2,3\}\), and \( 0\leq i < d\).

◆ axis_direction_cartesian_end()

template<typename Traits >
Axes_direction_coordinate_iterator CGAL::Approximate_min_ellipsoid_d< Traits >::axis_direction_cartesian_end ( int  i)

returns the past-the-end iterator corresponding to axis_direction_cartesian_begin().

Precondition
ame.is_full_dimensional() == true, and \( d\in\{2,3\}\), and \( 0\leq i < d\).

◆ center_cartesian_begin()

template<typename Traits >
Center_coordinate_iterator CGAL::Approximate_min_ellipsoid_d< Traits >::center_cartesian_begin ( )

returns an iterator pointing to the first of the \( d\) Cartesian coordinates of the computed ellipsoid's center.

The returned point is a floating-point approximation to the ellipsoid's exact center; no guarantee is given w.r.t. the involved relative error.

Precondition
ame.is_full_dimensional() == true.

◆ center_cartesian_end()

template<typename Traits >
Center_coordinate_iterator CGAL::Approximate_min_ellipsoid_d< Traits >::center_cartesian_end ( )

returns the past-the-end iterator corresponding to center_cartesian_begin().

Precondition
ame.is_full_dimensional() == true.

◆ defining_matrix()

template<typename Traits >
double CGAL::Approximate_min_ellipsoid_d< Traits >::defining_matrix ( int  i,
int  j 
) const

gives access to the \( (i,j)\)th entry of the matrix \( E\) in the representation \( \{ x\in\E^d \mid x^T E x + x^T e + \eta\leq0 \}\) of the computed approximation ellipsoid \( {\cal E}\).

The number returned by this routine is \( (1+\epsilon')(d+1)\,E_{ij}\), where \( \epsilon'\) is the number returned by achieved_epsilon().

Precondition
\( 0\leq i,j\leq d\), where \( d\) is the dimension of the points \( P\), and ame.is_full_dimensional() == true.

◆ defining_scalar()

template<typename Traits >
double CGAL::Approximate_min_ellipsoid_d< Traits >::defining_scalar ( ) const

gives access to the scalar \( \eta\) from the representation \( \{ x\in\E^d \mid x^T E x + x^T e + \eta\leq0 \}\) of the computed approximation ellipsoid \( {\cal E}\).

The number returned by this routine is \( (1+\epsilon')(d+1)\,(\eta+1)\), where \( \epsilon'\) is the number returned by achieved_epsilon().

Precondition
ame.is_full_dimensional() == true.

◆ defining_vector()

template<typename Traits >
double CGAL::Approximate_min_ellipsoid_d< Traits >::defining_vector ( int  i) const

gives access to the \( i\)th entry of the vector \( e\) in the representation \( \{ x\in\E^d \mid x^T E x + x^T e + \eta\leq0 \}\) of the computed approximation ellipsoid \( {\cal E}\).

The number returned by this routine is \( (1+\epsilon')(d+1)\,e_{i}\), where \( \epsilon'\) is the number returned by achieved_epsilon().

Precondition
\( 0\leq i\leq d\), where \( d\) is the dimension of the points \( P\), and ame.is_full_dimensional() == true.

◆ is_full_dimensional()

template<typename Traits >
bool CGAL::Approximate_min_ellipsoid_d< Traits >::is_full_dimensional ( ) const

returns whether \( P\) is full-dimensional or not, i.e., returns true if and only if \( P\) is full-dimensional.

Note: due to the limited precision in the algorithm's underlying arithmetic, the result of this method is not always correct. Rather, a return value of false means that the points \( P\) are contained in a "very thin" linear subspace of \( \E^d\), and as a consequence, the algorithm cannot compute an approximation. More precisely, a return value of false means that the points \( P\) are contained between two parallel hyperplanes in \( \E^d\) that are very close to each other (possibly at distance zero) - so close, that the algorithm could not compute an approximation ellipsoid. Similarly, a return value of true does not guarantee \( P\) to be full-dimensional; but there exists an input pointset \( P'\) such that the points \( P'\) and \( P\) have almost identical coordinates and \( P'\) is full-dimensional.

◆ is_valid()

template<typename Traits >
bool CGAL::Approximate_min_ellipsoid_d< Traits >::is_valid ( bool  verbose = false) const

returns true iff ame is valid according to the above definition.

If verbose is true, some messages concerning the performed checks are written to the standard error stream.

◆ write_eps()

template<typename Traits >
void CGAL::Approximate_min_ellipsoid_d< Traits >::write_eps ( const std::string &  name) const

Writes the points \( P\) and the computed approximation to \( (P)\) as an EPS-file under pathname name.

Precondition
The dimension of points \( P\) must be \( 2\). Note: this routine is provided as a debugging routine; future version of CGAL might not provide it anymore.
ame.is_full_dimensional() == true.