G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsFitting< T > Class Template Reference

Detailed Description

template<class T>
class gismo::gsFitting< T >

Class for performing a fit of a parametrized point cloud with a gsGeometry.

+ Inheritance diagram for gsFitting< T >:
+ Collaboration diagram for gsFitting< T >:

Public Types

enum  tdm_method
 choose the method for computing the coefficients: TDM, PDM, HDM with different blending weights
 

Public Member Functions

void applySmoothing (T lambda, gsSparseMatrix< T > &A_mat)
 Adds to the matrix A_mat terms for minimization of second derivative, weighted with parameter lambda.
 
void assembleSystem (const gsMatrix< T > &points_int, const gsMatrix< T > &params_int, const gsMatrix< T > &points_bdy, const gsMatrix< T > &params_bdy, const index_t &num_basis, const gsSparseMatrix< T > &NNT, gsSparseMatrix< T > &A_tilde, gsMatrix< T > &rhs)
 assembleSystem: assembles the linear system for the Hybrid Distance Minimization method
 
void assembleSystem (gsSparseMatrix< T > &A_mat, gsMatrix< T > &B)
 Assembles system for the least square fit.
 
void blending_weights (const gsSparseMatrix< T > &N_int, const index_t &num_int, const T &mu, const T &sigma, const gsMatrix< T > &params_int, tdm_method method, gsSparseMatrix< T > &NNT)
 blending_weights: computes the blending weights mu and sigma for the balance mu * PDM + sigma * TDM in the HDM method
 
void compute (T lambda=0)
 compute: Computes the coefficients of the spline geometry via penalized least squares
 
void compute_normals (const index_t &num_int, const gsMatrix< T > &params_int, gsSparseMatrix< T > &N_int)
 compute_normals: Computes the normals of the fitted geometry at the input parameter values params_int
 
void compute_tdm (T lambda, T mu, T sigma, const std::vector< index_t > &interpIdx, tdm_method method=hybrid_curvature_pdm_tdm_boundary_pdm)
 compute_hdm: computes the coefficients of the spline geometry via Hybrid Distance Minimization (HDM)
 
void computeApproxError (T &error, int type=0) const
 
void computeErrors ()
 
std::vector< T > computeErrors (const gsMatrix<> &param_values, const gsMatrix<> &points)
 
void computeMaxNormErrors ()
 
gsMatrix< T > fill_pointWiseErrors (const index_t &num_int, T &max_err_int)
 vector of length num_int containing all the point-wise errors; store also the max err value in max_err_int
 
void get_Error (std::vector< T > &errors, int type=0) const
 
const gsBasis< T > & getBasis () const
 Returns the basis of the approximation.
 
 gsFitting ()
 default constructor
 
 gsFitting (gsMatrix< T > const &param_values, gsMatrix< T > const &points, gsBasis< T > &basis)
 gsFitting: Main constructor of the fitting class
 
 gsFitting (gsMatrix< T > const &param_values, gsMatrix< T > const &points, gsVector< index_t > offset, gsMappedBasis< 2, T > &mbasis)
 gsFitting: Main constructor of the fitting class for multi-patch
 
void initializeGeometry (const gsMatrix< T > &coefficients, const gsMatrix< T > &parameters)
 initializeGeometry: Initializes the fitted geometry with given coefficients and parameters
 
void initParametricDomain ()
 Initialize the parametric domain of the point cloud.
 
gsMatrix< T > inverse_principal_curvatures (const index_t &num_int, const gsMatrix< T > &params_int)
 vector of length num_int containing rho = 1/max(c1, c2), where c1, c2 are the principal curvature values computed at every parametric point params_int
 
bool is_corner (gsMatrix< T > &parametric_domain, gsVector< T > &parameter)
 check if the given parameter parameter is a corner of the domain parametric_domain
 
bool is_point_inside_support (const gsMatrix< T > &parameter, const gsMatrix< T > &support)
 
bool is_point_inside_support (const T x, const T y, const gsMatrix< T > &support)
 check if the given parameter x, y is inside the support support; same as is_point_inside_support, but different input format
 
bool is_point_within_cell (const gsMatrix< T > &parameter, const gsMatrix< T > &element)
 check if the given parameter parameter is within the cell element
 
bool is_point_within_cell (const T x, const T y, const gsMatrix< T > &element)
 check if the given parameter x, y is within the cell element; same as is_point_within_cell, but different input format
 
