Classes | Typedefs | Functions | Variables

crv Namespace Reference

the curving functions are contained in this namespace More...

Classes

class  MeshCurver
 Base Mesh curving object. More...
class  InterpolatingCurver
 curves an already changed mesh More...
class  BezierCurver
 this curves a mesh with Bezier shapes More...
class  GregoryCurver
 this curves a mesh with 4th order G1 Patches More...
class  Quality
 class to store matrices used in quality assessment and validity checking More...
class  Adapt
 base crv::Adapt class, looks the same as ma::Adapt, but carries tag identifying validity (see crvShape.h) More...

Typedefs

typedef void(* bezierShape )(int P, apf::Vector3 const &xi, apf::NewArray< double > &values)
 typedef for table of shape functions
typedef void(* bezierShapeGrads )(int P, apf::Vector3 const &xi, apf::NewArray< apf::Vector3 > &grads)
 typedef for table of shape function gradients
typedef void(* SubdivisionFunction )(int P, apf::NewArray< double > &nodes, apf::NewArray< double > *subNodes)
 typedef for table of jacobian det subdivision functions

Functions

void setOrder (const int order)
 sets order used in bezier shape functions
int getOrder ()
 gets order used in bezier shape functions
void setBlendingOrder (const int type, const int b)
 sets the blending order, if shape blending is used
int getBlendingOrder (const int type)
 gets the blending order
int countNumberInvalidElements (apf::Mesh2 *m)
 count invalid elements of the mesh
ma::InputconfigureShapeCorrection (ma::Mesh *m, ma::SizeField *f=0, ma::SolutionTransfer *s=0)
 configure for fixing invalid elements
void adapt (ma::Input *in)
 crv adapt with custom configuration
void stats (ma::Input *in, std::vector< double > &edgeLengths, std::vector< double > &linearQualities, std::vector< double > &curvedQualities, bool inMetric=true)
 crv stats to get statistic information about the mesh
apf::FieldShapegetBezier (int order)
 Get the Bezier Curve or Shape of some order.
apf::FieldShapegetGregory ()
 Get the 4th order Gregory Surface.
double getQuality (apf::Mesh *m, apf::MeshEntity *e)
 computes min det Jacobian / max det Jacobian. Quality::getQuality should be used if multiple elements checked in a row
int checkValidity (apf::Mesh *m, apf::MeshEntity *e, int algorithm=2)
 checks validity of it and returns integer corresponding to invalid entity. Quality::checkValidity should be used if multiple elements checked in a row
QualitymakeQuality (apf::Mesh *m, int algorithm=2)
 use this to make a quality object with the correct dimension
double interpolationError (apf::Mesh *m, apf::MeshEntity *e, int n)
 computes interpolation error of a curved entity on a mesh
void writeCurvedVtuFiles (apf::Mesh *m, int type, int n, const char *prefix)
 Visualization, writes file for specified type, n is number of subdivisions, higher number -> better resolution, but bigger file.
void writeCurvedWireFrame (apf::Mesh *m, int n, const char *prefix)
 Visualization, writes wireframe of the curved mesh, n is number of subdivisions, higher number -> better resolution, but bigger file.
void writeControlPointVtuFiles (apf::Mesh *m, const char *prefix)
 Visualization, writes file of control nodes for each entity.
void writeInterpolationPointVtuFiles (apf::Mesh *m, const char *prefix)
 Visualization, writes file of shapes evaluated at node xi for each entity.
int getTriNodeIndex (int P, int i, int j)
 publically accessible functions
void fail (const char *why) __attribute__((noreturn))
 crv fail function
void changeMeshOrder (apf::Mesh2 *m, int newOrder)
 change the order of a Bezier Mesh
double Bij (const int i, const int j, const double u, const double v)
 polynomial part of bernstein polynomial, Bij, Bijk, Bijkl
double Bij (const int ij[], const double xi[])
 a different form of Bij, Bijk, Bijkl
int computeTriNodeIndex (int P, int i, int j)
 computes node index, use getTriNodeIndex to leverage tables
int computeTetNodeIndex (int P, int i, int j, int k)
 computes node index, use getTetNodeIndex to leverage tables
int getNumControlPoints (int type, int order)
 calculates total number of control points, use tables for smaller numbers, this is for quality
int getNumInternalControlPoints (int type, int order)
 calculates number of internal control points, use tables for smaller numbers, this is for quality
void getTriNodesFromTetNodes (int f, int P, apf::NewArray< apf::Vector3 > &tetNodes, apf::NewArray< apf::Vector3 > &triNodes)
 computes nodes of face f from tet
void getTriDetJacNodesFromTetDetJacNodes (int f, int P, apf::NewArray< double > &tetNodes, apf::NewArray< double > &triNodes)
 computes det(Jacobian) nodes of face f from tet
