SCOREC core
Parallel unstructured mesh tools
Classes | Functions
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...
 

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
 
void interpolatingToBezier (apf::Mesh2 *m)
 converts Interpolating nodes to Control points for a Bezier 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 More...
 
void adapt (const ma::Input *in)
 crv adapt with custom configuration More...
 
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 More...
 
apf::FieldShapegetBezier (int order)
 Get the Bezier Curve or Shape of some order. More...
 
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 More...
 
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 More...
 
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
 

Detailed Description

the curving functions are contained in this namespace

Function Documentation

◆ adapt() [1/2]

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

crv adapt with custom configuration

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

◆ adapt() [2/2]

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

◆ checkValidity()

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

◆ getBezier()

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

Get the Bezier Curve or Shape of some order.

goes from first to sixth order

◆ interpolationError()

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

◆ stats()

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)