MBA Class Reference

Main interaface for multilevel B-spline approximation More...

#include <MBA.h>

List of all members.

Public Member Functions

 MBA (boost::shared_ptr< dVec > U, boost::shared_ptr< dVec > V, boost::shared_ptr< dVec > Z)
void init (boost::shared_ptr< dVec > U, boost::shared_ptr< dVec > V, boost::shared_ptr< dVec > Z)
void init (UCBspl::SplineSurface &surf)
void setDomain (double umin, double vmin, double umax, double vmax)
void setBaseType (MBAbaseType baseType)
void setBaseValue (double base)
 Level of the base surface if MBA_CONSTVAL is used as base type.
void MBAalg (int m0, int n0, int h=0, int smoothing_steps=0)
const MBAdatagetData () const
 Get the data object (See class MBAdata).
void readScatteredData (char filename[])
 Read scattered data from file in the format u v z on each line.
void smoothZeros (int no_iter)
 Smooth the spline surface without affecting the approximation accuracy.
UCBspl::SplineSurface getSplineSurface () const
void getIndexDomain (int &m, int &n) const
void cleanup (int type=0)
boost::shared_ptr< GenMatrixTypePHI () const
void checkSparsity () const
void checkError () const

Static Public Member Functions

static void smoothMatrix (GenMatrixType &matrix, int no_iter)


Detailed Description

Main interaface for multilevel B-spline approximation

MBA - Multilevel B-spline approximation (and interpolation). This is the main interface to multilevel B-spline approximation of scattered data for functional surfaces. The method is based on B-spline refinement and produces one single B-spline surface, as opposed to MBAadaptive that produces a hierarchy of spline surfaces. The domain of the surface will be the rectangular domain of the given scattered data in the (x,y)-axes directions as default. (We use u and v for x and y in function headers) The function MBA::MBAalg produces the spline surface from the scattered data. The spline surface can be retrieved as a SplineSurface object with MBA::getSplineSurface. Note that this is a (mathematical) lower level interface; thus there is no check for the validity of given scattered data or of arguments passed to member functions.

Examle of use (minimal):

   MBA mba(x_arr, y_arr, z_arr);      // initialize with scattered data
   mba.MBAalg(1,1,7);                 // create spline surface

   UCBspl::SplineSurface surf = mba.getSplineSurface();
   double z = surf.f(x,y);            // evaluate (height of) surface in (x,y)
   surf.normalVector(x,y, nx,ny,nz);  // evaluate normal vector in (x,y)


Constructor & Destructor Documentation

MBA::MBA ( boost::shared_ptr< dVec U,
boost::shared_ptr< dVec V,
boost::shared_ptr< dVec Z 
) [inline]

Constructor with (standard) shared pointers to scattered data


Member Function Documentation

void MBA::init ( boost::shared_ptr< dVec U,
boost::shared_ptr< dVec V,
boost::shared_ptr< dVec Z 
) [inline]

Initialization that takes reference to arrays of scattered data. Note: see documenation of corresponding constructor above.

void MBA::init ( UCBspl::SplineSurface surf  ) 

Initialization of surface using data retrieved by getSplineSurface()

void MBA::setDomain ( double  umin,
double  vmin,
double  umax,
double  vmax 
) [inline]

Set the domain over which the surface is defined. The default is the xy-range of the scattered data.

Note:
This function can only be used to expand the domain boyond 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)

void MBA::setBaseType ( MBAbaseType  baseType  )  [inline]

Set surface base. Since the multilevel spline method is not affine-invariant, it is recommended to set a base surface over which the approximation is created.

Parameters:
MBA_ZERO,: A planar surface at level zero.
MBA_CONSTLS,: A konstant (least square) surface at the level of the mean z-value (default and recommended).
MBA_CONSTVAL,: A konstant surface at the level as given by setBaseValue.

void MBA::MBAalg ( int  m0,
int  n0,
int  h = 0,
int  smoothing_steps = 0 
)

Create a B-spline approximation to the scattered data.

Parameters:
m0,n0 (>=1): The initial size of the spline space in the hierarchical construction. If the rectangular domain is a square, m0=n0=1 is recommended. If the rectangular domain in the y-direction is twice of that the x-direction, m0=1, n0=2 is recommended. In general, if the rectangular domain in the y-direction is k times the length in the x-direction, m0=1, n0=k is recommended.
h,: Number of levels in the hierarchical construction. If, e.g., m0=n0=1 and h=8, The resulting spline surface has a coefficient grid of size 2^h + 3 = 259 in each direction if a cubic C^2 spline surface is created. (Currently a compile option UNIFORM_CUBIC_C1_SPLINES can be set in MBAunifBsplines.h to create a cubic C^1 surface. In this case the coefficient grid will be approx. twice the size of that in the C^2 case when using the same h value.)

UCBspl::SplineSurface MBA::getSplineSurface (  )  const [inline]

Retrieve the spline surface.

void MBA::getIndexDomain ( int &  m,
int &  n 
) const [inline]

Index domain of spline coefficient matrix. The surface can also be evaluated by f(i,j) and other functions with 0 <= i <= m and 0 <= j <= n. Note that the size of the tensor product grid (in the C2 case) is (m+3) x (n+3)

void MBA::cleanup ( int  type = 0  ) 

Clean-up array structures and reduce memory usage. Depending on the given type described below, temporary array structures will be deleted. All options preserve information necessary to evaluate the surface.

Parameters:
0,: Temporary (mathematical) stuff not used further after MBA::MBAalg.
1,: The scattered data arrays and data derived from them.
2,: Everything, including that above, that is not needed by evaluators.

boost::shared_ptr<GenMatrixType> MBA::PHI (  )  const [inline]

Get the coefficient grid of the tensor product spline surface.


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