void getFullRepFromBlended (int type, apf::NewArray< double > &transformCoefficients, apf::NewArray< apf::Vector3 > &elemNodes)
 gets full set of bezier control points given blended points
double computeTriJacobianDetFromBezierFormulation (apf::Mesh *m, apf::MeshEntity *e, apf::Vector3 &xi)
 computes det(Jacobian) for tri from the Bezier conversion
double computeTetJacobianDetFromBezierFormulation (apf::Mesh *m, apf::MeshEntity *e, apf::Vector3 &xi)
 computes det(Jacobian) for tri from the Bezier conversion
void getBezierNodeXi (int type, int P, int node, apf::Vector3 &xi)
 get one bezier node location in parameter space
void collectNodeXi (int parentType, int childType, int P, const apf::Vector3 *range, apf::NewArray< apf::Vector3 > &xi)
 get all bezier node locations in parameter space
void elevateBezierEdge (int P, int r, apf::NewArray< apf::Vector3 > &nodes, apf::NewArray< apf::Vector3 > &elevatedNodes)
 elevation functions for beziers
void elevateBezierCurve (apf::Mesh2 *m, apf::MeshEntity *edge, int n, int r)
 Elevate a bezier curve to a higher order.
void subdivideBezierEdge (int P, double t, apf::NewArray< apf::Vector3 > &nodes, apf::NewArray< apf::Vector3 >(&subNodes)[2])
 subdivision functions for beziers
void convertInterpolationPoints (int n, int ne, apf::NewArray< apf::Vector3 > &nodes, apf::NewArray< double > &c, apf::NewArray< apf::Vector3 > &newNodes)
 converts interpolating points to control points
void getBezierTransformationCoefficients (int P, int type, apf::NewArray< double > &c)
 get coefficients for interpolating points to control points
void snapToInterpolate (apf::Mesh2 *m, apf::MeshEntity *e)
 a per entity version of above
void getTransformationMatrix (apf::Mesh *m, apf::MeshEntity *e, mth::Matrix< double > &A, const apf::Vector3 *range)
 compute the matrix to transform between Bezier and Lagrange Points
void BlendedTriangleGetValues (apf::Mesh *m, apf::MeshEntity *e, apf::Vector3 const &xi, apf::NewArray< double > &values)
 shape blending functions
void getBezierTransformationMatrix (int type, int P, mth::Matrix< double > &A, const apf::Vector3 *range)
 Get transformation matrix corresponding to a parametric range.
void getBezierTransformationMatrix (int parentType, int childType, int P, mth::Matrix< double > &A, const apf::Vector3 *childRange)
 Get transformation matrix of a lower entity in a higher entity over a parametric range.
int binomial (int n, int i)
 binomial function n!/(i!(n-i)!)
int trinomial (int n, int i, int j)
 trinomial function n!/(i!j!(n-i-j)!)
int quadnomial (int n, int i, int j, int k)
 "quadnomial" function n!/(i!j!k!(n-i-j-k)!)
double intpow (const double b, const int e)
 faster power for integers
void invertMatrixWithQR (int n, mth::Matrix< double > &A, mth::Matrix< double > &Ai)
 invert a matrix using QR factorization
void invertMatrixWithPLU (int n, mth::Matrix< double > &A, mth::Matrix< double > &Ai)
 invert a matrix using Pivoting and LU decomposition
void subdivideBezierEntityJacobianDet (int P, int type, apf::NewArray< double > &c, apf::NewArray< double > &nodes, apf::NewArray< double > *subNodes)
 subdivide jacobian det using subdivision matrices
void getBezierJacobianDetSubdivisionCoefficients (int P, int type, apf::NewArray< double > &c)
 get matrices used for uniform subdivision, 2^dim matrices, unrolled into a double
void elevateBezierJacobianDet (int type, int P, int r, apf::NewArray< double > &nodes, apf::NewArray< double > &elevatedNodes)
 elevate jacobian det to higher order, used in getQuality
bool isBoundaryEntity (apf::Mesh *m, apf::MeshEntity *e)
 checks if is a boundary entity
void repositionInteriorWithBlended (ma::Mesh *m, ma::Entity *e)
 uses blending to position interior points, based on edge locations
void splitEdges (ma::Adapt *a)
 Split edges marked with ma::SPLIT and place high order nodes using subdivision, see ma::refine.
int markInvalidEntities (Adapt *a)
 mark invalid entities with validity tag
int getTag (Adapt *a, ma::Entity *e)
 get validityTag
void setTag (Adapt *a, ma::Entity *e, int tag)
 set validityTag
void clearTag (Adapt *a, ma::Entity *e)
 reset validityTag
