iRoCS Toolbox  1.1.0
ATBSpline.hh
Go to the documentation of this file.
1 /**************************************************************************
2  *
3  * Copyright (C) 2015 Thorsten Falk
4  *
5  * Image Analysis Lab, University of Freiburg, Germany
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 3 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software Foundation,
19  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
20  *
21  **************************************************************************/
22 
23 #ifndef ATBSPLINE_HH
24 #define ATBSPLINE_HH
25 
26 #ifdef HAVE_CONFIG_H
27 #include <config.hh>
28 #endif
29 
30 #include <cstdlib>
31 #include <vector>
32 #include <map>
33 #include <omp.h>
34 
35 #include <blitz/array.h>
36 
38 
39 #include "ATBPolynomial.hh"
40 #include "RuntimeError.hh"
41 #include "ATBLinAlg.hh"
42 
43 // #define ATBSPLINE_DEBUG
44 // #define ATBSPLINE_PROFILE
45 
46 #ifdef ATBSPLINE_PROFILE
47 #include "ATBTiming.hh"
48 #endif
49 
50 namespace atb
51 {
52 
53 /*======================================================================*/
66 /*======================================================================*/
67  template<typename ControlPointT>
68  class BSpline
69  {
70 
71  public:
72 
73 /*======================================================================*/
79 /*======================================================================*/
81  {
82  BaseCacheIndex(double u, size_t index, int degree, int derivative);
83  BaseCacheIndex(BaseCacheIndex const &idx);
84  bool operator<(BaseCacheIndex const &i) const;
85 
86  double u;
87  size_t j;
88  int p, l;
89  };
90 
91 /*======================================================================*/
98 /*======================================================================*/
100  {
102  size_t segment, size_t index, int degree, int derivative);
104  bool operator<(BasePolynomialCacheIndex const &i) const;
105 
106  size_t m, j;
107  int p, l;
108  };
109 
110 /*======================================================================*/
120 /*======================================================================*/
121  BSpline(int degree = 3, size_t nControlPoints = 0);
122 
123 /*======================================================================*/
134 /*======================================================================*/
135  BSpline(int degree, std::vector<ControlPointT> const &controlPoints);
136 
137 /*======================================================================*/
145 /*======================================================================*/
146  BSpline(BSpline<ControlPointT> const &spline);
147 
148 /*======================================================================*/
152 /*======================================================================*/
153  ~BSpline();
154 
155 /*======================================================================*/
163 /*======================================================================*/
165 
166 /*======================================================================*/
172 /*======================================================================*/
173  int degree() const;
174 
175 /*======================================================================*/
182 /*======================================================================*/
183  void setDegree(int degree);
184 
185 /*======================================================================*/
191 /*======================================================================*/
192  size_t nControlPoints() const;
193 
194 /*======================================================================*/
202 /*======================================================================*/
203  void setNControlPoints(size_t nControlPoints);
204 
205 /*======================================================================*/
214 /*======================================================================*/
215  void setControlPoints(std::vector<ControlPointT> const &controlPoints);
216 
217 /*======================================================================*/
223 /*======================================================================*/
224  std::vector<ControlPointT> const &controlPoints() const;
225 
226 /*======================================================================*/
235 /*======================================================================*/
236  ControlPointT const &controlPoint(size_t index) const;
237 
238 /*======================================================================*/
245 /*======================================================================*/
246  void setControlPoint(size_t index, ControlPointT const &point);
247 
248 /*======================================================================*/
256 /*======================================================================*/
257  size_t nKnots() const;
258 
259 /*======================================================================*/
273 /*======================================================================*/
274  void setKnots(std::vector<double> const &knots);
275 
276 /*======================================================================*/
284 /*======================================================================*/
285  double knot(size_t index) const;
286 
287 /*======================================================================*/
293 /*======================================================================*/
294  std::vector<double> const &knots() const;
295 
296 /*======================================================================*/
308 /*======================================================================*/
309  void setKnot(size_t index, double value);
310 
311 /*======================================================================*/
319 /*======================================================================*/
320  void setOpenUniformKnots(double uMin = 0.0, double uMax = 1.0);
321 
322 /*======================================================================*/
339 /*======================================================================*/
340  double base(double u, size_t index, int derivative = 0) const;
341 
342 /*======================================================================*/
353 /*======================================================================*/
354  double base(BaseCacheIndex const &i) const;
355 
356 /*======================================================================*/
378 /*======================================================================*/
380  size_t segment, size_t index, int derivative = 0) const;
381 
382 /*======================================================================*/
393 /*======================================================================*/
395  BasePolynomialCacheIndex const &idx) const;
396 
397 /*======================================================================*/
408 /*======================================================================*/
409  std::string printBasePolynomial(BasePolynomialCacheIndex idx) const;
410 
411 /*======================================================================*/
425 /*======================================================================*/
426  size_t curvePositionToSegment(double u) const;
427 
428 /*======================================================================*/
448 /*======================================================================*/
449  size_t curvePositionToSegment(double u, double tolerance) const;
450 
451 /*======================================================================*/
459 /*======================================================================*/
460  ControlPointT operator()(double u) const;
461 
462 /*======================================================================*/
484 /*======================================================================*/
485  ControlPointT operator()(double u, double lBound, double uBound) const;
486 
487 /*======================================================================*/
496 /*======================================================================*/
497  ControlPointT derivative(double u, size_t derivative = 1) const;
498 
499 /*======================================================================*/
508 /*======================================================================*/
509  BSpline<ControlPointT> derivative(size_t derivative = 1) const;
510 
511 /*======================================================================*/
528 /*======================================================================*/
529  double curveIntegral(double uFrom, double uTo, double uStep = 0.001) const;
530 
531 /*======================================================================*/
542 /*======================================================================*/
543  double curveLengthToU(double origin, double distance) const;
544 
545 /*======================================================================*/
556 /*======================================================================*/
557  double extendedCurveLengthToU(double origin, double distance) const;
558 
559 /*======================================================================*/
576 /*======================================================================*/
577  ControlPointT extendedDerivative(
578  double u, size_t derivative = 1, double lBound = 0.0,
579  double uBound = 1.0) const;
580 
581 /*======================================================================*/
602 /*======================================================================*/
603  double extendedCurveIntegral(
604  double uFrom, double uTo, double lBound = 0.0, double uBound = 1.0,
605  double uStep = 0.001) const;
606 
607 /*======================================================================*/
624 /*======================================================================*/
625  void save(std::string const &fileName, std::string const &groupName,
626  bool throwExceptions = false) const;
627 
628 /*======================================================================*/
643 /*======================================================================*/
644  void load(std::string const &fileName, std::string const &groupName,
645  bool throwExceptions = false);
646 
647 /*======================================================================*/
658 /*======================================================================*/
659  void save(BlitzH5File &outFile, std::string const &groupName) const;
660 
661 /*======================================================================*/
671 /*======================================================================*/
672  void load(BlitzH5File const &inFile, std::string const &groupName);
673 
674  private:
675 
676 /*======================================================================*/
684 /*======================================================================*/
685  void _setNKnots(size_t nKnots);
686 
687 /*======================================================================*/
691 /*======================================================================*/
692  void _clearCache() const;
693 
694 /*======================================================================*/
698 /*======================================================================*/
699  void _copyCache(BSpline<ControlPointT> const &spline) const;
700 
701 /*======================================================================*/
708 /*======================================================================*/
709  void _updateCache() const;
710 
711  int _degree;
712  std::vector<ControlPointT> _controlPoints;
713  mutable std::vector<double> _knots;
714 
715  mutable std::vector<
716  std::vector< std::vector< std::vector<Polynomial<double>* > > > >
717  _basePolynomialDerivativeCache;
718  mutable std::vector<
719  std::vector< std::vector< std::vector<Polynomial<double>* > > > >
720  _basePolynomialIntegralCache;
721 
722  static Polynomial<double> _zeroPolynomial;
723 
724  friend void fitSplineToSpline(
725  BSpline<double> &spline, BSpline<double> const &reference);
726 
727  template<int Dim>
728  friend void fitSplineToSpline(
729  BSpline< blitz::TinyVector<double,Dim> > &spline,
730  BSpline< blitz::TinyVector<double,Dim> > const &reference);
731 
732  friend void fitSplineToPointCloud(
733  BSpline<double> &spline, const std::vector<double> &u,
734  std::vector<double> const &points);
735 
736  template<int Dim>
737  friend void fitSplineToPointCloud(
738  BSpline< blitz::TinyVector<double,Dim> > &spline,
739  std::vector<double> const &u,
740  std::vector< blitz::TinyVector<double,Dim> > const &points);
741 
742  friend double distance(BSpline<double> const &spline, double x, double &u);
743 
744  template<int Dim>
745  friend double distance(
746  BSpline< blitz::TinyVector<double,Dim> > const &spline,
747  blitz::TinyVector<double,Dim> const &x, double &u);
748 
749  template<int Dim>
750  friend double extendedDistance(
751  BSpline< blitz::TinyVector<double,Dim> > const &spline,
752  blitz::TinyVector<double,Dim> const &x, double &u,
753  double lBound, double uBound);
754 
755 };
756 
757 /*======================================================================*/
767 /*======================================================================*/
768  template<typename ControlPointT>
769  std::ostream &operator<<(
770  std::ostream &os, BSpline<ControlPointT> const &spline);
771 
772 /*======================================================================*/
781 /*======================================================================*/
782  void fitSplineToSpline(
783  BSpline<double> &spline, BSpline<double> const &reference);
784 
785 /*======================================================================*/
794 /*======================================================================*/
795  template<int Dim>
796  void fitSplineToSpline(
797  BSpline< blitz::TinyVector<double,Dim> > &spline,
798  BSpline< blitz::TinyVector<double,Dim> > const &reference);
799 
800 /*======================================================================*/
817 /*======================================================================*/
819  BSpline<double> &spline, std::vector<double> const &u,
820  std::vector<double> const &points);
821 
822 /*======================================================================*/
839 /*======================================================================*/
840  template<int Dim>
842  BSpline< blitz::TinyVector<double,Dim> > &spline,
843  std::vector<double> const &u,
844  std::vector< blitz::TinyVector<double,Dim> > const &points);
845 
846 /*======================================================================*/
858 /*======================================================================*/
859  double distance(BSpline<double> const &spline, double x, double &u);
860 
861 /*======================================================================*/
873 /*======================================================================*/
874  template<int Dim>
875  double distance(
876  BSpline< blitz::TinyVector<double,Dim> > const &spline,
877  blitz::TinyVector<double,Dim> const &x, double &u);
878 
879 /*======================================================================*/
896 /*======================================================================*/
897  double extendedDistance(
898  BSpline<double> const &spline, double x, double &u,
899  double lBound = 0.0, double uBound = 1.0);
900 
901 /*======================================================================*/
918 /*======================================================================*/
919  template<int Dim>
920  double extendedDistance(
921  BSpline< blitz::TinyVector<double,Dim> > const &spline,
922  blitz::TinyVector<double,Dim> const &x, double &u,
923  double lBound = 0.0, double uBound = 1.0);
924 
925 }
926 
927 #include "ATBSpline.icc"
928 
929 #endif
ATBPolynomial.hh provides the ATB::Polynomial class to store and do arithmetics within the Polynomial...
BSpline< ControlPointT > & operator=(BSpline< ControlPointT > const &spline)
Copy assignment operator.
The BSpline class provides functions for fitting B-Splines to point clouds and evaluating them at arb...
Definition: ATBSpline.hh:68
Exception specialization for error handling within libArrayToolbox.
std::ostream & operator<<(std::ostream &os, BasicTreeNode< KeyT, ContentT > const &n)
std::vector< double > const & knots() const
Read accessor to the knot vector of the B-Spline.
ControlPointT derivative(double u, size_t derivative=1) const
Computation of the derivative of the spline function at position u.
friend double extendedDistance(BSpline< blitz::TinyVector< double, Dim > > const &spline, blitz::TinyVector< double, Dim > const &x, double &u, double lBound, double uBound)
Compute the distance and the corresponding curve position of the given point to the spline...
void setDegree(int degree)
Set the degree of the B-Spline.
friend double distance(BSpline< double > const &spline, double x, double &u)
Compute the distance and the corresponding curve position of the given point to the spline...
void setNControlPoints(size_t nControlPoints)
Set the number of control points of the B-Spline.
The BaseCacheIndex struct provides a sortable quadrupel for uniquely identifying the evaluation of th...
Definition: ATBSpline.hh:80
double extendedCurveLengthToU(double origin, double distance) const
Get the curve parameter u at the specified distance from the specified origin.
size_t curvePositionToSegment(double u) const
Compute the index of the spline segment, that corresponds to the curve position u.
Polynomial< double > const & basePolynomial(size_t segment, size_t index, int derivative=0) const
Get the polynomial representation of the (ell&#39;th derivative/integral) of the j&#39;th basis function on t...
Lightweight alternative to libBlitzHDF5 providing its basic functionality.
void setKnot(size_t index, double value)
Set the value of the knot at the specified index.
size_t nControlPoints() const
Returns the number of control points of the B-Spline.
std::vector< ControlPointT > const & controlPoints() const
Read accessor to the control point vector of the B-Spline.
BSpline(int degree=3, size_t nControlPoints=0)
Default constructor.
void save(std::string const &fileName, std::string const &groupName, bool throwExceptions=false) const
Save the spline with corresponding meta-data to an hdf5 file with given name, using the groupName pro...
size_t nKnots() const
Returns the number of knots of the B-Spline.
BaseCacheIndex(double u, size_t index, int degree, int derivative)
double knot(size_t index) const
Read accessor to the k&#39;th knot of the B-Spline.
ATBTiming.hh provides the MyDateTime and Timer classes for high accuracy profiling.
The BasePolynomialCacheIndex struct provides a sortable quadrupel to uniquely identify the polynomial...
Definition: ATBSpline.hh:99
friend void fitSplineToSpline(BSpline< double > &spline, BSpline< double > const &reference)
Fitting of a spline curve to another spline curve.
ControlPointT operator()(double u) const
Evaluation of the spline curve at position u.
friend void fitSplineToPointCloud(BSpline< double > &spline, const std::vector< double > &u, std::vector< double > const &points)
Fitting of a spline curve to scalar data.
bool operator<(BaseCacheIndex const &i) const
~BSpline()
Destructor.
void setControlPoints(std::vector< ControlPointT > const &controlPoints)
Replace the current control points by the provided ones.
int degree() const
Returns the B-Spline degree.
ControlPointT const & controlPoint(size_t index) const
Read the jth control point.
void setKnots(std::vector< double > const &knots)
With this function the knot Array can be set manually.
std::string printBasePolynomial(BasePolynomialCacheIndex idx) const
Get a string representation of the (ell&#39;th derivative/integral) of the j&#39;th basis function polynomial...
void load(std::string const &fileName, std::string const &groupName, bool throwExceptions=false)
Load the spline with corresponding meta-data from an hdf5 file with given name, using the groupName p...
double base(double u, size_t index, int derivative=0) const
Evalation of the B-Spline basis function index at position u along the curve.
double curveLengthToU(double origin, double distance) const
Get the curve parameter u at the specified distance from the specified origin.
void setOpenUniformKnots(double uMin=0.0, double uMax=1.0)
Generation of an open uniform B-Spline knot configuration.
ControlPointT extendedDerivative(double u, size_t derivative=1, double lBound=0.0, double uBound=1.0) const
Computation of the derivative of the spline function at position u.
double extendedCurveIntegral(double uFrom, double uTo, double lBound=0.0, double uBound=1.0, double uStep=0.001) const
Estimate the curve length from curve position uFrom to curve position uTo.
double curveIntegral(double uFrom, double uTo, double uStep=0.001) const
Estimate the curve length from curve position uFrom to curve position uTo.
void setControlPoint(size_t index, ControlPointT const &point)
Set the ControlPoint at the specified index.