BezierTriangle Class Reference

The Face type for BezierMesh. More...

#include </usr1/tp517/Tumble/trunk/src/tumble/cell.h>

Inheritance diagram for BezierTriangle:

[legend]
Collaboration diagram for BezierTriangle:
[legend]
List of all members.

Public Types

typedef FaceCell super

Public Member Functions

 BezierTriangle (PersistantStore &, BezierEdge *const es[3], const DataStore &, bool inverse_handed=false)
 Create a BezierTriangle.
void set_boundary_face (BoundaryFace *f)
virtual ~BezierTriangle ()
 Destructor.
int get_num_cp () const
const Point2Dget_cp (int i) const
 Get the geometry of the ith control point.
Point2Daccess_cp (int i)
 Read-write access to the geometry.
void set_cp (int i, const Point2D &p)
 Write access to the geometry.
ControlPoint get_control_point (int i) const
int get_num_dp () const
const LinearDataget_dp (int i) const
 Get the data at the ith control point.
LinearDataaccess_dp (int i)
 Read-write access to the data.
void set_dp (int i, const LinearData &p)
 Write access to the data.
DataPoint get_data_point (int i) const
int get_num_edges () const
BezierEdgeget_edge (int i)
 Since we always have exactly three, allow index-based retrieval of edges.
const BezierEdgeget_edge (int i) const
bool has_bdry_face () const
const BoundaryFaceget_bdry_face () const
BoundaryFaceget_bdry_face ()
BoundaryFaceget_bdry_face_or_null ()
double get_min_angle () const
int get_color () const
const BezierEdgeget_edge (int p0, int p1) const
 Get the BezierEdge that is a sub-cell of this Face, which shares the given control points.
BezierEdgeget_edge (int p0, int p1)
Point2D eval (double alpha, double beta) const
 Given parameters u, v, interpolate the geometric coordinates in the Bezier Basis.
LinearData dataeval (double alpha, double beta) const
 Given parameters u, v, interpolate the functional data in the Bezier Basis.
void evalIntermediate (double alpha, double beta, Point2D &p, Point2D &pc0, Point2D &pc1, Point2D &pc2) const
 Given a parameters u and v, interpolate the geometric coordinates in the Bezier Basis.
void dataevalIntermediate (double alpha, double beta, LinearData &d, LinearData &dc0, LinearData &dc1, LinearData &dc2) const
 Given a parameters u and v, interpolate the functional data in the Bezier Basis.
double jacobian (double alpha, double beta) const
 Gives the Jacobian determinate at the given barycentric coordinates.
void jacobian_matrix (double alpha, double beta, Point2D &pdu, Point2D &pdv) const
 Gives the full Jacobian at the given barycentric coordinates.
void bound_jacobian (double &min, double &max, double &area) const
 Approximates the jacobian determinate and the area of the triangle.
bool newton_find (const Point2D &p, double &u, double &v, int &niter, int &edge_num)
double area () const
 Gives the area of the BezierTriangle.
bool small_angle () const
 Determine if angle is small.
bool is_linear_inverted () const
int is_inverted () const
bool is_definitely_inverted () const
bool is_definitely_not_inverted () const
Point2D circumcenter () const
 Get the circumcenter of this triangle.
Point2D offcenter () const
 Get the offcenter of this triangle using the angle bound enforced by the BoundaryFace that contains this triangle. Error if no such face or bound exists.
Point2D offcenter (double beta) const
 Get the offcenter of this triangle (as defined by Ungor et al.).
void print () const
 Print info.
void print_mathematica () const
 Print info in a Mathematica friendly format.
 declare_iterators_canon (edge, edges, BezierEdge)

Private Member Functions

int get_edge_index_by_cps (int p0, int p1) const
 Given the index of two of the triangle's cps, find the edge that spans them.
void replace_cp (const ControlPoint &oldp, const ControlPoint &newp)
 BezierTriangle (const BezierTriangle &o)