int getValidityTag (ma::Mesh *m, ma::Entity *e, ma::Entity *bdry)
 get validityTag
int fixLargeBoundaryAngles (Adapt *a)
 Take boundary triangles where two edges on the boundary form an angle of 180 (or greater) at a vertex and split the edge opposite them.
int fixInvalidEdges (Adapt *a)
 If an edge is flagged as invalid, try and collapse or swap it away.
void fixCrvElementShapes (Adapt *a)
 attempts to fix the shape of the elements in a same manner as ma::fixElementShape
ma::ShapeHandler * getShapeHandler (ma::Adapt *a)
 get bezier shape handler
void transferParametricOnEdgeSplit (apf::Mesh *m, apf::MeshEntity *e, double t, apf::Vector3 &p)
 gets parametric location on geometry given t in [0,1] on edge
void transferParametricOnTriSplit (apf::Mesh2 *m, apf::MeshEntity *e, apf::Vector3 &t, apf::Vector3 &p)
 gets parametric location on geometry given barycentric coordinate t on tri
void transferParametricOnGeometricEdgeSplit (apf::Mesh2 *m, apf::MeshEntity *e, double t, apf::Vector3 &p)
 gets parametric location on geometry given t in [0,1] on edge, but splitting in geometric space first and then projecting
void transferParametricOnGeometricTriSplit (apf::Mesh2 *m, apf::MeshEntity *e, apf::Vector3 &t, apf::Vector3 &p)
 gets parametric location on geometry given barycentric coordinate t on tri, but splitting in geometric space first and then projecting

Variables

const bezierShape bezier [apf::Mesh::TYPES]
 table of shape functions
const bezierShapeGrads bezierGrads [apf::Mesh::TYPES]
 table of shape function gradients
const SubdivisionFunction subdivideBezierJacobianDet [apf::Mesh::TYPES]
 table of jacobian det subdivision functions
unsigned const *const *const b2 [11]
 table of indices for triangles, b2[order][i][j], only up to 10th order is stored, higher can be generated on the fly
unsigned const *const *const *const b3 [5]
 table of indices for tets, b3[order][i][j][k], only up to 4th order is stored, higher can be generated on the fly
unsigned const *const *const *const tet_tri [7]
 table of alignment used in alignSharedNodes, tet_tri[order][flip][rotate][node];
apf::Vector3 const *const elem_vert_xi [apf::Mesh::TYPES]
 parametric locations of midpoint nodes given a vertex number, elem_vert_xi[type][vertex_index]
apf::Vector3 const *const elem_edge_xi [apf::Mesh::TYPES]
 parametric locations of midpoint nodes given an edge number, elem_edge_xi[type][edge_index]

Detailed Description

the curving functions are contained in this namespace


Function Documentation

void crv::adapt ( ma::Input in  ) 

crv adapt with custom configuration

see maInput.h for details. note that this function will delete the Input object

void crv::BlendedTriangleGetValues ( apf::Mesh m,
apf::MeshEntity *  e,
apf::Vector3 const &  xi,
apf::NewArray< double > &  values 
)

shape blending functions

see bezier.tex in SCOREC/docs repo

void crv::changeMeshOrder ( apf::Mesh2 m,
int  newOrder 
)

change the order of a Bezier Mesh

going up in order is exact, except for boundary elements, where snapping changes things Going down in order is approximate everywhere

int crv::checkValidity ( apf::Mesh m,
apf::MeshEntity *  e,
int  algorithm = 2 
)

checks validity of it and returns integer corresponding to invalid entity. Quality::checkValidity should be used if multiple elements checked in a row

Use an integer to determine the vuality tag 0 -> Not checked 1 -> Okay Quality 2-7 -> Vertex of index+2 is bad 8-13 -> Edge of index+6 is bad 14-17 -> Face of index+12 bad 20 -> Tet itself is bad, this one is the worst

6*dim + 2 + index

double crv::computeTetJacobianDetFromBezierFormulation ( apf::Mesh m,
apf::MeshEntity *  e,
apf::Vector3 xi 
)

computes det(Jacobian) for tri from the Bezier conversion

this evaluates an order d(P-1) bezier to compute it and is much, much slower than the direct method, but exists for comparison

double crv::computeTriJacobianDetFromBezierFormulation ( apf::Mesh m,
apf::MeshEntity *  e,
apf::Vector3 xi 
)

computes det(Jacobian) for tri from the Bezier conversion

this evaluates an order d(P-1) bezier to compute it and is much, much slower than the direct method, but exists for comparison

void crv::convertInterpolationPoints ( int  n,
int  ne,
apf::NewArray< apf::Vector3 > &  nodes,
apf::NewArray< double > &  c,
apf::NewArray< apf::Vector3 > &  newNodes 
)

converts interpolating points to control points

