DOLFIN
DOLFIN C++ interface
|
21 #ifndef __MESH_QUALITY_H
22 #define __MESH_QUALITY_H
27 #include <boost/multi_array.hpp>
80 static std::pair<std::vector<double>, std::vector<double> >
82 std::size_t num_bins = 50);
90 std::size_t num_intervals = 50);
96 static std::pair<double, double>
101 static std::pair<std::vector<double>, std::vector<double> >
103 std::size_t num_bins = 100);
108 std::size_t num_intervals = 100);
static std::string dihedral_angles_matplotlib_histogram(const Mesh &mesh, std::size_t num_intervals=100)
Create Matplotlib string to plot dihedral angles quality histogram.
Definition: MeshQuality.cpp:278
static std::pair< double, double > dihedral_angles_min_max(const Mesh &mesh)
Get internal minimum and maximum dihedral angles of a 3D mesh.
Definition: MeshQuality.cpp:200
static void dihedral_angles(const Cell &cell, std::vector< double > &dh_angle)
Get internal dihedral angles of a tetrahedral cell.
Definition: MeshQuality.cpp:167
A Cell is a MeshEntity of topological codimension 0.
Definition: Cell.h:43
The class provides functions to quantify mesh quality.
Definition: MeshQuality.h:39
static MeshFunction< double > radius_ratios(std::shared_ptr< const Mesh > mesh)
Definition: MeshQuality.cpp:33
static MeshFunction< double > aspect_ratio_gamma(std::shared_ptr< const Mesh > mesh)
Definition: MeshQuality.cpp:46
static std::string radius_ratio_matplotlib_histogram(const Mesh &mesh, std::size_t num_intervals=50)
Definition: MeshQuality.cpp:116
static std::pair< std::vector< double >, std::vector< double > > dihedral_angles_histogram_data(const Mesh &mesh, std::size_t num_bins=100)
Definition: MeshQuality.cpp:234
Definition: MeshFunction.h:58
static std::pair< std::vector< double >, std::vector< double > > radius_ratio_histogram_data(const Mesh &mesh, std::size_t num_bins=50)
Definition: MeshQuality.cpp:85
static std::pair< double, double > radius_ratio_min_max(const Mesh &mesh)
Definition: MeshQuality.cpp:68