BezierTriangleoperator= (const BezierTriangle &o)

Private Attributes

ControlPoint cp_ [6]
 An array of iterators, pointing to the geometric data (Point2D).
DataPoint dp_ [6]
 An array of iterators, pointing to the functional data (LinearData).
BoundaryFacebdry_face_
 A pointer to the BoundaryFace containing this BezierTriangle.

Friends

class BezierMesh
std::ostream & operator<< (std::ostream &stream, const BezierTriangle &t)

Detailed Description

The Face type for BezierMesh.

As we are using a second order BezierBasis for our geometric and functional basis, each BezierTriangle will have six ControlPoints and six DataPoints. Interpolation is done with the basis functions as:

$ \forall u,v \in [0,1]$ such that $ u+v <= 1 $

$ B(u, v) = (cp0)\cdot (1-u-v)^2 + (cp1)\cdot 2 v (1-u-v) + (cp2)\cdot v^2 + (cp3)\cdot 2 u (1-u-v) +(cp4)\cdot 2 u v + (cp5)\cdot u^2 $

btri.png

BezierTriangle

A BezierTriangle contains ControlPoints which point to the geometric Point2D information, and a DataPoints which point to the functional (LinearData) information. To get direct access to the information use the get_cp(int) and get_dp(int) methods.

Use get_edge(int) to get the edges of this face.

Also contains a pointer to the BoundaryFace containing this BezierTriangle.

Definition at line 719 of file cell.h.


Member Typedef Documentation

typedef FaceCell BezierTriangle::super

Definition at line 722 of file cell.h.


Constructor & Destructor Documentation

BezierTriangle::BezierTriangle ( PersistantStore ,
BezierEdge *const  es[3],
const DataStore ,
bool  inverse_handed = false 
)

Create a BezierTriangle.

Definition at line 1472 of file cell.C.

References cp_, dp_, BezierEdge::get_control_point(), get_control_point(), BezierEdge::get_vertex(), and map_idx().

Here is the call graph for this function:

virtual BezierTriangle::~BezierTriangle (  )  [inline, virtual]

Destructor.

Definition at line 731 of file cell.h.

BezierTriangle::BezierTriangle ( const BezierTriangle o  )  [private]


Member Function Documentation

void BezierTriangle::set_boundary_face ( BoundaryFace f  ) 

Definition at line 1558 of file cell.C.

References bdry_face_.

Referenced by BezierMesh::add_bezier_triangle().

int BezierTriangle::get_num_cp (  )  const [inline]

Definition at line 734 of file cell.h.

Referenced by access_cp(), access_dp(), get_control_point(), get_cp(), get_num_dp(), and replace_cp().

const Point2D & BezierTriangle::get_cp ( int  i  )  const

Get the geometry of the ith control point.

Definition at line 1566 of file cell.C.

References cp_, and get_num_cp().

Referenced by bound_jacobian(), circumcenter(), Visualization::draw_control_net(), BezierMesh::draw_insert_error(), eval(), evalIntermediate(), BezierMesh::find_start_triangle(), GhostTriangle::GhostTriangle(), is_definitely_not_inverted(), is_linear_inverted(), Visualization::is_visible(), jacobian_matrix(), offcenter(), operator<<(), and print_mathematica().

Here is the call graph for this function:

Point2D & BezierTriangle::access_cp ( int  i  ) 

Read-write access to the geometry.

Definition at line 1573 of file cell.C.

References PersistantList< Value >::iterator::access(), cp_, and get_num_cp().

Referenced by set_cp().

Here is the call graph for this function:

void BezierTriangle::set_cp ( int  i,
const Point2D p 
)

Write access to the geometry.

Definition at line 1580 of file cell.C.

References access_cp().

Here is the call graph for this function:

ControlPoint BezierTriangle::get_control_point ( int  i  )  const

Definition at line 1585 of file cell.C.

