/*=====================================================================
File: nurbs.cpp
Purpose:
Revision: $Id: nurbs.cpp,v 1.3 2002/05/24 17:27:24 philosophil Exp $
Author: Philippe Lavoie (3 Oct, 1996)
Modified by:
Copyright notice:
Copyright (C) 1996-1997 Philippe Lavoie
This library is free software; you can redistribute it and/or
modify it under the terms of the GNU Library General Public
License as published by the Free Software Foundation; either
version 2 of the License, or (at your option) any later version.
This library is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
Library General Public License for more details.
You should have received a copy of the GNU Library General Public
License along with this library; if not, write to the Free
Software Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
=====================================================================*/
#include <nurbs.h>
#include <fstream>
#include <string.h>
#include <nurbsS.h>
#include "integrate.h"
#include <malloc.h>
/*!
*/
namespace PLib {
/*!
\brief default constructor
\author Philippe Lavoie
\date 24 January 1997
*/
template <class T, int N>
NurbsCurve<T,N>::NurbsCurve(): P(1),U(1),deg_(0)
{
}
/*!
\brief A copy constructor.
\param nurb the NURBS curve to copy
\author Philippe Lavoie
\date 24 January 1997
*/
template <class T, int N>
NurbsCurve<T,N>::NurbsCurve(const NurbsCurve<T,N>& nurb):
ParaCurve<T,N>(), P(nurb.P),U(nurb.U),deg_(nurb.deg_)
{
}
/*!
\brief Resets a NURBS curve to new values
\param P1 the new values for the control points
\param U1 the new values for the knot vector
\param Degree the new degree of the curve
\warning The size of P1,U1 and Degree must agree: P.n()+degree+1=U.n()
\author Philippe Lavoie
\date 24 January 1997
*/
template <class T, int N>
void NurbsCurve<T,N>::reset(const Vector< HPoint_nD<T,N> >& P1, const Vector<T> &U1, int Degree) {
int nSize = P1.n() ;
int mSize = U1.n() ;
deg_ = Degree ;
if(nSize != mSize-deg_-1){
#ifdef USE_EXCEPTION
throw NurbsSizeError(P1.n(),U1.n(),Degree) ;
#else
Error err("reset");
err << "Invalid input size for the control points and the knot vector when reseting a Nurbs Curve.\n";
err << nSize << " control points and " << mSize << " knots\n" ;
err.fatal() ;
#endif
}
P.resize(P1.n()) ;
U.resize(U1.n()) ;
P = P1 ;
U = U1 ;
}
/*!
\brief Constructor with control points in 4D
\param P1 the control points
\param U1 the knot vector
\param Degree the degree of the curve
\warning The size of P1,U1 and Degree must agree: P.n()+degree+1=U.n()
\author Philippe Lavoie
\date 24 January 1997
*/
template <class T, int N>
NurbsCurve<T,N>::NurbsCurve(const Vector< HPoint_nD<T,N> >& P1, const Vector<T> &U1, int Degree): P(P1), U(U1), deg_(Degree)
{
if(P.n() != U.n()-deg_-1){
#ifdef USE_EXCEPTION
throw NurbsSizeError(P.n(),U.n(),deg_) ;
#else
Error err("NurbsCurve(P1,U1,Degree)");
err << "Invalid input size for the control points and the knot vector.\n";
err << P.n() << " control points and " << U.n() << " knots\n" ;
err.fatal() ;
#endif
}
}
/*!
\brief Constructor with control points in 3D
\param P1 --> the control point vector
\param W --> the weight for each control points
\param U1 --> the knot vector
\param Degree --> the degree of the curve
\warning The size of P1,U1 and Degree must agree: P.n()+degree+1=U.n()
\author Philippe Lavoie
\date 24 January 1997
*/
template <class T, int N>
NurbsCurve<T,N>::NurbsCurve(const Vector< Point_nD<T,N> >& P1, const Vector<T>& W, const Vector<T>& U1, int Degree): P(P1.n()), U(U1), deg_(Degree)
{
int nSize = P1.n() ;
int mSize = U1.n() ;
if(nSize != mSize-deg_-1){
#ifdef USE_EXCEPTION
throw NurbsSizeError(P.n(),U.n(),deg_) ;
#else
Error err("NurbsCurve(P1,W,U1,Degree)") ;
err << "Invalid input size for the control points and the knot vector.\n" ;
err << nSize << " control points and " << mSize << " knots\n" ;
err.fatal() ;
#endif
}
if(nSize != W.n()){
#ifdef USE_EXCEPTION
throw NurbsInputError(nSize,W.n()) ;
#else
Error err("NurbsCurve(P1,W,U1,Degree)") ;
err << "Size mismatched between the control points and the weights\n" ;
err << "ControlPoints size = " << nSize << ", Weight size = " << W.n() << endl ;
err.fatal() ;
#endif
}
for(int i = 0 ;i<nSize;i++){
const Point_nD<T,N>& pt = P1[i] ; // This makes the SGI compiler happy
for(int j=0;j<N;j++)
P[i].data[j] = pt.data[j] * W[i] ;
P[i].w() = W[i] ;
}
}
/*!
\brief The assignment operator for a NURBS curve
\param curve the NURBS curve to copy
\return A reference to itself
\warning The curve being copied must be valid, otherwise strange
results might occur.
\author Philippe Lavoie
\date 24 January 1997
*/
template <class T, int N>
NurbsCurve<T,N>& NurbsCurve<T,N>::operator=(const NurbsCurve<T,N>& curve) {
if(curve.U.n() != curve.P.n()+curve.deg_+1){
#ifdef USE_EXCEPTION
throw NurbsSizeError(curve.P.n(),curve.U.n(),curve.deg_) ;
#else
Error err("operator=") ;
err << "Invalid assignment... the curve being assigned to isn't valid\n" ;
err.fatal() ;
#endif
}
deg_ = curve.deg_ ;
U = curve.U ;
P = curve.P ;
if(U.n()!=P.n()+deg_+1){
#ifdef USE_EXCEPTION
throw NurbsSizeError(P.n(),U.n(),deg_) ;
#else
Error err("operator=") ;
err << "Error in assignment... couldn't assign properly the vectors\n" ;
err.fatal() ;
#endif
}
return *this ;
}
/*!
\brief draws a NURBS curve on an image
This will draw very primitively the NURBS curve on an image.
The drawing assumes the line is only in the xy plane (the z
is not used for now).
The algorithm finds the points on the curve at a \a step
parametric intervall between them and join them by a line.
No fancy stuff.
\param Img <-- draws the nurbs curve to this Image
\param color --> the line is drawn in this color
\param step --> the parametric distance between two computed points.
\author Philippe Lavoie
\date 24 January 1997
*/
template <class T, int N>
void NurbsCurve<T,N>::drawImg(Image_UBYTE& Img,unsigned char color,T step){
Point_nD<T,N> a1,a2 ;
T u_max = U[U.n()-1-deg_] ;
if(step<=0)
step = 0.01 ;
a1 = this->pointAt(U[deg_]) ;
T u ;
int i1,j1,i2,j2 ;
getCoordinates(a1,i1,j1,Img.rows(),Img.cols()) ;
for(u=U[deg_]+step ; u < u_max+(step/2.0) ; u+=step){ // the <= u_max doesn't work
a2 = this->pointAt(u) ;
if(!getCoordinates(a2,i2,j2,Img.rows(),Img.cols()))
continue ;
Img.drawLine(i1,j1,i2,j2,color) ;
i1 = i2 ;
j1 = j2 ;
}
a2 = this->pointAt(U[P.n()]) ;
if(getCoordinates(a2,i2,j2,Img.rows(),Img.cols()))
Img.drawLine(i1,j1,i2,j2,color) ;
}
/*!
\brief Draws a NURBS curve on an image
This will draw very primitively the NURBS curve on an image.
The drawing assumes the line is only in the xy plane (the z
is not used for now).
The algorithm finds the points on the curve at a \a step
parametric intervall between them and join them by a line.
No fancy stuff.
\param Img draws the nurbs curve to this Image
\param color the line is drawn in this color
\param step the parametric distance between two computed points.
\author Philippe Lavoie
\date 24 January 1997
*/
template <class T, int N>
void NurbsCurve<T,N>::drawImg(Image_Color& Img,const Color& color,T step){
Point_nD<T,N> a1,a2 ;
T u_max = U[U.n()-1-deg_] ;
if(step<=0)
step = 0.01 ;
a1 = this->pointAt(U[deg_]) ;
int i1,j1,i2,j2 ;
getCoordinates(a1,i1,j1,Img.rows(),Img.cols()) ;
T u ;
for(u=U[deg_]+step ; u < u_max+(step/2.0
没有合适的资源?快使用搜索试试~ 我知道了~
nurbs++-3.0.11 nurbs类库
共254个文件
cpp:131个
h:40个
dsp:21个
3星 · 超过75%的资源 需积分: 14 38 下载量 54 浏览量
2008-09-17
19:33:54
上传
评论
收藏 675KB ZIP 举报
温馨提示
nurbs++-3.0.11 nurbs类库包含了nurbs曲线的生成的代码,是图形学的重要工具
资源推荐
资源详情
资源评论
收起资源包目录
nurbs++-3.0.11 nurbs类库 (254个子文件)
nurbs++-config.1 2KB
Makefile.am 2KB
Makefile.am 2KB
Makefile.am 1KB
Makefile.am 590B
Makefile.am 572B
Makefile.am 501B
Makefile.am 426B
Makefile.am 394B
Makefile.am 356B
Makefile.am 353B
Makefile.am 350B
Makefile.am 108B
Makefile.am 108B
Makefile.am 87B
Makefile.am 23B
AUTHORS 258B
config_mvc.bat 79B
ChangeLog 3KB
configure 349KB
COPYING 0B
nurbs.cpp 180KB
nurbsS.cpp 127KB
nurbsGL.cpp 90KB
nurbsSub.cpp 50KB
hnurbsS.cpp 34KB
surface.cpp 34KB
matrixMat.cpp 31KB
matrix.cpp 21KB
intccq.cpp 18KB
fft.cpp 18KB
nurbsArray.cpp 17KB
rec_filter.cpp 17KB
f_nurbsSub.cpp 15KB
hnurbsS_sp.cpp 14KB
vector.cpp 14KB
matrixRT.cpp 13KB
d_nurbsSub.cpp 12KB
tri_spline.cpp 12KB
curve.cpp 11KB
hnurbs.cpp 11KB
barray2d.cpp 10KB
barray.cpp 10KB
image.cpp 9KB
nurbsS_sp.cpp 8KB
chebexp.cpp 7KB
barray2d_hpoint.cpp 7KB
statistic.cpp 7KB
matrix_hpoint.cpp 7KB
test_barray2d.cpp 6KB
matrix_point.cpp 6KB
image_.cpp 6KB
color.cpp 5KB
f_nurbs.cpp 5KB
d_nurbs.cpp 5KB
error.cpp 5KB
barray_hpoint.cpp 4KB
topengl.cpp 4KB
test_matrix.cpp 4KB
vector_hpoint.cpp 4KB
vector_point.cpp 4KB
tnApprox.cpp 4KB
nurbs_sp.cpp 3KB
matrix_complex.cpp 3KB
filter.cpp 3KB
test_nurbs.cpp 3KB
thnurbsS.cpp 3KB
intccq_.cpp 3KB
matrix_uchar.cpp 3KB
tnlength.cpp 3KB
vector_complex.cpp 3KB
matrix_int.cpp 3KB
barray_point.cpp 3KB
matrix_char.cpp 3KB
tnMovePoint.cpp 3KB
barray2d_point.cpp 3KB
tmatrixMat.cpp 3KB
statistic_.cpp 3KB
trispline.cpp 3KB
tclose.cpp 2KB
tnurbsS.cpp 2KB
list.cpp 2KB
tnsMovePoint.cpp 2KB
tnInterp.cpp 2KB
filter_.cpp 2KB
barray_complex.cpp 2KB
barray2d_uchar.cpp 2KB
matrix_double.cpp 2KB
matrix_float.cpp 2KB
tnurbs.cpp 2KB
vector_double.cpp 2KB
vector_float.cpp 2KB
vector_uchar.cpp 2KB
tnurbsSub.cpp 2KB
vector_int.cpp 2KB
tnsSweep.cpp 2KB
vector_char.cpp 2KB
tnurbs_sp.cpp 2KB
barray_uchar.cpp 2KB
barray_coordinate.cpp 2KB
共 254 条
- 1
- 2
- 3
资源评论
- jinhaiqiang2012-08-02这个,能附上说明吗,不好理解啊。
p4bailinp4
- 粉丝: 0
- 资源: 2
上传资源 快速赚钱
- 我的内容管理 展开
- 我的资源 快来上传第一个资源
- 我的收益 登录查看自己的收益
- 我的积分 登录查看自己的积分
- 我的C币 登录后查看C币余额
- 我的收藏
- 我的下载
- 下载帮助
最新资源
- Screenshot_20240427_031602.jpg
- 网页PDF_2024年04月26日 23-46-14_QQ浏览器网页保存_QQ浏览器转格式(6).docx
- 直接插入排序,冒泡排序,直接选择排序.zip
- 在排序2的基础上,再次对快排进行优化,其次增加快排非递归,归并排序,归并排序非递归版.zip
- 实现了7种排序算法.三种复杂度排序.三种nlogn复杂度排序(堆排序,归并排序,快速排序)一种线性复杂度的排序.zip
- 冒泡排序 直接选择排序 直接插入排序 随机快速排序 归并排序 堆排序.zip
- 课设-内部排序算法比较 包括冒泡排序、直接插入排序、简单选择排序、快速排序、希尔排序、归并排序和堆排序.zip
- Python排序算法.zip
- C语言实现直接插入排序、希尔排序、选择排序、冒泡排序、堆排序、快速排序、归并排序、计数排序,并带图详解.zip
- 常用工具集参考用于图像等数据处理
资源上传下载、课程学习等过程中有任何疑问或建议,欢迎提出宝贵意见哦~我们会及时处理!
点击此处反馈
安全验证
文档复制为VIP权益,开通VIP直接复制
信息提交成功