Hermes2D  3.0
mesh.h
1 // This file is part of Hermes2D.
2 //
3 // Hermes2D is free software: you can redistribute it and/or modify
4 // it under the terms of the GNU General Public License as published by
5 // the Free Software Foundation, either version 2 of the License, or
6 // (at your option) any later version.
7 //
8 // Hermes2D is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11 // GNU General Public License for more details.
12 //
13 // You should have received a copy of the GNU General Public License
14 // along with Hermes2D. If not, see <http://www.gnu.org/licenses/>.
15 
16 #ifndef __H2D_MESH_H
17 #define __H2D_MESH_H
18 
19 #include "element.h"
20 #include "mesh_util.h"
21 #include "hash.h"
22 #include "../mixins2d.h"
23 
24 namespace Hermes
25 {
26  namespace Hermes2D
27  {
28  class Element;
29  class HashTable;
30 
31  template<typename Scalar> class Space;
32  template<typename Scalar> class KellyTypeAdapt;
33  struct MItem;
34  struct Rect;
35  extern HERMES_API unsigned g_mesh_seq;
36 
37  namespace RefinementSelectors
38  {
39  template<typename Scalar> class Selector;
40  template<typename Scalar> class POnlySelector;
41  template<typename Scalar> class HOnlySelector;
42  template<typename Scalar> class OptimumSelector;
43  template<typename Scalar> class ProjBasedSelector;
44  template<typename Scalar> class H1ProjBasedSelector;
45  template<typename Scalar> class L2ProjBasedSelector;
46  template<typename Scalar> class HcurlProjBasedSelector;
47  }
61  class HERMES_API Mesh : public HashTable,
62  public Hermes::Mixins::StateQueryable
63  {
64  public:
65  Mesh();
66  virtual ~Mesh();
67 
70  void init(int size = H2D_DEFAULT_HASH_SIZE);
71 
73  virtual bool isOkay() const;
74  inline std::string getClassName() const { return "Mesh"; }
75 
77  bool rescale(double x_ref, double y_ref);
78 
80  void copy(MeshSharedPtr mesh);
81 
83  void copy_base(MeshSharedPtr mesh);
84 
86  void copy_converted(MeshSharedPtr mesh);
87 
89  void create(int nv, double2* verts, int nt, int3* tris, std::string* tri_markers,
90  int nq, int4* quads, std::string* quad_markers, int nm, int2* mark, std::string* boundary_markers);
91 
92 #pragma region MeshHashGrid
93  Element* element_on_physical_coordinates(double x, double y);
99 
100  MeshHashGrid* meshHashGrid;
101 #pragma endregion
102 
103 #pragma region MarkerArea
104  double get_marker_area(int marker);
105 
106  std::map<int, MarkerArea*> marker_areas;
107 #pragma endregion
108 
109 #pragma region getters
110  Element* get_element(int id) const;
112 
114  int get_num_elements() const;
115 
117  int get_num_base_elements() const;
118 
120  int get_num_used_base_elements() const;
121 
123  int get_num_active_elements() const;
124 
126  int get_max_element_id() const;
127 
129  int get_num_vertex_nodes() const;
130 
132  int get_num_edge_nodes() const;
133 
139  void get_bounding_box(double& bottom_left_x, double& bottom_left_y, double& top_right_x, double& top_right_y);
140 
142  Element* get_element_fast(int id) const;
143 
145  unsigned get_seq() const;
146 #pragma endregion
147 
148 #pragma region refinements
149  void refine_element_id(int id, int refinement = 0);
156 
159  void refine_all_elements(int refinement = 0, bool mark_as_initial = false);
160 
168  void refine_by_criterion(int(*criterion)(Element* e), int depth = 1, bool mark_as_initial = false);
169 
172  void refine_towards_vertex(int vertex_id, int depth = 1, bool mark_as_initial = false);
173 
178  void refine_towards_boundary(std::vector<std::string> markers, int depth = 1, bool aniso = true, bool mark_as_initial = false);
179  void refine_towards_boundary(std::string marker, int depth = 1, bool aniso = true, bool mark_as_initial = false);
180 
182  void refine_in_area(std::string marker, int depth = 1, int refinement = 0, bool mark_as_initial = false);
184  void refine_in_areas(std::vector<std::string> markers, int depth = 1, int refinement = 0, bool mark_as_initial = false);
185 
194  int* regularize(int n);
195 
198  void unrefine_element_id(int id);
199 
205  void unrefine_all_elements(bool keep_initial_refinements = true);
206 #pragma endregion
207 
209  class HERMES_API ReferenceMeshCreator
210  {
211  public:
218  ReferenceMeshCreator(MeshSharedPtr coarse_mesh, int refinement = 0);
219 
222  virtual MeshSharedPtr create_ref_mesh();
223 
224  private:
226  MeshSharedPtr coarse_mesh;
227  int refinement;
228  };
229 
230  class HERMES_API MarkersConversion
231  {
232  public:
234 
239  int insert_marker(std::string user_marker);
240 
242  struct StringValid
243  {
244  StringValid();
245  StringValid(std::string marker, bool valid);
246  std::string marker;
247  bool valid;
248  };
249 
251  struct IntValid
252  {
253  IntValid();
254  IntValid(int marker, bool valid);
255  int marker;
256  bool valid;
257  };
258 
261  StringValid get_user_marker(int internal_marker) const;
262 
264  IntValid get_internal_marker(std::string user_marker) const;
265 
266  enum MarkersConversionType {
267  HERMES_ELEMENT_MARKERS_CONVERSION = 0,
268  HERMES_BOUNDARY_MARKERS_CONVERSION = 1
269  };
270 
271  virtual MarkersConversionType get_type() const = 0;
272 
273  int size() const;
274 
275  protected:
280 
283  std::map<int, std::string> conversion_table;
284 
287  std::map<std::string, int> conversion_table_inverse;
288 
289  friend class Space < double > ;
290  friend class Space < std::complex<double> > ;
291  friend class Mesh;
292  };
293 
295  void free();
296 
298  void set_seq(unsigned seq);
299 
300  int nbase, ntopvert, ninitial, nactive;
301 
303  {
304  public:
306  virtual MarkersConversionType get_type() const;
307  };
308 
310  {
311  public:
313  virtual MarkersConversionType get_type() const;
314  };
315 
316  ElementMarkersConversion element_markers_conversion;
317  BoundaryMarkersConversion boundary_markers_conversion;
318  Array<Element> elements;
319 
320  unsigned seq;
321 
323  void initial_single_check();
324 
325  private:
329  void convert_quads_to_triangles();
330 
332  void convert_to_base();
333 
334  void refine_element_to_quads_id(int id);
335 
336  void refine_triangle_to_quads(Element* e, Element** elems_out = nullptr);
337 
338  void refine_element_to_triangles_id(int id);
339 
340  void refine_quad_to_triangles(Element* e);
341 
346  void refine_quad_to_quads(Element* e, int refinement = 0);
347 
348  void convert_element_to_base_id(int id);
349  void convert_triangles_to_base(Element* e);
350  void convert_quads_to_base(Element* e);
351 
353  double bottom_left_x, bottom_left_y, top_right_x, top_right_y;
355  bool bounding_box_calculated;
357  void calc_bounding_box();
358 
359  void unrefine_element_internal(Element* e);
360 
361  int* parents;
362  int parents_size;
363 
364  int get_edge_degree(Node* v1, Node* v2);
365  void assign_parent(Element* e, int i);
366  void regularize_triangle(Element* e);
367  void regularize_quad(Element* e);
368  void flatten();
369 
372  class HERMES_API CurvedException : public Hermes::Exceptions::Exception
373  {
374  public:
377  CurvedException(int elementId);
378  CurvedException(const CurvedException & e);
379 
380  int getElementId() const;
381  private:
382  int elementId;
383  };
384 
385  friend class MeshHashGrid;
386  friend class MeshReaderH2D;
387  friend class MeshReaderH2DBSON;
388  friend class MeshReaderH2DXML;
389  friend class MeshReaderH1DXML;
390  friend class MeshReaderExodusII;
391  friend class DiscreteProblem < double > ;
392  friend class DiscreteProblem < std::complex<double> > ;
393  template<typename Scalar> friend class DiscreteProblemDGAssembler;
394  friend class WeakForm < double > ;
395  friend class WeakForm < std::complex<double> > ;
396  template<typename Scalar> friend class Adapt;
397  friend class KellyTypeAdapt < double > ;
398  friend class KellyTypeAdapt < std::complex<double> > ;
399  template<typename Scalar> friend class Solution;
400  template<typename Scalar> friend class Filter;
401  template<typename Scalar> friend class MeshFunction;
402  friend class RefMap;
403  friend class Traverse;
404  friend class Transformable;
405  friend class Curved;
406  friend class Views::MeshView;
407  template<typename Scalar> friend class RefinementSelectors::Selector;
408  template<typename Scalar> friend class RefinementSelectors::POnlySelector;
409  template<typename Scalar> friend class RefinementSelectors::HOnlySelector;
410  template<typename Scalar> friend class RefinementSelectors::OptimumSelector;
411  template<typename Scalar> friend class RefinementSelectors::ProjBasedSelector;
412  template<typename Scalar> friend class RefinementSelectors::H1ProjBasedSelector;
413  template<typename Scalar> friend class RefinementSelectors::L2ProjBasedSelector;
414  template<typename Scalar> friend class RefinementSelectors::HcurlProjBasedSelector;
415  friend class PrecalcShapeset;
416  template<typename Scalar> friend class Space;
417  template<typename Scalar> friend class H1Space;
418  template<typename Scalar> friend class HcurlSpace;
419  template<typename Scalar> friend class HdivSpace;
420  template<typename Scalar> friend class L2Space;
421  friend class Views::ScalarView;
422  friend class Views::Orderizer;
423  public:
424  const ElementMarkersConversion &get_element_markers_conversion() const;
425  const BoundaryMarkersConversion &get_boundary_markers_conversion() const;
426  ElementMarkersConversion &get_element_markers_conversion();
427  BoundaryMarkersConversion &get_boundary_markers_conversion();
428 
429  Element* create_quad(int marker, Node* v0, Node* v1, Node* v2, Node* v3, CurvMap* cm, int id = -1);
430  Element* create_triangle(int marker, Node* v0, Node* v1, Node* v2, CurvMap* cm, int id = -1);
431  void refine_element(Element* e, int refinement);
432 
435  std::vector<std::pair<unsigned int, int> > refinements;
436 
442  void refine_quad(Element* e, int refinement, Element** sons_out = nullptr);
443  void refine_triangle_to_triangles(Element* e, Element** sons = nullptr);
444 
446  static double vector_length(double a_1, double a_2);
447 
449  static bool same_line(double p_1, double p_2, double q_1, double q_2, double r_1, double r_2);
450 
452  static bool is_convex(double a_1, double a_2, double b_1, double b_2);
453  static void check_triangle(int i, Node *&v0, Node *&v1, Node *&v2);
454  static void check_quad(int i, Node *&v0, Node *&v1, Node *&v2, Node *&v3);
455 
456  friend class EggShell;
457  };
458 
459  class HERMES_API EggShell
460  {
461  public:
468  static MeshSharedPtr get_egg_shell(MeshSharedPtr mesh, std::string marker, unsigned int levels);
469 
475  static MeshSharedPtr get_egg_shell(MeshSharedPtr mesh, std::vector<std::string> markers, unsigned int levels);
476 
481 
484 
487  static bool egg_shell_verbose;
488 
489  private:
498  static void get_egg_shell_structures(MeshSharedPtr target_mesh, std::vector<std::string> markers, unsigned int levels);
499 
502  static void handle_vertex_on_target_mesh(Element* e, int vertex, MeshSharedPtr target_mesh, std::vector<int> markers, int* neighbor_targets_local);
503 
506  static void fix_hanging_nodes(MeshSharedPtr target_mesh);
507 
509  static void mark_elements_down_used(int eggShell_marker_volume, Element* element);
510 
512  static const std::string eggShellInnerMarker;
513 
516  static void fix_markers(MeshSharedPtr target_mesh, MeshSharedPtr original_mesh);
517  };
518  }
519 }
520 #endif
std::map< std::string, int > conversion_table_inverse
Definition: mesh.h:287
Definition: adapt.h:24
static bool egg_shell_verbose
Definition: mesh.h:487
A projection-based selector for Hcurl space.
Definition: function.h:36
Stores one element of a mesh.
Definition: element.h:107
Caches precalculated shape function values.
Definition: precalc.h:32
Represents a finite element mesh. Typical usage: MeshSharedPtr mesh; Hermes::Hermes2D::MeshReaderH2DX...
Definition: mesh.h:61
Represents a function defined on a mesh.
Definition: mesh_function.h:56
A projection-based selector for H1 space.
Definition: function.h:35
std::map< int, std::string > conversion_table
Definition: mesh.h:283
Stores one node of a mesh.
Definition: element.h:45
A selector that selects H-refinements only.
Definition: function.h:30
Struct for return type of get_user_marker().
Definition: mesh.h:242
A general projection-based selector.
Definition: function.h:33
Represents the weak formulation of a PDE problem.
Definition: global.h:86
static const std::string eggShell0Marker
The mesh returned from get_egg_shell has this marker on the "0" boundary.
Definition: mesh.h:480
A selector that increases order (i.e., it selects P-refinements only).
Definition: function.h:31
Struct for return type of get_internal_marker().
Definition: mesh.h:251
static const std::string eggShell1Marker
The mesh returned from get_egg_shell has this marker on the "1" boundary.
Definition: mesh.h:478
Visualizes a Scalar PDE solution.
Definition: scalar_view.h:37
::xsd::cxx::tree::string< char, simple_type > string
C++ type corresponding to the string XML Schema built-in type.
A selector that chooses an optimal candidates based on a score.
Definition: function.h:32
static const std::string eggShellMarker
Internal marker for eggshell elements.
Definition: mesh.h:483
Represents a finite element space over a domain.
Definition: api2d.h:34
Represents the reference mapping.
Definition: refmap.h:40
A parent of all refinement selectors. Abstract class.
Definition: function.h:29
Class for creating reference mesh.
Definition: mesh.h:209
A projection-based selector for L2 space.
Definition: function.h:34
Represents the solution of a PDE.
Definition: api2d.h:35
Stores and searches node tables.
Definition: hash.h:39