n is total number of nodes on the shape ne is nodes on the entity, that belong to it c is a coefficient matrix in vector form corresponding to the matrix

void crv::elevateBezierCurve ( apf::Mesh2 m,
apf::MeshEntity *  edge,
int  n,
int  r 
)

Elevate a bezier curve to a higher order.

This elevates from nth order to n+rth order requires the curve be order n+r in memory already, and that the first n points correspond to the lower order curve

apf::FieldShape* crv::getBezier ( int  order  ) 

Get the Bezier Curve or Shape of some order.

goes from first to sixth order

void crv::getBezierTransformationCoefficients ( int  P,
int  type,
apf::NewArray< double > &  c 
)

get coefficients for interpolating points to control points

works only for prescribed optimal point locations up to 6th order in 2D and

void crv::getBezierTransformationMatrix ( int  type,
int  P,
mth::Matrix< double > &  A,
const apf::Vector3 range 
)

Get transformation matrix corresponding to a parametric range.

Range is an array of size(num vertices), this is used for subdivision, refinement. It is the element transformation matrix, see notes

void crv::getBezierTransformationMatrix ( int  parentType,
int  childType,
int  P,
mth::Matrix< double > &  A,
const apf::Vector3 childRange 
)

Get transformation matrix of a lower entity in a higher entity over a parametric range.

Range is an array of size(num vertices), this is used for refinement, It is the lower dimensional component, as part of the higher array, see notes

void crv::getFullRepFromBlended ( int  type,
apf::NewArray< double > &  transformCoefficients,
apf::NewArray< apf::Vector3 > &  elemNodes 
)

gets full set of bezier control points given blended points

this is used for quality assessment of blended shapes, and reallocates elemNodes

int crv::getNumControlPoints ( int  type,
int  order 
)

calculates total number of control points, use tables for smaller numbers, this is for quality

This gives the numbers for full bezier shapes, this is not accurate for blended shapes

int crv::getNumInternalControlPoints ( int  type,
int  order 
)

calculates number of internal control points, use tables for smaller numbers, this is for quality

This gives the numbers for full bezier shapes, this is not accurate for blended shapes

void crv::getTransformationMatrix ( apf::Mesh m,
apf::MeshEntity *  e,
mth::Matrix< double > &  A,
const apf::Vector3 range 
)

compute the matrix to transform between Bezier and Lagrange Points

this is a support function, not actual ever needed. Bezier control points, C, can be written as L = A*C where A is a matrix of Bernstein polynomials and binomial coefficients. To compute the control points from Lagrange points, A^{-1} is used. crvBezierPoints.cc contains A^{-1}, precomputed, for the nodeXi locations in getBezierNodeXi.

void crv::getTriNodesFromTetNodes ( int  f,
int  P,
apf::NewArray< apf::Vector3 > &  tetNodes,
apf::NewArray< apf::Vector3 > &  triNodes 
)

computes nodes of face f from tet

does not consider alignment, extra care needed

int crv::getValidityTag ( ma::Mesh m,
ma::Entity e,
ma::Entity bdry 
)

get validityTag

Use an integer to determine the validity tag 0 -> Not checked 1 -> Okay Quality 2-7 -> Vertices are bad 8-13 -> Edge is are bad 14-17 -> Face is are bad 20 -> Tet itself is bad, this one is the worst

6*dim + 2 + index

double crv::interpolationError ( apf::Mesh m,
apf::MeshEntity *  e,
int  n 
)

computes interpolation error of a curved entity on a mesh

this computes the Hausdorff distance by sampling n points per dimension of the entity through uniform sampling locations in parameter space

int crv::markInvalidEntities ( Adapt *  a  ) 

mark invalid entities with validity tag

since validity checking is expensive, do this as little as possible and keep information about the check

void crv::stats ( ma::Input in,
std::vector< double > &  edgeLengths,
std::vector< double > &  linearQualities,
std::vector< double > &  curvedQualities,
bool  inMetric = true 
)

crv stats to get statistic information about the mesh

statistic considered are (1)final/desired edge-lengths (2) linear quality (3) curved quality (minJ/maxJ)

void crv::subdivideBezierEntityJacobianDet ( int  P,
int  type,
apf::NewArray< double > &  c,
apf::NewArray< double > &  nodes,
apf::NewArray< double > *  subNodes 
)

subdivide jacobian det using subdivision matrices

see getBezierJacobianDetSubdivisionCoefficients

void crv::transferParametricOnTriSplit ( apf::Mesh2 m,
apf::MeshEntity *  e,
apf::Vector3 t,
apf::Vector3 p 
)

gets parametric location on geometry given barycentric coordinate t on tri

uses two linear splits to find location on a triangle, which is more stable and handles degeneracy better than using pure barycentrics

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Defines