SCOREC core
Parallel unstructured mesh tools
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
crvBezier.h
Go to the documentation of this file.
1 /*
2  * Copyright 2015 Scientific Computation Research Center
3  *
4  * This work is open source software, licensed under the terms of the
5  * BSD license as described in the LICENSE file in the top-level directory.
6  */
7 
8 #ifndef CRVBEZIER_H
9 #define CRVBEZIER_H
10 
11 #include "crv.h"
12 #include "mth.h"
13 
17 namespace crv {
18 
20 double Bij(const int i, const int j,const double u, const double v);
21 double Bijk(const int i, const int j, const int k, const double u,
22  const double v, const double w);
23 double Bijkl(const int i, const int j, const int k, const int l,
24  const double u, const double v, const double w, const double t);
25 
27 double Bij(const int ij[], const double xi[]);
28 double Bijk(const int ijk[], const double xi[]);
29 double Bijkl(const int ijkl[], const double xi[]);
30 
32 int computeTriNodeIndex(int P, int i, int j);
34 int computeTetNodeIndex(int P, int i, int j, int k);
35 
40 int getNumControlPoints(int type, int order);
45 int getNumInternalControlPoints(int type, int order);
46 
49 void getTriNodesFromTetNodes(int f, int P,
50  apf::NewArray<apf::Vector3>& tetNodes,
51  apf::NewArray<apf::Vector3>& triNodes);
53 void getTriDetJacNodesFromTetDetJacNodes(int f, int P,
54  apf::NewArray<double>& tetNodes,
55  apf::NewArray<double>& triNodes);
59 void getFullRepFromBlended(int type,
60  apf::NewArray<double>& transformCoefficients,
61  apf::NewArray<apf::Vector3>& elemNodes);
67  apf::MeshEntity* e, apf::Vector3& xi);
73  apf::MeshEntity* e, apf::Vector3& xi);
74 
76 void getBezierNodeXi(int type, int P, int node, apf::Vector3& xi);
78 void collectNodeXi(int parentType, int childType, int P,
79  const apf::Vector3* range, apf::NewArray<apf::Vector3>& xi);
81 void elevateBezierEdge(int P, int r, apf::NewArray<apf::Vector3>& nodes,
82  apf::NewArray<apf::Vector3>& elevatedNodes);
83 void elevateBezierTriangle(int P, int r, apf::NewArray<apf::Vector3>& nodes,
84  apf::NewArray<apf::Vector3>& elevatedNodes);
85 void elevateBezierTet(int P, int r, apf::NewArray<apf::Vector3>& nodes,
86  apf::NewArray<apf::Vector3>& elevatedNodes);
87 void elevateBezier(int type, int P, int r,
88  apf::NewArray<apf::Vector3>& nodes,
89  apf::NewArray<apf::Vector3>& elevatedNodes);
90 
95 void elevateBezierCurve(apf::Mesh2* m, apf::MeshEntity* edge, int n, int r);
96 
98 void subdivideBezierEdge(int P, double t, apf::NewArray<apf::Vector3>& nodes,
99  apf::NewArray<apf::Vector3> (&subNodes)[2]);
100 void subdivideBezierTriangle(int P, apf::Vector3& p,
101  apf::NewArray<apf::Vector3>& nodes,
102  apf::NewArray<apf::Vector3> (&subNodes)[3]);
103 void subdivideBezierTriangle(int P, apf::NewArray<apf::Vector3>& nodes,
104  apf::NewArray<apf::Vector3> (&subNodes)[4]);
105 void subdivideBezierTet(int P, apf::Vector3& p,
106  apf::NewArray<apf::Vector3>& nodes,
107  apf::NewArray<apf::Vector3> (&subNodes)[4]);
108 
114 void convertInterpolationPoints(int n, int ne,
115  apf::NewArray<apf::Vector3>& nodes,
116  apf::NewArray<double>& c,
117  apf::NewArray<apf::Vector3>& newNodes);
118 
119 void convertInterpolationPoints(apf::Mesh2* m, apf::MeshEntity* e,
120  int n, int ne, apf::NewArray<double>& c);
121 
125 void getBezierTransformationCoefficients(int P, int type,
126  apf::NewArray<double>& c);
127 void getInternalBezierTransformationCoefficients(apf::Mesh* m, int P, int blend,
128  int type, apf::NewArray<double>& c);
129 void getGregoryTransformationCoefficients(int type, apf::NewArray<double>& c);
130 void getGregoryBlendedTransformationCoefficients(int blend, int type,
131  apf::NewArray<double>& c);
132 
134 void snapToInterpolate(apf::Mesh2* m, apf::MeshEntity* e, bool isNew = false);
135 
143 void getTransformationMatrix(apf::Mesh* m, apf::MeshEntity* e,
144  mth::Matrix<double>& A, const apf::Vector3 *range);
145 
146 }
147 
148 #endif
templated math functions
compile-time (static) matrix
Definition: mthMatrix.h:23
double Bij(const int i, const int j, const double u, const double v)
polynomial part of bernstein polynomial, Bij, Bijk, Bijkl
void getBezierTransformationCoefficients(int P, int type, apf::NewArray< double > &c)
get coefficients for interpolating points to control points
int computeTriNodeIndex(int P, int i, int j)
computes node index, use getTriNodeIndex to leverage tables
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 getBezierNodeXi(int type, int P, int node, apf::Vector3 &xi)
get one bezier node location in parameter space
int getNumInternalControlPoints(int type, int order)
calculates number of internal control points, use tables for smaller numbers, this is for quality ...
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
Interface to a mesh part.
Definition: apfMesh.h:103
void elevateBezierCurve(apf::Mesh2 *m, apf::MeshEntity *edge, int n, int r)
Elevate a bezier curve to a higher order.
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
int getNumControlPoints(int type, int order)
calculates total number of control points, use tables for smaller numbers, this is for quality ...
convenience wrapper over apf::Vector&lt;3&gt;
Definition: apfVector.h:150
void subdivideBezierEdge(int P, double t, apf::NewArray< apf::Vector3 > &nodes, apf::NewArray< apf::Vector3 >(&subNodes)[2])
subdivision functions for beziers
void snapToInterpolate(apf::Mesh2 *m, apf::MeshEntity *e, bool isNew=false)
a per entity version of above
Extended mesh interface for modification.
Definition: apfMesh2.h:29
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
int computeTetNodeIndex(int P, int i, int j, int k)
computes node index, use getTetNodeIndex to leverage tables
void getTriDetJacNodesFromTetDetJacNodes(int f, int P, apf::NewArray< double > &tetNodes, apf::NewArray< double > &triNodes)
computes det(Jacobian) nodes of face f from tet
void elevateBezierEdge(int P, int r, apf::NewArray< apf::Vector3 > &nodes, apf::NewArray< apf::Vector3 > &elevatedNodes)
elevation functions for beziers
void getTriNodesFromTetNodes(int f, int P, apf::NewArray< apf::Vector3 > &tetNodes, apf::NewArray< apf::Vector3 > &triNodes)
computes 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
main file for curved element support