#include <nurbs.hpp>
|
| KnotVector () |
| Create KnotVector. More...
|
|
| KnotVector (std::istream &input) |
|
| KnotVector (int Order_, int NCP) |
|
| KnotVector (const KnotVector &kv) |
|
KnotVector & | operator= (const KnotVector &kv) |
|
int | GetNE () const |
|
int | GetNKS () const |
|
int | GetNCP () const |
|
int | GetOrder () const |
|
int | Size () const |
|
void | GetElements () |
| Count the number of elements. More...
|
|
bool | isElement (int i) const |
|
double | getKnotLocation (double xi, int ni) const |
|
int | findKnotSpan (double u) const |
|
void | CalcShape (Vector &shape, int i, double xi) const |
|
void | CalcDShape (Vector &grad, int i, double xi) const |
|
void | CalcDnShape (Vector &gradn, int n, int i, double xi) const |
|
void | CalcD2Shape (Vector &grad2, int i, double xi) const |
|
void | FindMaxima (Array< int > &ks, Vector &xi, Vector &u) |
|
void | FindInterpolant (Array< Vector *> &x) |
|
void | Difference (const KnotVector &kv, Vector &diff) const |
|
void | UniformRefinement (Vector &newknots) const |
|
KnotVector * | DegreeElevate (int t) const |
|
void | Flip () |
|
void | Print (std::ostream &out) const |
|
void | PrintFunctions (std::ostream &out, int samples=11) const |
|
| ~KnotVector () |
| Destroys KnotVector. More...
|
|
double & | operator[] (int i) |
|
const double & | operator[] (int i) const |
|
Definition at line 32 of file nurbs.hpp.
◆ KnotVector() [1/4]
mfem::KnotVector::KnotVector |
( |
| ) |
|
|
inline |
◆ KnotVector() [2/4]
mfem::KnotVector::KnotVector |
( |
std::istream & |
input | ) |
|
◆ KnotVector() [3/4]
mfem::KnotVector::KnotVector |
( |
int |
Order_, |
|
|
int |
NCP |
|
) |
| |
◆ KnotVector() [4/4]
mfem::KnotVector::KnotVector |
( |
const KnotVector & |
kv | ) |
|
|
inline |
◆ ~KnotVector()
mfem::KnotVector::~KnotVector |
( |
| ) |
|
|
inline |
◆ CalcD2Shape()
void mfem::KnotVector::CalcD2Shape |
( |
Vector & |
grad2, |
|
|
int |
i, |
|
|
double |
xi |
|
) |
| const |
|
inline |
◆ CalcDnShape()
void mfem::KnotVector::CalcDnShape |
( |
Vector & |
gradn, |
|
|
int |
n, |
|
|
int |
i, |
|
|
double |
xi |
|
) |
| const |
◆ CalcDShape()
void mfem::KnotVector::CalcDShape |
( |
Vector & |
grad, |
|
|
int |
i, |
|
|
double |
xi |
|
) |
| const |
◆ CalcShape()
void mfem::KnotVector::CalcShape |
( |
Vector & |
shape, |
|
|
int |
i, |
|
|
double |
xi |
|
) |
| const |
◆ DegreeElevate()
KnotVector * mfem::KnotVector::DegreeElevate |
( |
int |
t | ) |
const |
Return a new KnotVector with elevated degree by repeating the endpoints of the knot vector.
- Note
- The returned object should be deleted by the caller.
Definition at line 58 of file nurbs.cpp.
◆ Difference()
void mfem::KnotVector::Difference |
( |
const KnotVector & |
kv, |
|
|
Vector & |
diff |
|
) |
| const |
◆ FindInterpolant()
void mfem::KnotVector::FindInterpolant |
( |
Array< Vector *> & |
x | ) |
|
Global curve interpolation through the points x. x is an array with the length of the spatial dimension containing vectors with spatial coordinates. The controlpoints of the interpolated curve are given in x in the same form.
Definition at line 409 of file nurbs.cpp.
◆ findKnotSpan()
int mfem::KnotVector::findKnotSpan |
( |
double |
u | ) |
const |
◆ FindMaxima()
Gives the locations of the maxima of the knotvector in reference space. The function gives the knotspan ks, the coordinate in the knotspan xi and the coordinate of the maximum in parameter space u
Definition at line 346 of file nurbs.cpp.
◆ Flip()
void mfem::KnotVector::Flip |
( |
| ) |
|
◆ GetElements()
void mfem::KnotVector::GetElements |
( |
| ) |
|
Count the number of elements.
Definition at line 101 of file nurbs.cpp.
◆ getKnotLocation()
double mfem::KnotVector::getKnotLocation |
( |
double |
xi, |
|
|
int |
ni |
|
) |
| const |
|
inline |
◆ GetNCP()
int mfem::KnotVector::GetNCP |
( |
| ) |
const |
|
inline |
◆ GetNE()
int mfem::KnotVector::GetNE |
( |
| ) |
const |
|
inline |
◆ GetNKS()
int mfem::KnotVector::GetNKS |
( |
| ) |
const |
|
inline |
◆ GetOrder()
int mfem::KnotVector::GetOrder |
( |
| ) |
const |
|
inline |
◆ isElement()
bool mfem::KnotVector::isElement |
( |
int |
i | ) |
const |
|
inline |
◆ operator=()
◆ operator[]() [1/2]
double& mfem::KnotVector::operator[] |
( |
int |
i | ) |
|
|
inline |
◆ operator[]() [2/2]
const double& mfem::KnotVector::operator[] |
( |
int |
i | ) |
const |
|
inline |
◆ Print()
void mfem::KnotVector::Print |
( |
std::ostream & |
out | ) |
const |
◆ PrintFunctions()
void mfem::KnotVector::PrintFunctions |
( |
std::ostream & |
out, |
|
|
int |
samples = 11 |
|
) |
| const |
◆ Size()
int mfem::KnotVector::Size |
( |
| ) |
const |
|
inline |
◆ UniformRefinement()
void mfem::KnotVector::UniformRefinement |
( |
Vector & |
newknots | ) |
const |
◆ knot
◆ MaxOrder
const int mfem::KnotVector::MaxOrder = 10 |
|
staticprotected |
◆ NumOfControlPoints
int mfem::KnotVector::NumOfControlPoints |
|
protected |
◆ NumOfElements
int mfem::KnotVector::NumOfElements |
|
protected |
◆ Order
int mfem::KnotVector::Order |
|
protected |
The documentation for this class was generated from the following files: