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);
 
  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)
 
  582   { 
return add_convex(bgeot::prism_geotrans(di, 1), ipts); }
 
  586     (dim_type di, 
const ITER &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); }
 
Inversion of geometric transformations.
 
generic definition of a convex ( bgeot::convex_structure + vertices coordinates )
 
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.
 
base class for static stored objects
 
Deal with interdependencies of objects.
 
structure used to hold a set of convexes and/or convex faces.
 
static mesh_region all_convexes()
provide a default value for the mesh_region parameters of assembly procedures etc.
 
Describe a mesh (collection of convexes (elements) and points).
 
const std::vector< size_type > & cuthill_mckee_ordering() const
Return the list of convex IDs for a Cuthill-McKee ordering.
 
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.
 
ref_convex convex(size_type ic) const
return a bgeot::convex object for the convex number ic.
 
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.
 
void update_from_context() const
this function has to be defined and should update the object when the context is modified.
 
size_type search_point(const base_node &pt, const scalar_type radius=0) const
Search a point given its coordinates.
 
void sup_region(size_type b)
Remove the region of index b.
 
gmm::uint64_type convex_version_number(size_type ic) const
return the version number of the convex ic.
 
const mesh_region region(size_type id) const
Return the region of index 'id'.
 
size_type add_convex(bgeot::pgeometric_trans pgt, ITER ipts)
Add a convex to the mesh.
 
void sup_point(size_type i)
Delete the point of index i from the mesh if it is not linked to a convex.
 
bool has_region(size_type s) const
Return true if the region of index 's' exists in the mesh.
 
void swap_points(size_type i, size_type j)
Swap the indexes of points of index i and j in the whole structure.
 
size_type add_simplex(dim_type di, ITER ipts)
Add a simplex to the mesh, given its dimension and point numbers.
 
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.
 
indexed array reference (given a container X, and a set of indexes I, this class provides a pseudo-co...
 
Deal with interdependencies of objects (getfem::context_dependencies).
 
region objects (set of convexes and/or convex faces)
 
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...
 
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...
 
size_type add_simplex_by_points(dim_type dim, ITER ipts)
Add a simplex to the mesh, given its dimension and point coordinates.
 
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_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.
 
mesh_region APIDECL select_faces_in_ball(const mesh &m, const mesh_region &mr, const base_node ¢er, scalar_type radius)
Select in the region mr the faces of the mesh m lying entirely in the ball delimated by pt1 and radiu...
 
scalar_type APIDECL convex_area_estimate(bgeot::pgeometric_trans pgt, const base_matrix &pts, pintegration_method pim)
rough estimate of the convex area.
 
size_type add_prism(dim_type di, const ITER &ipts)
Add a prism to the mesh.
 
size_type add_parallelepiped(dim_type di, const ITER &ipts)
Add a parallelepiped to the mesh.
 
mesh_region APIDECL inner_faces_of_mesh(const mesh &m, const mesh_region &mr=mesh_region::all_convexes())
Select all the faces sharing at least two element of the given mesh region.
 
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.
 
size_type add_parallelepiped_by_points(dim_type di, const ITER &ps)
Add a parallelepiped to the mesh.
 
size_type add_prism_by_points(dim_type di, const ITER &ps)
Add a prism to the mesh.
 
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.
 
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.
 
mesh_region APIDECL all_faces_of_mesh(const mesh &m, const mesh_region &mr=mesh_region::all_convexes())
Select all the faces of the given mesh region.
 
gmm::uint16_type short_type
used as the common short type integer in the library
 
size_t size_type
used as the common size type in the library
 
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
 
GEneric Tool for Finite Element Methods.