LSsystem Class Reference

Least squares approximation to scattered data with B-splines More...

#include <LSsystem.h>

List of all members.

Public Member Functions

 LSsystem (boost::shared_ptr< std::vector< double > > U, boost::shared_ptr< std::vector< double > > V, boost::shared_ptr< std::vector< double > > Zvals, int noX, int noY, boost::shared_ptr< GenMatrix< UCBspl_real > > PHI, double smoothingFac=1.0)
void setWeights (boost::shared_ptr< std::map< int, double, std::less< int > > > weights)
void buildEqSystem (bool buildRHS=true)
void relaxGaussSeidel (int noIterations)
void relaxCG (int noIterations)
void relaxCGPrecond (int noIterations)
void relaxJacobi (int noIterations, double omega=2./3.)
void setDomain (double umin, double vmin, double umax, double vmax)
void getDomain (double &umin, double &vmin, double &umax, double &vmax)
int numberOfUnknowns () const
void setSmoothingFactor (double smoothingFac)
bool addPoint (double u, double v, double z, double weight=1.0, bool addToRHS=true)
double currentSolutionNorm () const
double residual_l2 (bool scaled=false) const
double residual_linf () const

Static Public Member Functions

static double rowMatVecMult (int row_no, const MatSparse &A, const GenMatrix< UCBspl_real > &x)

Friends

class MGsystem


Detailed Description

Least squares approximation to scattered data with B-splines

LSsystem - The class takes a set of scattered data and produces a B-spline surface that approximates the scattered data. The approximation scheme is a least square fit with a thin plate spline smoothing term. A set of iterative equation solvers can be used to produce the set of B-spline coefficients.

No stopping criteria for the equation solvers are provided - the user control is only through the given number of iterations. Thus, we call it relaxation in the documentation. But the relaxation step may be repeated many times.

Note that the format of the scattered data and the solution vector is exactly the same as in the SINTEF MBA (Multilevel B-spline Approximation) Library. Thus, a good starting vector can be found by running MBA first. This is also necessary for surfaces with large coefficient grids (and unevenly distributed data), which otherwise will converge very slow.

This class can also be (and should be) considered as a relaxation step for the multigrid methods in class MGsystem. MGsystem contains a number of multigrid schemes which produces better results than LSsystem alone.

Examle of use (minimal):

   LSsystem lssys(xArr, yArr, zArr, noX, noY, PHI);  // initialize with scattered data and initial solution vector
   lssys.buildEqSystem(); // Build the equation system
   lssys.relaxCG(noIterations); // Run the given number of iterations with conjugate gradients
Author:
Øyvind Hjelle <Oyvind.Hjelle@math.sintef.no>


Constructor & Destructor Documentation

LSsystem::LSsystem ( boost::shared_ptr< std::vector< double > >  U,
boost::shared_ptr< std::vector< double > >  V,
boost::shared_ptr< std::vector< double > >  Zvals,
int  noX,
int  noY,
boost::shared_ptr< GenMatrix< UCBspl_real > >  PHI,
double  smoothingFac = 1.0 
)

Constructor with (standard) shared pointers to scattered data

Parameters:
U,: x-values
V,: y-values
Z,: z-values
noX,noY,: Number of coefficients in each direction of the spline surface
PHI,: The spline coefficients
smoothingFac,: Smoothing factor; a weight that determines the smoothness of the surface. See LSsystem::setSmoothingFactor.


Member Function Documentation

void LSsystem::buildEqSystem ( bool  buildRHS = true  ) 

Builds the system matrix and the right hand side (uses solution vector given by constructor)

void LSsystem::relaxGaussSeidel ( int  noIterations  ) 

Gauss-Seidel relaxation

void LSsystem::relaxCG ( int  noIterations  ) 

Conjugate Gradient relaxation

void LSsystem::relaxCGPrecond ( int  noIterations  ) 

Conjugate Gradient relaxation with multigrid preconditioning

void LSsystem::relaxJacobi ( int  noIterations,
double  omega = 2./3. 
)

Jacobi relaxation

void LSsystem::setDomain ( double  umin,
double  vmin,
double  umax,
double  vmax 
)

Set the domain over which the surface is to be defined. The default is the xy-range of the scattered data. If used, this must be done before creating the surface.

Note:
This function can only be used to expand the domain beyond the xy-range of the scattered data. It is the users responsibility to check that no scattered data falls outside the domain. (use std::min_element and std::max_element to find range of data for std::vector)
See also:
MGsystem::setDomain

void LSsystem::getDomain ( double &  umin,
double &  vmin,
double &  umax,
double &  vmax 
)

Get the domain over which the surface is defined

int LSsystem::numberOfUnknowns (  )  const [inline]

Number of unknowns in the equation system, i.e. B-spline coefficients.

void LSsystem::setSmoothingFactor ( double  smoothingFac  ) 

Adjust the thin plate spline energy. By default this value is 1.0. To get a surface that is smoother (but less accurate) than the default, a value grater than 1.0 should be given and vice versa. This will typically be done after relaxation has been run and one want a a surface that is smoother, or less smooth and more accurate, than obtained in from previous relaxation. Relaxation must then be run afterwards.

See also:
MGsystem::setSmoothingFactor

bool LSsystem::addPoint ( double  u,
double  v,
double  z,
double  weight = 1.0,
bool  addToRHS = true 
)

Add a point to the data set. This can be done after relaxation has been run as in the example above, but relaxation must be run again afterwards. A weight of importance different fro 1.0 can be given.

Return values:
bool false if (u,v) is outside the domain.

double LSsystem::currentSolutionNorm (  )  const

Calculate the l2 norm of the current solution (Frobenius)

double LSsystem::residual_l2 ( bool  scaled = false  )  const

The discrete l_2 norm of the residual ||b-Ax||_2 possibly scaled by area of grid cell

double LSsystem::residual_linf (  )  const

The max norm of the residual (unscaled) : ||b-Ax||_inf


The documentation for this class was generated from the following files:
Generated on Wed Nov 28 12:27:29 2007 for LSMG by  doxygen 1.5.1