References cp_, and get_num_cp().

Referenced by BezierTriangle(), BezierMesh::consistency_check_helper(), get_edge_index_by_cps(), BezierMesh::locate_point_in_quadratic_triangle(), BezierMesh::orient_edge_triangle(), replace_cp(), and BezierMesh::tri_cps_off_edge().

Here is the call graph for this function:

int BezierTriangle::get_num_dp (  )  const [inline]

Definition at line 741 of file cell.h.

References get_num_cp().

Referenced by get_data_point(), and get_dp().

Here is the call graph for this function:

const LinearData & BezierTriangle::get_dp ( int  i  )  const

Get the data at the ith control point.

Definition at line 1592 of file cell.C.

References dp_, and get_num_dp().

Referenced by dataeval(), dataevalIntermediate(), GhostTriangle::GhostTriangle(), and operator<<().

Here is the call graph for this function:

LinearData & BezierTriangle::access_dp ( int  i  ) 

Read-write access to the data.

Definition at line 1599 of file cell.C.

References PersistantList< Value >::iterator::access(), dp_, and get_num_cp().

Referenced by set_dp().

Here is the call graph for this function:

void BezierTriangle::set_dp ( int  i,
const LinearData p 
)

Write access to the data.

Definition at line 1606 of file cell.C.

References access_dp().

Here is the call graph for this function:

DataPoint BezierTriangle::get_data_point ( int  i  )  const

Definition at line 1611 of file cell.C.

References dp_, and get_num_dp().

Referenced by BezierMesh::consistency_check_helper().

Here is the call graph for this function:

int BezierTriangle::get_num_edges (  )  const [inline]

Definition at line 749 of file cell.h.

Referenced by get_edge(), and get_edge_index_by_cps().

BezierEdge * BezierTriangle::get_edge ( int  i  ) 

Since we always have exactly three, allow index-based retrieval of edges.

Definition at line 1619 of file cell.C.

References FaceCell::begin_edges(), and get_num_edges().

Referenced by BezierMesh::check_consistency(), MeshOutput::ele_to_file(), get_edge(), get_edge_index_by_cps(), BezierMesh::locate_point_in_quadratic_triangle(), small_angle(), and BezierMesh::smooth_mesh().

Here is the call graph for this function:

const BezierEdge * BezierTriangle::get_edge ( int  i  )  const

Definition at line 1623 of file cell.C.

References FaceCell::begin_edges(), and get_num_edges().

Here is the call graph for this function:

bool BezierTriangle::has_bdry_face (  )  const

Definition at line 1628 of file cell.C.

References bdry_face_.

Referenced by get_bdry_face(), and operator<<().

const BoundaryFace * BezierTriangle::get_bdry_face (  )  const

Definition at line 1631 of file cell.C.

References bdry_face_, and has_bdry_face().

Referenced by get_color(), get_min_angle(), and operator<<().

Here is the call graph for this function:

BoundaryFace * BezierTriangle::get_bdry_face (  ) 

Definition at line 1635 of file cell.C.

References bdry_face_, and has_bdry_face().

Here is the call graph for this function:

BoundaryFace * BezierTriangle::get_bdry_face_or_null (  ) 

Definition at line 1639 of file cell.C.

References bdry_face_.

Referenced by MeshOutput::ele_to_file(), BezierMesh::flip(), BezierMesh::insert_edge_midpoint(), and BezierMesh::insert_point().

double BezierTriangle::get_min_angle (  )  const

Definition at line 1643 of file cell.C.

References get_bdry_face(), and BoundaryFace::get_min_angle().

Referenced by offcenter(), and small_angle().

Here is the call graph for this function:

int BezierTriangle::get_color (  )  const

Definition at line 1646 of file cell.C.

References get_bdry_face(), and BoundaryFace::get_color().

Referenced by Visualization::draw_bezier_triangle(), and EPSWrite::shade_triangle().

