G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsCurvatureSmoothing.h
Go to the documentation of this file.
1 
19 #pragma once
20 
21 #include <iostream>
22 #include <stdexcept>
23 
24 #include <gsCore/gsGeometry.h>
25 #include <gsNurbs/gsBSpline.h>
26 #include <gsCore/gsLinearAlgebra.h>
27 
28 namespace gismo
29 {
39 template<class T>
41 {
42 public:
45  {
46  m_curve_original = NULL;
47  m_curve_smooth = NULL;
48  }
49 
51  gsCurvatureSmoothing(const gsBSpline<T> & init_curve, const gsMatrix<T>& param_values, const gsMatrix<T>& points)
52  {
53  m_curve_original = &init_curve;
54  m_curve_smooth = init_curve.clone().release();
55  m_param_values = param_values;
56  m_points = points;
57  }
58 
61  {
62  delete m_curve_smooth; // deleting the B-spline curve
63  }
64 
65 
67  const gsBSpline<T> & curveOriginal() const { return *m_curve_original; }
69  const gsBSpline<T> & curveSmooth() const { return *m_curve_smooth; }
70 
71 
74  void smoothTotalVariation(const T omega1, const T omega2, const T lamda, const T tau, const unsigned iter=50);
75 
78  void smoothTotalVariationSelectLamda(const T omega1, const T omega2, const gsMatrix<T> listlamdas, const unsigned iter=50);
79 
82  void smoothTotalVariationSelectLamda(const T omega1, const T omega2, const T lamda, const unsigned iter=50);
83 
85  void smoothHadenfeld(const unsigned smooth_degree, const T delta, const index_t iter_step, const index_t iter_total, gsVector<index_t> &iterated, const bool original=true);
86 
87 
91  void smoothAllHadenfeld(const unsigned smooth_degree=4, const unsigned iter=500);
92 
94  void write(std::ostream &os);
95 
97  void computeApproxError(T & error);
98 
100  void computeApproxErrorL2(T & error);
101 
103  void computeApproxErrorLMax(T & error);
104 
106  void computeApproxErrorCoef(T & error);
107 
109  void computeCurvatureError(T & error);
110 
111 private:
112 
121 
123  void compute_AllValues(gsBSplineBasis<T> * basis, gsMatrix<T> u, gsMatrix<T> *coefs, gsMatrix<T> & values0, gsMatrix<T> & values1, gsMatrix<T> & values2, gsMatrix<T> & values3);
124 
126  void compute_ObjectiveFunction(gsBSplineBasis<T> * basis, gsMatrix<T> *coefs, const T omega1, const T omega2, T &value);
127 
129  // TODO: What is exactly the purpose of this function; why do we want to change output?
130  void reset(gsBSpline<T> * newCurve)
131  {
132  if (m_curve_smooth)
133  delete m_curve_smooth;
134  m_curve_smooth = newCurve;
135  }
136 
137 }; // class gsCurvatureSmoothing
138 
139 } // namespace gismo
140 
141 #ifndef GISMO_BUILD_LIB
142 #include GISMO_HPP_HEADER(gsCurvatureSmoothing.hpp)
143 #endif
gsMatrix< T > m_points
the points of the original point cloud
Definition: gsCurvatureSmoothing.h:120
gsMatrix< T > m_param_values
the parameter values of the original point cloud
Definition: gsCurvatureSmoothing.h:118
void smoothTotalVariationSelectLamda(const T omega1, const T omega2, const gsMatrix< T > listlamdas, const unsigned iter=50)
Definition: gsCurvatureSmoothing.hpp:203
const gsBSpline< T > * m_curve_original
the original B-spline curve
Definition: gsCurvatureSmoothing.h:114
void reset(gsBSpline< T > *newCurve)
set the smooth curve to the the original curve
Definition: gsCurvatureSmoothing.h:130
void computeCurvatureError(T &error)
computes the curvature error of the smoother curve
Definition: gsCurvatureSmoothing.hpp:632
void computeApproxErrorCoef(T &error)
computes the max approximation error between the coeffeicientsof the original and the smoother curve ...
Definition: gsCurvatureSmoothing.hpp:616
Class for computing a closed B-spline curve with a smaller number of curvature extrema compared to a ...
Definition: gsCurvatureSmoothing.h:40
Provides declaration of Geometry abstract interface.
void compute_ObjectiveFunction(gsBSplineBasis< T > *basis, gsMatrix< T > *coefs, const T omega1, const T omega2, T &value)
computes the objective function for given coefs and omega1 and omega2 – objective function = omega1*A...
Definition: gsCurvatureSmoothing.hpp:686
#define index_t
Definition: gsConfig.h:32
const gsBSpline< T > & curveSmooth() const
gives back the smoother B-spline curve
Definition: gsCurvatureSmoothing.h:69
void computeApproxErrorL2(T &error)
computes the L^{2}-norm approximation error of the smoother curve
Definition: gsCurvatureSmoothing.hpp:591
A B-spline function of one argument, with arbitrary target dimension.
Definition: gsBSpline.h:50
void smoothTotalVariation(const T omega1, const T omega2, const T lamda, const T tau, const unsigned iter=50)
Definition: gsCurvatureSmoothing.hpp:27
A univariate B-spline basis.
Definition: gsBSplineBasis.h:694
const gsBSpline< T > & curveOriginal() const
gives back the original B-spline curve
Definition: gsCurvatureSmoothing.h:67
void computeApproxErrorLMax(T &error)
computes the L-max-norm approximation error of the smoother curve
Definition: gsCurvatureSmoothing.hpp:601
~gsCurvatureSmoothing()
Destructor.
Definition: gsCurvatureSmoothing.h:60
void smoothHadenfeld(const unsigned smooth_degree, const T delta, const index_t iter_step, const index_t iter_total, gsVector< index_t > &iterated, const bool original=true)
smooth the curve by smoothing only one cofficient in each step using the Hadenfeld algorithm — the us...
Definition: gsCurvatureSmoothing.hpp:388
void compute_AllValues(gsBSplineBasis< T > *basis, gsMatrix< T > u, gsMatrix< T > *coefs, gsMatrix< T > &values0, gsMatrix< T > &values1, gsMatrix< T > &values2, gsMatrix< T > &values3)
computes all values and derivatives (up to three) at the parameter values u for the given coefs ...
Definition: gsCurvatureSmoothing.hpp:650
uPtr clone()
Clone methode. Produceds a deep copy inside a uPtr.
void smoothAllHadenfeld(const unsigned smooth_degree=4, const unsigned iter=500)
Definition: gsCurvatureSmoothing.hpp:500
Represents a B-spline curve/function with one parameter.
This is the main header file that collects wrappers of Eigen for linear algebra.
void write(std::ostream &os)
writes the smooth curve to a file, which can be visualized in Mathematica (Name of Mathematica File V...
Definition: gsCurvatureSmoothing.hpp:565
gsCurvatureSmoothing()
default constructor
Definition: gsCurvatureSmoothing.h:44
gsCurvatureSmoothing(const gsBSpline< T > &init_curve, const gsMatrix< T > &param_values, const gsMatrix< T > &points)
constructor
Definition: gsCurvatureSmoothing.h:51
gsBSpline< T > * m_curve_smooth
the smoother B-spline curve
Definition: gsCurvatureSmoothing.h:116
void computeApproxError(T &error)
computes approximation error of the smoother curve to the original point cloud
Definition: gsCurvatureSmoothing.hpp:578