void iterativeCompute (T const &tolerance, unsigned const &num_iters=10)
 Computes the least squares fit for a gsBasis.
 
lambda () const
 Returns the smoothing weight used in the last fitting.
 
maxPointError () const
 Returns the maximum point-wise error from the pount cloud (or zero if not fitted)
 
minPointError () const
 Returns the minimum point-wise error from the pount cloud (or zero if not fitted)
 
const gsMappedSpline< 2, T > & mresult () const
 gives back the computed approximation for multipatch geometry
 
size_t numPointsBelow (T threshold) const
 Computes the number of points below the error threshold (or zero if not fitted)
 
void parameterCorrection (T accuracy=1e-8, index_t maxIter=10, T tolOrth=1e-6)
 parameterCorrection: globally apply maxIter steps of parameter correction to the least squares fitted geometry
 
void parameterCorrectionSepBoundary_pdm (T accuracy, index_t maxIter, const std::vector< index_t > &sepIndex)
 parameterCorrectionSepBoundary_pdm: apply maxIter steps of parameter correction for PDM method, separating interior and boundary points
 
void parameterCorrectionSepBoundary_tdm (T accuracy, index_t maxIter, T mu, T sigma, const std::vector< index_t > &sepIndex, tdm_method method=hybrid_curvature_pdm_tdm_boundary_pdm)
 parameterCorrectionSepBoundary_tdm: apply maxIter steps of parameter correction for HDM method, separating interior and boundary points
 
void parameterProjectionSepBoundary (T accuracy, const std::vector< index_t > &interpIdx)
 parameterProjectionSepBoundary: project the points onto the fitted geometry, separating interior and boundary points
 
const std::vector< T > & pointWiseErrors () const
 Return the errors for each point.
 
gsMatrix< T > pointWiseErrors (const gsMatrix<> &parameters, const gsMatrix<> &points)
 
gsMatrix< T > principal_curvatures (const gsMatrix< T > &params)
 compute the principal curvatures (c1, c2) at the given parameters params
 
gsGeometry< T > * result () const
 gives back the computed approximation
 
gsMatrix< T > & returnParamValues ()
 returns the parameter values
 
gsMatrix< T > returnPoints () const
 returns the points
 
void setConstraints (const gsSparseMatrix< T > &lhs, const gsMatrix< T > &rhs)
 
void setConstraints (const std::vector< boxSide > &fixedSides)
 
void setConstraints (const std::vector< boxSide > &fixedSides, const std::vector< gsBSpline< T > > &fixedCurves)
 
void setConstraints (const std::vector< index_t > &indices, const std::vector< gsMatrix< T > > &coefs)
 
void updateGeometry (gsMatrix< T > coefficients, gsMatrix< T > parameters)
 updateGeometry: Updates the fitted geometry with new coefficients and parameters
 
virtual ~gsFitting ()
 Destructor.
 

Protected Member Functions

void assembleBlockB (const gsMatrix< T > &points, const gsMatrix< T > &params, index_t num_basis, gsSparseMatrix< T > &result) const
 Assembles 3xblock collocation matrix.
 
void assembleBlockX (const gsMatrix< T > &points, gsMatrix< T > &result) const
 Assembles the right hand side vectors for PDM/TDM.
 

Protected Attributes

gsFunctionSet< T > * m_basis
 Pointer keeping the basis.
 
gsSparseMatrix< T > m_constraintsLHS
 
gsMatrix< T > m_constraintsRHS
 
m_max_error
 Maximum point-wise error.
 
m_min_error
 Minimum point-wise error.
 
gsMappedSpline< 2, T > m_mresult
 Pointer keeping the resulting multipatch geometry.
 
gsMatrix< T > m_param_values
 the parameter values of the point cloud
 
gsMatrix< T > m_points
 the points of the point cloud
 
gsGeometry< T > * m_result
 Pointer keeping the resulting geometry.
 

Private Member Functions

void extendSystem (gsSparseMatrix< T > &A_mat, gsMatrix< T > &m_B)
 Extends the system of equations by taking constraints into account.
 

Constructor & Destructor Documentation

◆ gsFitting()

template<class T >
gsFitting ( gsMatrix< T > const &  param_values,
gsMatrix< T > const &  points,
gsBasis< T > &  basis 
)

gsFitting: Main constructor of the fitting class

Parameters
param_valuesa matrix containing the parameter values that parametrize the points
pointsmatrix containing the points to be fitted
basisbasis to use for fitting

Member Function Documentation

◆ assembleSystem()