Here is the call graph for this function:

const BezierEdge * BezierTriangle::get_edge ( int  p0,
int  p1 
) const

Get the BezierEdge that is a sub-cell of this Face, which shares the given control points.

Parameters:
p0 The index in cps of the edge's first control point
p1 The index in cps of the edge's second control point
Returns:
The edge with the given control points

Definition at line 1655 of file cell.C.

References get_edge(), and get_edge_index_by_cps().

Here is the call graph for this function:

BezierEdge * BezierTriangle::get_edge ( int  p0,
int  p1 
)

Definition at line 1659 of file cell.C.

References get_edge(), and get_edge_index_by_cps().

Here is the call graph for this function:

Point2D BezierTriangle::eval ( double  u,
double  v 
) const [virtual]

Given parameters u, v, interpolate the geometric coordinates in the Bezier Basis.

Assumes that $ u,v,u+v \in [0,1] $.

returns:$ B(u, v) = (cp0)\cdot (1-u-v)^2 + (cp1)\cdot 2 v (1-u-v) + (cp2)\cdot v^2 + (cp3)\cdot 2 u (1-u-v) +(cp4)\cdot 2 u v + (cp5)\cdot u^2 $

Parameters:
u The first parameter (u in [0,1] )
v The second parameter (v in [0,1] )
Returns:
The interpolated coordinates

Implements CurvedTriangle.

Definition at line 1675 of file cell.C.

References get_cp().

Referenced by BezierMesh::clean_insert_point(), Visualization::draw_bezier_triangle(), Visualization::draw_point(), and BezierMesh::interpolant_error().

Here is the call graph for this function:

LinearData BezierTriangle::dataeval ( double  u,
double  v 
) const

Given parameters u, v, interpolate the functional data in the Bezier Basis.

Assumes that $ u,v,u+v \in [0,1] $.

returns:$ B(u, v) = (dp0)\cdot (1-u-v)^2 + (dp1)\cdot 2 v (1-u-v) + (dp2)\cdot v^2 + (dp3)\cdot 2 u (1-u-v) +(dp4)\cdot 2 u v + (dp5)\cdot u^2 $

Parameters:
u The first parameter (u in [0,1] )
v The second parameter (v in [0,1] )
Returns:
The interpolated functional value

Definition at line 1699 of file cell.C.

References get_dp().

Referenced by Visualization::draw_bezier_triangle(), and BezierMesh::interpolant_error().

Here is the call graph for this function:

void BezierTriangle::evalIntermediate ( double  u,
double  v,
Point2D p,
Point2D pc0,
Point2D pc1,
Point2D pc2 
) const

Given a parameters u and v, interpolate the geometric coordinates in the Bezier Basis.

Runs the de Casteljau algorithm, and returns intermediates as well as the final interpolated value.

Parameters:
u The first parameter (u in [0, 1])
v The second parameter (v in [0, 1])
[out] p The interpolated geometric point
[out] pc0 The first intermediate
[out] pc1 The second intermediate
[out] pc2 The third intermediate

Definition at line 1723 of file cell.C.

References get_cp().

Referenced by BezierMesh::insert_edge_midpoint(), and BezierMesh::insert_point().

Here is the call graph for this function:

void BezierTriangle::dataevalIntermediate ( double  u,
double  v,
LinearData d,
LinearData dc0,
LinearData dc1,
LinearData dc2 
) const

Given a parameters u and v, interpolate the functional data in the Bezier Basis.

Runs the de Casteljau algorithm, and returns intermediates as well as the final interpolated value.

Parameters:
u The first parameter (u in [0, 1])
v The second parameter (v in [0, 1])
[out] d The interpolated functional data
[out] dc0 The first intermediate
[out] dc1 The second intermediate
[out] dc2 The third intermediate

Definition at line 1755 of file cell.C.

References get_dp().

Referenced by BezierMesh::insert_edge_midpoint(), and BezierMesh::insert_point().

