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;
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"); }
336 scalar_type tol=scalar_type(0));
339 void sup_convex(
size_type ic,
bool sup_points =
false);
356 const base_node &pt)
const;
376 base_small_vector mean_normal_of_face_of_convex(
size_type ic,
386 const base_node &pt)
const;
407 scalar_type minimal_convex_radius_estimate()
const;
409 scalar_type maximal_convex_radius_estimate()
const;
411 void translation(
const base_small_vector &);
413 void transformation(
const base_matrix &);
415 void bounding_box(base_node& Pmin, base_node &Pmax)
const;
424 else if (!has_region(
id)) {
425 valid_cvf_sets.add(
id);
432 if (!has_region(
id)) {
433 valid_cvf_sets.add(
id);
441 const dal::bit_vector ®ions_index()
const {
return valid_cvf_sets; }
448 if (valid_cvf_sets[b])
449 { cvf_sets[b].clear(); valid_cvf_sets.sup(b); touch(); }
453 void sup_convex_from_regions(
size_type cv);
456 void optimize_structure(
bool with_renumbering =
true);
458 const std::vector<size_type> &cuthill_mckee_ordering()
const;
464 void write_to_file(
const std::string &name)
const;
468 void write_to_file(std::ostream &ost)
const;
473 void read_from_file(
const std::string &name);
478 void read_from_file(std::istream &ist);
480 void copy_from(
const mesh& m);
486 void touch_from_region(
size_type ) { touch(); }
497 struct green_simplex {
499 std::vector<size_type> sub_simplices;
501 std::vector<size_type> ipt_loc;
505 bool operator <(
const edge &e)
const;
509 typedef std::set<edge> edge_set;
511 struct Bank_info_struct {
512 dal::bit_vector is_green_simplex;
513 std::map<size_type, size_type> num_green_simplex;
514 dal::dynamic_tas<green_simplex> green_simplices;
518 std::unique_ptr<Bank_info_struct> Bank_info;
521 void set_name(
const std::string&);
524 std::vector<size_type> &);
525 bool Bank_is_convex_having_points(
size_type,
526 const std::vector<size_type> &);
527 void Bank_sup_convex_from_green(
size_type);
529 void Bank_build_first_mesh(mesh &,
size_type);
530 void Bank_basic_refine_convex(
size_type);
531 void Bank_refine_normal_convex(
size_type);
534 void Bank_build_green_simplexes(
size_type, std::vector<size_type> &);
539 void Bank_refine(dal::bit_vector);
544 const mesh &dummy_mesh();
559 ITER ipts,
const scalar_type tol)
562 std::vector<size_type> ind(nb);
563 for (
short_type i = 0; i < nb; ++ipts, ++i) ind[i] = add_point(*ipts, tol);
573 {
return add_convex(bgeot::parallelepiped_geotrans(di, 1), ipts); }
577 (dim_type di,
const ITER &ps)
578 {
return add_convex_by_points(bgeot::parallelepiped_geotrans(di, 1), ps); }
582 {
return add_convex(bgeot::prism_geotrans(di, 1), ipts); }
586 (dim_type di,
const ITER &ps)
587 {
return add_convex_by_points(bgeot::prism_geotrans(di, 1), ps); }
596 const base_matrix& pts,
597 pintegration_method pim);
602 const base_matrix& pts);
607 const base_matrix& pts);
611 typedef std::vector<convex_face> convex_face_ct;
612 struct APIDECL convex_face {
615 inline bool operator < (
const convex_face &e)
const {
616 if (cv < e.cv)
return true;
617 if (cv > e.cv)
return false;
618 if (f < e.f)
return true;
else if (f > e.f)
return false;
621 bool is_face()
const {
return f !=
short_type(-1); }
631 convex_face_ct &flist);
669 const mesh_region &mr,
670 const base_small_vector &V,
678 const base_node &pt1,
679 const base_node &pt2);
686 const base_node ¢er,
690 select_convexes_in_box(
const mesh &m,
const mesh_region &mr,
691 const base_node &pt1,
692 const base_node &pt2);
694 inline mesh_region APIDECL
695 select_convexes_in_box(
const mesh &m,
696 const base_node &pt1,
697 const base_node &pt2)
698 {
return select_convexes_in_box(m, m.region(-1), pt1, pt2); }