SCOREC core
Parallel unstructured mesh tools
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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 
46 {
47  public:
48  MeshCurver(apf::Mesh2* m, int P) : m_mesh(m), m_order(P) {};
49  virtual ~MeshCurver() {};
50  virtual bool run() = 0;
51 
53  void snapToInterpolate(int dim);
54 
56  void synchronize();
57 
58  protected:
59  apf::Mesh2* m_mesh;
60  int m_order;
61 };
62 
68 {
69  public:
70  InterpolatingCurver(apf::Mesh2* m, int P) : MeshCurver(m,P) {};
71  virtual ~InterpolatingCurver() {};
72  virtual bool run();
73 
74 };
75 
80 class BezierCurver : public MeshCurver
81 {
82  public:
83  BezierCurver(apf::Mesh2* m, int P, int B) : MeshCurver(m,P)
84  {
86  };
87 
91  virtual bool run();
94 };
95 
100 {
101  public:
102  GregoryCurver(apf::Mesh2* m, int P, int B)
103  : BezierCurver(m,P,B) {};
105  virtual bool run();
110 };
111 
114  ma::Mesh* m, ma::SizeField* f=0,
115  ma::SolutionTransfer* s=0);
116 
120 void adapt(ma::Input* in);
121 
126 void stats(ma::Input* in,
127  std::vector<double> &edgeLengths,
128  std::vector<double> &linearQualities,
129  std::vector<double> &curvedQualities,
130  bool inMetric=true);
131 
132 
135 apf::FieldShape* getBezier(int order);
138 
141 double getQuality(apf::Mesh* m,apf::MeshEntity* e);
142 
155 int checkValidity(apf::Mesh* m, apf::MeshEntity* e,
156  int algorithm = 2);
159 class Quality
160 {
161 public:
162 /* \brief three options for algorithm:
163  * 0 - subdivision
164  * 1 - elevation
165  * 2 - subdivision, using matrices */
166  Quality(apf::Mesh* m, int algorithm_);
167  virtual ~Quality() {};
169  virtual double getQuality(apf::MeshEntity* e) = 0;
171  virtual int checkValidity(apf::MeshEntity* e) = 0;
172 protected:
173  apf::Mesh* mesh;
174  int algorithm;
175  int order;
176 };
178 Quality* makeQuality(apf::Mesh* m, int algorithm = 2);
179 
184 double interpolationError(apf::Mesh* m, apf::MeshEntity* e, int n);
185 
189 void writeCurvedVtuFiles(apf::Mesh* m, int type, int n, const char* prefix);
190 
194 void writeCurvedWireFrame(apf::Mesh* m, int n, const char* prefix);
195 
197 void writeControlPointVtuFiles(apf::Mesh* m, const char* prefix);
200 void writeInterpolationPointVtuFiles(apf::Mesh* m, const char* prefix);
201 
203 int getTriNodeIndex(int P, int i, int j);
204 int getTetNodeIndex(int P, int i, int j, int k);
205 
207 void fail(const char* why) __attribute__((noreturn));
208 
209 }
210 
211 #endif
The MeshAdapt interface.
templated math functions
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 -&gt; better r...
user-defined solution transfer base
int getTriNodeIndex(int P, int i, int j)
publically accessible functions
void adapt(ma::Input *in)
crv adapt with custom configuration
double getQuality(apf::Mesh *m, apf::MeshEntity *e)
computes min det Jacobian / max det Jacobian. Quality::getQuality should be used if multiple elements...
placeholder to set array sizes
Definition: apfMesh.h:166
void setCubicEdgePointsUsingNormals()
sets cubic edge points using normals
this curves a mesh with Bezier shapes
Definition: crv.h:80
void writeCurvedWireFrame(apf::Mesh *m, int n, const char *prefix)
Visualization, writes wireframe of the curved mesh, n is number of subdivisions, higher number -&gt; bet...
Describes field distribution and shape functions.
Definition: apfShape.h:73
ma::Input * configureShapeCorrection(ma::Mesh *m, ma::SizeField *f=0, ma::SolutionTransfer *s=0)
configure for fixing invalid elements
virtual double getQuality(apf::MeshEntity *e)=0
get scaled jacobian, a quality measure
int getBlendingOrder(const int type)
gets the blending order
void setBlendingOrder(const int type, const int b)
sets the blending order, if shape blending is used
The APF Mesh modification interface.
virtual bool run()
curves a mesh using G1 gregory surfaces, see crvBezier.cc
int getOrder()
gets order used in bezier shape functions
User configuration for a MeshAdapt run.
Definition: maInput.h:34
void setInternalPointsLocally()
sets internal points locally
double interpolationError(apf::Mesh *m, apf::MeshEntity *e, int n)
computes interpolation error of a curved entity on a mesh
Interface to a mesh part.
Definition: apfMesh.h:103
void setOrder(const int order)
sets order used in bezier shape functions
void writeControlPointVtuFiles(apf::Mesh *m, const char *prefix)
Visualization, writes file of control nodes for each entity.
Quality * makeQuality(apf::Mesh *m, int algorithm=2)
use this to make a quality object with the correct dimension
int countNumberInvalidElements(apf::Mesh2 *m)
count invalid elements of the mesh
Field node distribution and shape functions.
Base Mesh curving object.
Definition: crv.h:45
void fail(const char *why) __attribute__((noreturn))
crv fail function
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...
class to store matrices used in quality assessment and validity checking
Definition: crv.h:159
Extended mesh interface for modification.
Definition: apfMesh2.h:29
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
virtual bool run()
curves a mesh using bezier curves of chosen order
apf::FieldShape * getBezier(int order)
Get the Bezier Curve or Shape of some order.
virtual int checkValidity(apf::MeshEntity *e)=0
check the validity (det(Jacobian) &gt; eps) of an element
void snapToInterpolate(int dim)
snaps points to interpolating locations
apf::FieldShape * getGregory()
Get the 4th order Gregory Surface.
this curves a mesh with 4th order G1 Patches
Definition: crv.h:99
curves an already changed mesh
Definition: crv.h:67
void synchronize()
wrapper around synchronizeFieldData
void writeInterpolationPointVtuFiles(apf::Mesh *m, const char *prefix)
Visualization, writes file of shapes evaluated at node xi for each entity.
void convertInterpolatingToBezier()
converts interpolating points to bezier control points