Here is the call graph for this function:

double BezierTriangle::jacobian ( double  u,
double  v 
) const

Gives the Jacobian determinate at the given barycentric coordinates.

Parameters:
u The first parameter
v The second parameter
Returns:
The Jacobian determinate

Definition at line 1781 of file cell.C.

References Point2D::cross(), and jacobian_matrix().

Referenced by is_definitely_inverted(), and BezierMesh::reinterpolate().

Here is the call graph for this function:

void BezierTriangle::jacobian_matrix ( double  u,
double  v,
Point2D pdu,
Point2D pdv 
) const [virtual]

Gives the full Jacobian at the given barycentric coordinates.

This returns the matrix as a pair of column vectors, where each column vector represents the partial with repsect to one of the two parameters.

Parameters:
u The first parameter
v The second parameter
[out] pdu A vector with coordinates representing the partial derivative with respect to u
[out] pdv A vector with coordinates representing the partial derivative with respect to v

Implements CurvedTriangle.

Definition at line 1798 of file cell.C.

References get_cp().

Referenced by jacobian().

Here is the call graph for this function:

void BezierTriangle::bound_jacobian ( double &  min,
double &  max,
double &  area 
) const

Approximates the jacobian determinate and the area of the triangle.

Returns a range [min, max ] for the Jacobian Determinate over the whole triangle. Thus, at no pair (u , v ) will the Jacobian determinate be less than min or more than max.

Also calculates the exact area as a side effect.

Parameters:
min The lower bound for the Jacobian determinate
max The upper bound for the Jacobian determinate
area The area of the BezierTriangle

Definition at line 1825 of file cell.C.

References get_cp().

Referenced by area(), and BezierMesh::smooth_mesh().

Here is the call graph for this function:

bool BezierTriangle::newton_find ( const Point2D p,
double &  u,
double &  v,
int &  niter,
int &  edge_num 
)

Definition at line 1855 of file cell.C.

References get_edge_index_by_cps(), and CurvedTriangle::newton_find().

Referenced by BezierMesh::locate_point_in_quadratic_triangle().

Here is the call graph for this function:

double BezierTriangle::area (  )  const

Gives the area of the BezierTriangle.

Definition at line 1873 of file cell.C.

References bound_jacobian().

Referenced by BezierMesh::remove_large_triangles().

Here is the call graph for this function:

bool BezierTriangle::small_angle (  )  const

Determine if angle is small.

This function looks up the angle bound from the Triangle's face. It then examines each of the three angles in the triangle, treating the edges as linear. If two edges are both on a boundary, then we ignore that angle, so as not to recurse on small angles.

Returns:
true if the angle is small

Definition at line 1888 of file cell.C.

References Point2D::angle(), BezierEdge::get_bdry_edge(), BezierVertex::get_cp(), get_edge(), get_min_angle(), BezierEdge::get_vertex(), BezierEdge::is_boundary(), and PI.

Referenced by BezierMesh::remove_small_angles().

Here is the call graph for this function:

bool BezierTriangle::is_linear_inverted (  )  const

Definition at line 1940 of file cell.C.

References get_cp(), and Point2D::is_left_of().

Here is the call graph for this function:

int BezierTriangle::is_inverted (  )  const

Definition at line 1949 of file cell.C.

References is_definitely_inverted(), and is_definitely_not_inverted().

Referenced by MeshOutput::ele_to_file(), BezierMesh::insert_point(), and BezierMesh::print_check_inverted().

Here is the call graph for this function:

bool BezierTriangle::is_definitely_inverted (  )  const

Definition at line 1971 of file cell.C.

References ALPHA, BETA, and jacobian().

Referenced by is_inverted().

Here is the call graph for this function:

bool BezierTriangle::is_definitely_not_inverted (  )  const

Definition at line 1984 of file cell.C.

References get_cp().

Referenced by is_inverted().