template<class T >
void assembleSystem ( const gsMatrix< T > &  points_int,
const gsMatrix< T > &  params_int,
const gsMatrix< T > &  points_bdy,
const gsMatrix< T > &  params_bdy,
const index_t num_basis,
const gsSparseMatrix< T > &  NNT,
gsSparseMatrix< T > &  A_tilde,
gsMatrix< T > &  rhs 
)

assembleSystem: assembles the linear system for the Hybrid Distance Minimization method

Parameters
points_intinterior points
params_intinterior parameters
points_bdyboundary points
params_bdyboundary parameters
num_basisdimension of the basis
NNTmatrix containing the normals and the blending weights
A_tildeoutput system matrix
rhsoutput right-hand side vector

◆ blending_weights()

template<class T >
void blending_weights ( const gsSparseMatrix< T > &  N_int,
const index_t num_int,
const T &  mu,
const T &  sigma,
const gsMatrix< T > &  params_int,
tdm_method  method,
gsSparseMatrix< T > &  NNT 
)

blending_weights: computes the blending weights mu and sigma for the balance mu * PDM + sigma * TDM in the HDM method

Parameters
N_intmatrix containing the nomals at the parameters params_int
num_intindeces of the interior parameters
muweight for PDM
sigmaweight for TDM
params_intinput parameter values
methodmethod for computing the blending weights: constant, based on point-wise error, based on curvature
NNToutput matrix containing the normals and the blending weights

◆ compute()

template<class T >
void compute ( lambda = 0)

compute: Computes the coefficients of the spline geometry via penalized least squares

Parameters
lambdasmoothing weight

◆ compute_normals()

template<class T >
void compute_normals ( const index_t num_int,
const gsMatrix< T > &  params_int,
gsSparseMatrix< T > &  N_int 
)

compute_normals: Computes the normals of the fitted geometry at the input parameter values params_int

Parameters
num_intindex of the input parameter values
params_intinput parameter values
N_intmatrix containing the normals

◆ compute_tdm()

template<class T >
void compute_tdm ( lambda,
mu,
sigma,
const std::vector< index_t > &  interpIdx,
tdm_method  method = hybrid_curvature_pdm_tdm_boundary_pdm 
)

compute_hdm: computes the coefficients of the spline geometry via Hybrid Distance Minimization (HDM)

Parameters
lambdasmoothing weight
muweight for the PDM
sigmaweight for the TDM
interpIdxvector containing the number of interior points and the indices of the boundary points
methodmethod for computing the blending weights

◆ computeApproxError()

template<class T >
void computeApproxError ( T &  error,
int  type = 0 
) const

Computes the approximation error error of the fitted geometry to the original point cloud type = 0: sum of squares, type = 1: sum of absolute values (l_1 norm)

◆ computeErrors() [1/2]

template<class T >
void computeErrors ( )

Computes the point-wise errors in euclidean norm as well as the max and min errors, and updates the member variables m_pointErrors, m_max_error, m_min_error; different from computeMaxNormErrors(), where the error is computed in inifinity/maximum norm

◆ computeErrors() [2/2]

template<class T >
std::vector< T > computeErrors ( const gsMatrix<> &  param_values,
const gsMatrix<> &  points 
)

Computes min, max and mse errors in euclidean norms between the fitted geometry at the parameters param_values and the input point cloud points it does not update the member variables

◆ computeMaxNormErrors()

template<class T >
void computeMaxNormErrors ( )

Computes the point-wise errors in infinity/maximum norm as well as the max and min errors, and updates the member variables m_pointErrors, m_max_error, m_min_error; different from computeErrors(), where the error is computed in euclidean norm

◆ get_Error()

template<class T >
void get_Error ( std::vector< T > &  errors,
int  type = 0 
) const

Compute the point-wise errors for each point type = 0: point-wise infinity/maximum norm

◆ initializeGeometry()

template<class T >
void initializeGeometry ( const gsMatrix< T > &  coefficients,
const gsMatrix< T > &  parameters 
)

initializeGeometry: Initializes the fitted geometry with given coefficients and parameters

Parameters
coefficientsthe input coefficients
parametersthe input parameters

◆ is_point_inside_support()

template<class T >
bool is_point_inside_support ( const gsMatrix< T > &  parameter,
const gsMatrix< T > &  support 
)

check if the given parameter parameter is inside the support support difference with is_point_inside_cell in the inclusion of the left and right interval extremes.

◆ parameterCorrection()

