38 #ifndef GETFEM_MESH_H__ 39 #define GETFEM_MESH_H__ 53 gmm::uint64_type APIDECL act_counter();
55 class integration_method;
56 typedef std::shared_ptr<const integration_method> pintegration_method;
95 class APIDECL
mesh :
public bgeot::basic_mesh,
98 public std::enable_shared_from_this<mesh>
103 typedef bgeot::mesh_structure::ind_cv_ct ind_cv_ct;
104 typedef bgeot::mesh_structure::ind_set ind_set;
117 mutable std::map<size_type, mesh_region> cvf_sets;
118 mutable dal::bit_vector valid_cvf_sets;
119 void handle_region_refinement(
size_type,
const std::vector<size_type> &,
122 mutable bool cuthill_mckee_uptodate;
124 mutable std::vector<size_type> cmk_order;
127 #if GETFEM_PARA_LEVEL > 1 128 mutable bool modified;
130 mutable std::map<size_type, mesh_region> mpi_sub_region;
132 mutable dal::bit_vector valid_sub_regions;
135 modified =
true; cuthill_mckee_uptodate =
false;
136 context_dependencies::touch();
138 void compute_mpi_region()
const ;
139 void compute_mpi_sub_region(
size_type)
const;
144 {
if (modified) compute_mpi_region();
return mpi_region; }
146 if (modified) compute_mpi_region();
147 if (n ==
size_type(-1))
return mpi_region;
148 if (!(valid_sub_regions.is_in(n))) compute_mpi_sub_region(n);
149 return mpi_sub_region[n];
151 void intersect_with_mpi_region(
mesh_region &rg)
const;
154 { cuthill_mckee_uptodate =
false; context_dependencies::touch(); }
159 if (n ==
size_type(-1))
return get_mpi_region();
162 void intersect_with_mpi_region(
mesh_region &)
const {}
166 explicit mesh(
const std::string name =
"");
167 mesh(
const bgeot::basic_mesh &m,
const std::string name =
"");
171 inline std::string get_name()
const {
return name_;}
174 using basic_mesh::dim;
176 using basic_mesh::points;
177 PT_TAB &points() {
return pts; }
180 using basic_mesh::points_of_convex;
184 ind_pt_face_ct rct = ind_points_of_face_of_convex(ic, f);
185 return ref_mesh_face_pt_ct(pts.begin(), rct.begin(), rct.end());
190 {
return ref_convex(structure_of_convex(ic), points_of_convex(ic)); }
192 using basic_mesh::add_point;
203 {
if (i != j) { pts.swap_points(i,j); mesh_structure::swap_points(i,j); } }
210 {
return pts.search_node(pt,radius); }
212 using basic_mesh::trans_of_convex;
216 {
return cvs_v_num[ic]; }
230 gtab[i] = pgt; trans_exists[i] =
true;
231 if (!present) { cvs_v_num[i] = act_counter(); touch(); }
247 const scalar_type tol=scalar_type(0));
253 {
return add_convex(bgeot::simplex_geotrans(di, 1), ipts); }
258 size_type add_simplex_by_points(dim_type dim, ITER ipts);
263 const base_node &pt2) {
266 return add_segment(ipt1, ipt2);
271 size_type add_triangle_by_points(
const base_node &p1,
273 const base_node &p3);
278 size_type add_tetrahedron_by_points(
const base_node &p1,
281 const base_node &p4);
291 size_type add_parallelepiped(dim_type di,
const ITER &ipts);
298 size_type add_parallelepiped_by_points(dim_type di,
const ITER &ps);
315 size_type add_prism(dim_type di,
const ITER &ipts);
323 size_type add_prism_by_points(dim_type di,
const ITER &ps);
326 { GMM_ASSERT1(
false,
"Sorry, to be done"); }
329 { GMM_ASSERT1(
false,
"Sorry, to be done"); }
332 void sup_convex(
size_type ic,
bool sup_points =
false);
349 const base_node &pt)
const;
369 base_small_vector mean_normal_of_face_of_convex(
size_type ic,
379 const base_node &pt)
const;
400 scalar_type minimal_convex_radius_estimate()
const;
402 scalar_type maximal_convex_radius_estimate()
const;
404 void translation(
const base_small_vector &);
406 void transformation(
const base_matrix &);
408 void bounding_box(base_node& Pmin, base_node &Pmax)
const;
417 else if (!has_region(
id)) {
418 valid_cvf_sets.add(
id);
419 cvf_sets[id] =
mesh_region(const_cast<mesh&>(*
this),
id);
425 if (!has_region(
id)) {
426 valid_cvf_sets.add(
id);
434 const dal::bit_vector ®ions_index()
const {
return valid_cvf_sets; }
441 if (valid_cvf_sets[b])
442 { cvf_sets[b].clear(); valid_cvf_sets.sup(b); touch(); }
446 void sup_convex_from_regions(
size_type cv);
449 void optimize_structure(
bool with_renumbering =
true);
451 const std::vector<size_type> &cuthill_mckee_ordering()
const;
457 void write_to_file(
const std::string &name)
const;
461 void write_to_file(std::ostream &ost)
const;
466 void read_from_file(
const std::string &name);
471 void read_from_file(std::istream &ist);
473 void copy_from(
const mesh& m);
479 void touch_from_region(
size_type ) { touch(); }
490 struct green_simplex {
492 std::vector<size_type> sub_simplices;
494 std::vector<size_type> ipt_loc;
498 bool operator <(
const edge &e)
const;
502 typedef std::set<edge> edge_set;
504 struct Bank_info_struct {
505 dal::bit_vector is_green_simplex;
506 std::map<size_type, size_type> num_green_simplex;
507 dal::dynamic_tas<green_simplex> green_simplices;
511 std::unique_ptr<Bank_info_struct> Bank_info;
514 void set_name(
const std::string&);
517 std::vector<size_type> &);
518 bool Bank_is_convex_having_points(
size_type,
519 const std::vector<size_type> &);
520 void Bank_sup_convex_from_green(
size_type);
523 void Bank_basic_refine_convex(
size_type);
524 void Bank_refine_normal_convex(
size_type);
527 void Bank_build_green_simplexes(
size_type, std::vector<size_type> &);
532 void Bank_refine(dal::bit_vector);
537 const mesh &dummy_mesh();
552 ITER ipts,
const scalar_type tol)
555 std::vector<size_type> ind(nb);
556 for (
short_type i = 0; i < nb; ++ipts, ++i) ind[i] = add_point(*ipts, tol);
557 return add_convex(pgt, ind.begin());
562 {
return add_convex_by_points(bgeot::simplex_geotrans(di, 1), ipts); }
566 {
return add_convex(bgeot::parallelepiped_geotrans(di, 1), ipts); }
570 (dim_type di,
const ITER &ps)
571 {
return add_convex_by_points(bgeot::parallelepiped_geotrans(di, 1), ps); }
575 {
return add_convex(bgeot::prism_geotrans(di, 1), ipts); }
579 (dim_type di,
const ITER &ps)
580 {
return add_convex_by_points(bgeot::prism_geotrans(di, 1), ps); }
589 const base_matrix& pts,
590 pintegration_method pim);
595 const base_matrix& pts);
600 const base_matrix& pts);
604 typedef std::vector<convex_face> convex_face_ct;
605 struct APIDECL convex_face {
608 inline bool operator < (
const convex_face &e)
const {
609 if (cv < e.cv)
return true;
610 if (cv > e.cv)
return false;
611 if (f < e.f)
return true;
else if (f > e.f)
return false;
614 bool is_face()
const {
return f !=
short_type(-1); }
624 convex_face_ct& flist);
653 const base_small_vector &V,
661 const base_node &pt1,
662 const base_node &pt2);
666 const base_node &pt1,
667 const base_node &pt2);
670 select_convexes_in_box(
const mesh &m,
671 const base_node &pt1,
672 const base_node &pt2)
673 {
return select_convexes_in_box(m, m.
region(-1), pt1, pt2); }
size_type add_parallelepiped(dim_type di, const ITER &ipts)
Add a parallelepiped to the mesh.
structure used to hold a set of convexes and/or convex faces.
mesh_region APIDECL select_faces_in_box(const mesh &m, const mesh_region &mr, const base_node &pt1, const base_node &pt2)
Select in the region mr the faces of the mesh m lying entirely in the box delimated by pt1 and pt2...
region objects (set of convexes and/or convex faces)
indexed array reference (given a container X, and a set of indexes I, this class provides a pseudo-co...
base class for static stored objects
static mesh_region all_convexes()
provide a default value for the mesh_region parameters of assembly procedures etc.
scalar_type APIDECL convex_area_estimate(bgeot::pgeometric_trans pgt, const base_matrix &pts, pintegration_method pim)
rough estimate of the convex area.
bool has_region(size_type s) const
Return true if the region of index 's' exists in the mesh.
ref_mesh_face_pt_ct points_of_face_of_convex(size_type ic, short_type f) const
Return a (pseudo)container of points of face of a given convex.
size_type add_simplex(dim_type di, ITER ipts)
Add a simplex to the mesh, given its dimension and point numbers.
Describe a mesh (collection of convexes (elements) and points).
size_type add_segment_by_points(const base_node &pt1, const base_node &pt2)
Add a segment to the mesh, given the coordinates of its vertices.
ref_convex convex(size_type ic) const
return a bgeot::convex object for the convex number ic.
size_type add_convex(bgeot::pgeometric_trans pgt, ITER ipts)
Add a convex to the mesh.
generic definition of a convex ( bgeot::convex_structure + vertices coordinates ) ...
void sup_point(size_type i)
Delete the point of index i from the mesh if it is not linked to a convex.
size_t size_type
used as the common size type in the library
Inversion of geometric transformations.
mesh_region APIDECL inner_faces_of_mesh(const mesh &m, mesh_region mr=mesh_region::all_convexes())
Select all the faces sharing at least two element of the given mesh region.
void update_from_context() const
this function has to be defined and should update the object when the context is modified.
size_type add_prism_by_points(dim_type di, const ITER &ps)
Add a prism to the mesh.
size_type search_point(const base_node &pt, const scalar_type radius=0) const
Search a point given its coordinates.
size_type add_prism(dim_type di, const ITER &ipts)
Add a prism to the mesh.
void APIDECL extrude(const mesh &in, mesh &out, size_type nb_layers, short_type degree=short_type(1))
build a N+1 dimensions mesh from a N-dimensions mesh by extrusion.
Deal with interdependencies of objects.
GEneric Tool for Finite Element Methods.
scalar_type APIDECL convex_radius_estimate(bgeot::pgeometric_trans pgt, const base_matrix &pts)
rough estimate of the radius of the convex using the largest eigenvalue of the jacobian of the geomet...
gmm::uint16_type short_type
used as the common short type integer in the library
scalar_type APIDECL convex_quality_estimate(bgeot::pgeometric_trans pgt, const base_matrix &pts)
rough estimate of the maximum value of the condition number of the jacobian of the geometric transfor...
Deal with interdependencies of objects (getfem::context_dependencies).
void sup_region(size_type b)
Remove the region of index b.
size_type add_convex(pconvex_structure cs, ITER ipts, bool *present=0)
Insert a new convex in the mesh_structure.
Store a set of points, identifying points that are nearer than a certain very small distance...
mesh_region APIDECL select_faces_of_normal(const mesh &m, const mesh_region &mr, const base_small_vector &V, scalar_type angle)
Select in the region mr the faces of the mesh m with their unit outward vector having a maximal angle...
size_type add_simplex_by_points(dim_type dim, ITER ipts)
Add a simplex to the mesh, given its dimension and point coordinates.
const dal::bit_vector & points_index() const
Return the points index.
size_type nb_points() const
Give the number of geometrical nodes in the mesh.
const mesh_region region(size_type id) const
Return the region of index 'id'.
size_type add_convex_by_points(bgeot::pgeometric_trans pgt, ITER ipts, const scalar_type tol=scalar_type(0))
Add a convex to the mesh, given a geometric transformation and a list of point coordinates.
gmm::uint64_type convex_version_number(size_type ic) const
return the version number of the convex ic.
void swap_points(size_type i, size_type j)
Swap the indexes of points of index i and j in the whole structure.
const dal::bit_vector & convex_index() const
Return the list of valid convex IDs.
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
size_type add_parallelepiped_by_points(dim_type di, const ITER &ps)
Add a parallelepiped to the mesh.
void APIDECL outer_faces_of_mesh(const mesh &m, const dal::bit_vector &cvlst, convex_face_ct &flist)
returns a list of "exterior" faces of a mesh (i.e.