SCOREC core
Parallel unstructured mesh tools
crv.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 CRV_H
9 #define CRV_H
10 
11 #include "apfMesh2.h"
12 #include "apfShape.h"
13 #include <ma.h>
14 #include <mth.h>
15 #include <stdio.h>
16 #include <vector>
17 
23 namespace crv {
24 
26 static unsigned const MAX_ORDER = 19;
27 
29 void setOrder(const int order);
31 int getOrder();
33 void setBlendingOrder(const int type, const int b);
35 int getBlendingOrder(const int type);
36 
39 
42 
49 {
50  public:
51  MeshCurver(apf::Mesh2* m, int P) : m_mesh(m), m_order(P) {};
52  virtual ~MeshCurver() {};
53  virtual bool run() = 0;
54 
56  void snapToInterpolate(int dim);
57 
59  void synchronize();
60 
61  protected:
62  apf::Mesh2* m_mesh;
63  int m_order;
64 };
65 
71 {
72  public:
73  InterpolatingCurver(apf::Mesh2* m, int P) : MeshCurver(m,P) {};
74  virtual ~InterpolatingCurver() {};
75  virtual bool run();
76 
77 };
78 
83 class BezierCurver : public MeshCurver
84 {
85  public:
86  BezierCurver(apf::Mesh2* m, int P, int B) : MeshCurver(m,P)
87  {
89  };
90 
94  virtual bool run();
97 };
98 
103 {
104  public:
105  GregoryCurver(apf::Mesh2* m, int P, int B)
106  : BezierCurver(m,P,B) {};
108  virtual bool run();
113 };
114 
117  ma::Mesh* m, ma::SizeField* f=0,
118  ma::SolutionTransfer* s=0);
119 
123 void adapt(ma::Input* in);
124 
128 void adapt(const ma::Input* in);
129 
134 void stats(ma::Input* in,
135  std::vector<double> &edgeLengths,
136  std::vector<double> &linearQualities,
137  std::vector<double> &curvedQualities,
138  bool inMetric=true);
139 
140 
146 
149 double getQuality(apf::Mesh* m,apf::MeshEntity* e);
150 
163 int checkValidity(apf::Mesh* m, apf::MeshEntity* e,
164  int algorithm = 2);
167 class Quality
168 {
169 public:
170 /* \brief three options for algorithm:
171  * 0 - subdivision
172  * 1 - elevation
173  * 2 - subdivision, using matrices */
174  Quality(apf::Mesh* m, int algorithm_);
175  virtual ~Quality() {};
177  virtual double getQuality(apf::MeshEntity* e) = 0;
179  virtual int checkValidity(apf::MeshEntity* e) = 0;
180 protected:
181  apf::Mesh* mesh;
182  int algorithm;
183  int order;
184 };
186 Quality* makeQuality(apf::Mesh* m, int algorithm = 2);
187 
192 double interpolationError(apf::Mesh* m, apf::MeshEntity* e, int n);
193 
197 void writeCurvedVtuFiles(apf::Mesh* m, int type, int n, const char* prefix);
198 
202 void writeCurvedWireFrame(apf::Mesh* m, int n, const char* prefix);
203 
205 void writeControlPointVtuFiles(apf::Mesh* m, const char* prefix);
208 void writeInterpolationPointVtuFiles(apf::Mesh* m, const char* prefix);
209 
211 int getTriNodeIndex(int P, int i, int j);
212 int getTetNodeIndex(int P, int i, int j, int k);
213 
215 void fail(const char* why) __attribute__((noreturn));
216 
217 }
218 
219 #endif
The APF Mesh modification interface.
Field node distribution and shape functions.
Describes field distribution and shape functions.
Definition: apfShape.h:92
Extended mesh interface for modification.
Definition: apfMesh2.h:30
Interface to a mesh part.
Definition: apfMesh.h:105
@ TYPES
placeholder to set array sizes
Definition: apfMesh.h:167
this curves a mesh with Bezier shapes
Definition: crv.h:84
void convertInterpolatingToBezier()
converts interpolating points to bezier control points
virtual bool run()
curves a mesh using bezier curves of chosen order
this curves a mesh with 4th order G1 Patches
Definition: crv.h:103
void setInternalPointsLocally()
sets internal points locally
virtual bool run()
curves a mesh using G1 gregory surfaces, see crvBezier.cc
void setCubicEdgePointsUsingNormals()
sets cubic edge points using normals
curves an already changed mesh
Definition: crv.h:71
Base Mesh curving object.
Definition: crv.h:49
void snapToInterpolate(int dim)
snaps points to interpolating locations
void synchronize()
wrapper around synchronizeFieldData
class to store matrices used in quality assessment and validity checking
Definition: crv.h:168
virtual double getQuality(apf::MeshEntity *e)=0
get scaled jacobian, a quality measure
virtual int checkValidity(apf::MeshEntity *e)=0
check the validity (det(Jacobian) > eps) of an element
User configuration for a MeshAdapt run.
Definition: maInput.h:35
user-defined solution transfer base
The MeshAdapt interface.
templated math functions
the curving functions are contained in this namespace
double interpolationError(apf::Mesh *m, apf::MeshEntity *e, int n)
computes interpolation error of a curved entity on a mesh
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
void setOrder(const int order)
sets order used in bezier shape functions
apf::FieldShape * getGregory()
Get the 4th order Gregory Surface.
int countNumberInvalidElements(apf::Mesh2 *m)
count invalid elements of the mesh
void setBlendingOrder(const int type, const int b)
sets the blending order, if shape blending is used
double getQuality(apf::Mesh *m, apf::MeshEntity *e)
computes min det Jacobian / max det Jacobian. Quality::getQuality should be used if multiple elements...
void fail(const char *why) __attribute__((noreturn))
crv fail function
void adapt(ma::Input *in)
crv adapt with custom configuration
int getTriNodeIndex(int P, int i, int j)
publically accessible functions
int getOrder()
gets order used in bezier shape functions
void writeInterpolationPointVtuFiles(apf::Mesh *m, const char *prefix)
Visualization, writes file of shapes evaluated at node xi for each entity.
apf::FieldShape * getBezier(int order)
Get the Bezier Curve or Shape of some order.
Quality * makeQuality(apf::Mesh *m, int algorithm=2)
use this to make a quality object with the correct dimension
int checkValidity(apf::Mesh *m, apf::MeshEntity *e, int algorithm=2)
checks validity of it and returns integer corresponding to invalid entity. Quality::checkValidity sho...
void interpolatingToBezier(apf::Mesh2 *m)
converts Interpolating nodes to Control points for a Bezier mesh
void writeCurvedWireFrame(apf::Mesh *m, int n, const char *prefix)
Visualization, writes wireframe of the curved mesh, n is number of subdivisions, higher number -> bet...
ma::Input * configureShapeCorrection(ma::Mesh *m, ma::SizeField *f=0, ma::SolutionTransfer *s=0)
configure for fixing invalid elements
int getBlendingOrder(const int type)
gets the blending order
void writeControlPointVtuFiles(apf::Mesh *m, const char *prefix)
Visualization, writes file of control nodes for each entity.
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 r...