20 #include "mesh_util.h"
22 #include "../mixins2d.h"
31 template<
typename Scalar>
class Space;
35 extern HERMES_API
unsigned g_mesh_seq;
37 namespace RefinementSelectors
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;
62 public Hermes::Mixins::StateQueryable
70 void init(
int size = H2D_DEFAULT_HASH_SIZE);
73 virtual bool isOkay()
const;
74 inline std::string getClassName()
const {
return "Mesh"; }
77 bool rescale(
double x_ref,
double y_ref);
80 void copy(MeshSharedPtr mesh);
83 void copy_base(MeshSharedPtr mesh);
86 void copy_converted(MeshSharedPtr mesh);
89 void create(
int nv, double2* verts,
int nt, int3* tris,
std::string* tri_markers,
92 #pragma region MeshHashGrid
93 Element* element_on_physical_coordinates(
double x,
double y);
103 #pragma region MarkerArea
104 double get_marker_area(
int marker);
106 std::map<int, MarkerArea*> marker_areas;
109 #pragma region getters
110 Element* get_element(
int id)
const;
114 int get_num_elements()
const;
117 int get_num_base_elements()
const;
120 int get_num_used_base_elements()
const;
123 int get_num_active_elements()
const;
126 int get_max_element_id()
const;
129 int get_num_vertex_nodes()
const;
132 int get_num_edge_nodes()
const;
139 void get_bounding_box(
double& bottom_left_x,
double& bottom_left_y,
double& top_right_x,
double& top_right_y);
142 Element* get_element_fast(
int id)
const;
145 unsigned get_seq()
const;
148 #pragma region refinements
149 void refine_element_id(
int id,
int refinement = 0);
159 void refine_all_elements(
int refinement = 0,
bool mark_as_initial =
false);
168 void refine_by_criterion(
int(*criterion)(
Element* e),
int depth = 1,
bool mark_as_initial =
false);
172 void refine_towards_vertex(
int vertex_id,
int depth = 1,
bool mark_as_initial =
false);
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);
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);
194 int* regularize(
int n);
198 void unrefine_element_id(
int id);
205 void unrefine_all_elements(
bool keep_initial_refinements =
true);
222 virtual MeshSharedPtr create_ref_mesh();
226 MeshSharedPtr coarse_mesh;
261 StringValid get_user_marker(
int internal_marker)
const;
266 enum MarkersConversionType {
267 HERMES_ELEMENT_MARKERS_CONVERSION = 0,
268 HERMES_BOUNDARY_MARKERS_CONVERSION = 1
271 virtual MarkersConversionType get_type()
const = 0;
289 friend class Space < double > ;
290 friend class Space <
std::complex<double> > ;
298 void set_seq(unsigned seq);
300 int nbase, ntopvert, ninitial, nactive;
306 virtual MarkersConversionType get_type()
const;
313 virtual MarkersConversionType get_type()
const;
318 Array<Element> elements;
323 void initial_single_check();
329 void convert_quads_to_triangles();
332 void convert_to_base();
334 void refine_element_to_quads_id(
int id);
336 void refine_triangle_to_quads(
Element* e,
Element** elems_out =
nullptr);
338 void refine_element_to_triangles_id(
int id);
340 void refine_quad_to_triangles(
Element* e);
346 void refine_quad_to_quads(
Element* e,
int refinement = 0);
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);
353 double bottom_left_x, bottom_left_y, top_right_x, top_right_y;
355 bool bounding_box_calculated;
357 void calc_bounding_box();
359 void unrefine_element_internal(
Element* e);
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);
372 class HERMES_API CurvedException :
public Hermes::Exceptions::Exception
377 CurvedException(
int elementId);
378 CurvedException(
const CurvedException & e);
380 int getElementId()
const;
387 friend class MeshReaderH2DBSON;
395 friend class
WeakForm < std::complex<double> > ;
396 template<typename Scalar> friend class
Adapt;
399 template<typename Scalar> friend class
Solution;
400 template<typename Scalar> friend class
Filter;
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;
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();
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);
435 std::vector<std::pair<unsigned int, int> > refinements;
442 void refine_quad(Element* e, int refinement, Element** sons_out = nullptr);
443 void refine_triangle_to_triangles(Element* e, Element** sons = nullptr);
446 static double vector_length(double a_1, double a_2);
449 static bool same_line(double p_1, double p_2, double q_1, double q_2, double r_1, double r_2);
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);
468 static MeshSharedPtr get_egg_shell(MeshSharedPtr mesh,
std::string marker,
unsigned int levels);
475 static MeshSharedPtr get_egg_shell(MeshSharedPtr mesh, std::vector<std::string> markers,
unsigned int levels);
498 static void get_egg_shell_structures(MeshSharedPtr target_mesh, std::vector<std::string> markers,
unsigned int levels);
502 static void handle_vertex_on_target_mesh(
Element* e,
int vertex, MeshSharedPtr target_mesh, std::vector<int> markers,
int* neighbor_targets_local);
506 static void fix_hanging_nodes(MeshSharedPtr target_mesh);
509 static void mark_elements_down_used(
int eggShell_marker_volume,
Element* element);
516 static void fix_markers(MeshSharedPtr target_mesh, MeshSharedPtr original_mesh);
std::map< std::string, int > conversion_table_inverse
static bool egg_shell_verbose
A projection-based selector for Hcurl space.
Stores one element of a mesh.
Caches precalculated shape function values.
Represents a finite element mesh. Typical usage: MeshSharedPtr mesh; Hermes::Hermes2D::MeshReaderH2DX...
Represents a function defined on a mesh.
A projection-based selector for H1 space.
std::map< int, std::string > conversion_table
Stores one node of a mesh.
A selector that selects H-refinements only.
Struct for return type of get_user_marker().
A general projection-based selector.
static const std::string eggShell0Marker
The mesh returned from get_egg_shell has this marker on the "0" boundary.
A selector that increases order (i.e., it selects P-refinements only).
Struct for return type of get_internal_marker().
static const std::string eggShell1Marker
The mesh returned from get_egg_shell has this marker on the "1" boundary.
Visualizes a Scalar PDE solution.
::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.
static const std::string eggShellMarker
Internal marker for eggshell elements.
Represents a finite element space over a domain.
Represents the reference mapping.
A parent of all refinement selectors. Abstract class.
Class for creating reference mesh.
A projection-based selector for L2 space.
Represents the solution of a PDE.
Stores and searches node tables.