DOLFIN
DOLFIN C++ interface
MeshQuality.h
1 // Copyright (C) 2013 Garth N. Wells
2 //
3 // This file is part of DOLFIN.
4 //
5 // DOLFIN is free software: you can redistribute it and/or modify
6 // it under the terms of the GNU Lesser General Public License as published by
7 // the Free Software Foundation, either version 3 of the License, or
8 // (at your option) any later version.
9 //
10 // DOLFIN is distributed in the hope that it will be useful,
11 // but WITHOUT ANY WARRANTY; without even the implied warranty of
12 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 // GNU Lesser General Public License for more details.
14 //
15 // You should have received a copy of the GNU Lesser General Public License
16 // along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
17 //
18 // First added: 2013-10-07
19 // Last changed:
20 
21 #ifndef __MESH_QUALITY_H
22 #define __MESH_QUALITY_H
23 
24 #include <string>
25 #include <utility>
26 #include <vector>
27 #include <boost/multi_array.hpp>
28 #include <memory>
29 #include "Cell.h"
30 
31 namespace dolfin
32 {
33 
34  class Mesh;
35 
37 
39  {
40  public:
41 
51  radius_ratios(std::shared_ptr<const Mesh> mesh);
52 
61  aspect_ratio_gamma(std::shared_ptr<const Mesh> mesh);
62 
72  static std::pair<double, double> radius_ratio_min_max(const Mesh& mesh);
73 
74 
80  static std::pair<std::vector<double>, std::vector<double> >
82  std::size_t num_bins = 50);
83 
88  static std::string
90  std::size_t num_intervals = 50);
91 
93  static void dihedral_angles(const Cell& cell, std::vector<double>& dh_angle);
94 
96  static std::pair<double, double>
97  dihedral_angles_min_max(const Mesh& mesh);
98 
101  static std::pair<std::vector<double>, std::vector<double> >
103  std::size_t num_bins = 100);
104 
106  static std::string
108  std::size_t num_intervals = 100);
109  };
110 
111 }
112 
113 #endif
Definition: adapt.h:41
Definition: adapt.h:29
static MeshFunction< double > aspect_ratio_gamma(std::shared_ptr< const Mesh > mesh)
Definition: MeshQuality.cpp:46
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::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
The class provides functions to quantify mesh quality.
Definition: MeshQuality.h:38
static std::string radius_ratio_matplotlib_histogram(const Mesh &mesh, std::size_t num_intervals=50)
Definition: MeshQuality.cpp:116
static MeshFunction< double > radius_ratios(std::shared_ptr< const Mesh > mesh)
Definition: MeshQuality.cpp:33
A Cell is a MeshEntity of topological codimension 0.
Definition: Cell.h:42
static std::pair< double, double > radius_ratio_min_max(const Mesh &mesh)
Definition: MeshQuality.cpp:68
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
static void dihedral_angles(const Cell &cell, std::vector< double > &dh_angle)
Get internal dihedral angles of a tetrahedral cell.
Definition: MeshQuality.cpp:167
Definition: Mesh.h:82
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