Here is the call graph for this function:

Point2D BezierTriangle::circumcenter (  )  const

Get the circumcenter of this triangle.

Returns:
The circumcenter of the BezierTriangle

Definition at line 2019 of file cell.C.

References get_cp().

Referenced by BezierMesh::clean_insert_circumcenter(), and offcenter().

Here is the call graph for this function:

Point2D BezierTriangle::offcenter (  )  const

Get the offcenter of this triangle using the angle bound enforced by the BoundaryFace that contains this triangle. Error if no such face or bound exists.

Returns:
the offcenter of the BezierTriangle

Definition at line 2039 of file cell.C.

References get_min_angle(), and PI.

Referenced by offcenter().

Here is the call graph for this function:

Point2D BezierTriangle::offcenter ( double  beta  )  const

Get the offcenter of this triangle (as defined by Ungor et al.).

Parameters:
beta : The radius-edge ratio for the offcenter
Returns:
the offcenter of the BezierTriangle

Definition at line 2052 of file cell.C.

References circumcenter(), get_cp(), Point2D::mag(), and offcenter().

Here is the call graph for this function:

void BezierTriangle::print (  )  const

Print info.

Definition at line 2093 of file cell.C.

void BezierTriangle::print_mathematica (  )  const

Print info in a Mathematica friendly format.

Definition at line 2099 of file cell.C.

References get_cp(), Point2D::x(), and Point2D::y().

Here is the call graph for this function:

BezierTriangle::declare_iterators_canon ( edge  ,
edges  ,
BezierEdge   
)

int BezierTriangle::get_edge_index_by_cps ( int  p0,
int  p1 
) const [private]

Given the index of two of the triangle's cps, find the edge that spans them.

p0 and p1 must correspond to the vertices of a single edge.

As the numbering of control points has been done in a counter-clockwise manner, this may not coorsepond to the ordering of edges. This is the case if the triangle is in the inverted set of the BezierMesh.

In all triangles, get_edge(0) spans control points 0-2; get_edge(1) however may span either 2-5 or 5-0 (get_edge(2) spans the other).

TODO: this should be abolished and the edge numbers made to agree with the control point numbers.

Definition at line 2165 of file cell.C.

References BezierVertex::get_control_point(), get_control_point(), get_edge(), get_num_edges(), and BezierEdge::get_vertex().

Referenced by get_edge(), and newton_find().

Here is the call graph for this function:

void BezierTriangle::replace_cp ( const ControlPoint oldp,
const ControlPoint newp 
) [private]

Definition at line 2187 of file cell.C.

References cp_, get_control_point(), and get_num_cp().

Here is the call graph for this function:

BezierTriangle& BezierTriangle::operator= ( const BezierTriangle o  )  [private]


Friends And Related Function Documentation

friend class BezierMesh [friend]

Definition at line 723 of file cell.h.

std::ostream& operator<< ( std::ostream &  stream,
const BezierTriangle t 
) [friend]


Member Data Documentation

ControlPoint BezierTriangle::cp_[6] [private]

An array of iterators, pointing to the geometric data (Point2D).

Definition at line 797 of file cell.h.

Referenced by access_cp(), BezierTriangle(), get_control_point(), get_cp(), and replace_cp().

DataPoint BezierTriangle::dp_[6] [private]

An array of iterators, pointing to the functional data (LinearData).

Definition at line 798 of file cell.h.

Referenced by access_dp(), BezierTriangle(), get_data_point(), and get_dp().

BoundaryFace* BezierTriangle::bdry_face_ [private]

A pointer to the BoundaryFace containing this BezierTriangle.

Definition at line 799 of file cell.h.

Referenced by get_bdry_face(), get_bdry_face_or_null(), has_bdry_face(), and set_boundary_face().


The documentation for this class was generated from the following files:
Generated on Mon May 24 09:53:33 2010 for TUMBLE by  doxygen 1.5.2