template<class T >
void parameterCorrection ( accuracy = 1e-8,
index_t  maxIter = 10,
tolOrth = 1e-6 
)

parameterCorrection: globally apply maxIter steps of parameter correction to the least squares fitted geometry

Parameters
accuracyaccuracy of the closest point computation
maxItermaximum number of parameter correction steps
tolOrthorthogonality tolerance

◆ parameterCorrectionSepBoundary_pdm()

template<class T >
void parameterCorrectionSepBoundary_pdm ( accuracy,
index_t  maxIter,
const std::vector< index_t > &  sepIndex 
)

parameterCorrectionSepBoundary_pdm: apply maxIter steps of parameter correction for PDM method, separating interior and boundary points

Parameters
accuracyaccuracy of the closest point computation
maxItermaximum number of parameter correction steps
sepIndexvector containing the number of interior points and the indices of the boundary points

◆ parameterCorrectionSepBoundary_tdm()

template<class T >
void parameterCorrectionSepBoundary_tdm ( accuracy,
index_t  maxIter,
mu,
sigma,
const std::vector< index_t > &  sepIndex,
tdm_method  method = hybrid_curvature_pdm_tdm_boundary_pdm 
)

parameterCorrectionSepBoundary_tdm: apply maxIter steps of parameter correction for HDM method, separating interior and boundary points

Parameters
accuracyaccuracy of the closest point computation
maxItermaximum number of parameter correction steps
muweight for PDM
sigmaweight for TDM
sepIndexvector containing the number of interior points and the indices of the boundary points
methodmethod for computing the blending weights

◆ parameterProjectionSepBoundary()

template<class T >
void parameterProjectionSepBoundary ( accuracy,
const std::vector< index_t > &  interpIdx 
)

parameterProjectionSepBoundary: project the points onto the fitted geometry, separating interior and boundary points

Parameters
accuracyaccuracy of the closest point computation, for the foot-point projection
interpIdxvector containing the number of interior points and the indices of the boundary points

◆ pointWiseErrors()

template<class T >
gsMatrix< T > pointWiseErrors ( const gsMatrix<> &  parameters,
const gsMatrix<> &  points 
)

Compute point-wise error in euclidean norm between the fitted geometry at the parameters parameters and the input point cloud points similar to computeErrors(), but different input and output format; it does not update the member variables

◆ setConstraints() [1/4]

template<class T >
void setConstraints ( const gsSparseMatrix< T > &  lhs,
const gsMatrix< T > &  rhs 
)
inline

Sets constraints that the coefficients of the resulting geometry have to conform to. More precisely, denoting the coefficient vector by x, it enforces lhs * x = rhs.

◆ setConstraints() [2/4]

template<class T >
void setConstraints ( const std::vector< boxSide > &  fixedSides)

Sets constraints in such a way that the previous values at fixedSides of the geometry remain intact.

◆ setConstraints() [3/4]

template<class T >
void setConstraints ( const std::vector< boxSide > &  fixedSides,
const std::vector< gsBSpline< T > > &  fixedCurves 
)

Set constraints in such a way that the resulting geometry on each of fixedSides will coincide with the corresponding curve in fixedCurves.

◆ setConstraints() [4/4]

template<class T >
void setConstraints ( const std::vector< index_t > &  indices,
const std::vector< gsMatrix< T > > &  coefs 
)

Sets constraints on that the coefficients of the resulting geometry have to conform to.

Parameters
indicesindices (in the coefficient vector) of the prescribed coefficients.
coefsprescribed coefficients.

◆ updateGeometry()

template<class T >
void updateGeometry ( gsMatrix< T >  coefficients,
gsMatrix< T >  parameters 
)

updateGeometry: Updates the fitted geometry with new coefficients and parameters

Parameters
coefficientsthe new coefficients
parametersthe new parameters

Member Data Documentation

◆ m_constraintsLHS

template<class T >
gsSparseMatrix<T> m_constraintsLHS
protected

Left hand-side of the constraints that the coefficients of the resulting geometry have to conform to. This corresponds to matrix D in Prautzch, Boehm, Paluszny: Bezier and B-spline techniques, Section 4.7.

◆ m_constraintsRHS

template<class T >
gsMatrix<T> m_constraintsRHS
protected

Right hand-side of the constraints that the coefficients of the resulting geometry have to conform to. This corresponds to vector q in Prautzch, Boehm, Paluszny: Bezier and B-spline techniques, Section 4.7.