GetFEM  5.5
getfem_mesh.h
Go to the documentation of this file.
1 /* -*- c++ -*- (enables emacs c++ mode) */
2 /*===========================================================================
3 
4  Copyright (C) 1999-2026 Yves Renard
5 
6  This file is a part of GetFEM
7 
8  GetFEM is free software; you can redistribute it and/or modify it
9  under the terms of the GNU Lesser General Public License as published
10  by the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version along with the GCC Runtime Library
12  Exception either version 3.1 or (at your option) any later version.
13  This program is distributed in the hope that it will be useful, but
14  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
15  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
16  License and GCC Runtime Library Exception for more details.
17  You should have received a copy of the GNU Lesser General Public License
18  along with this program. If not, see https://www.gnu.org/licenses/.
19 
20  As a special exception, you may use this file as it is a part of a free
21  software library without restriction. Specifically, if other files
22  instantiate templates or use macros or inline functions from this file,
23  or you compile this file and link it with other files to produce an
24  executable, this file does not by itself cause the resulting executable
25  to be covered by the GNU Lesser General Public License. This exception
26  does not however invalidate any other reasons why the executable file
27  might be covered by the GNU Lesser General Public License.
28 
29 ===========================================================================*/
30 
31 /**@file getfem_mesh.h
32  @author Yves Renard <[email protected]>
33  @date November 05, 1999.
34  @brief Define a getfem::getfem_mesh object.
35 */
36 
37 #ifndef GETFEM_MESH_H__
38 #define GETFEM_MESH_H__
39 
40 #include <bitset>
41 #include "bgeot_ftool.h"
42 #include "bgeot_mesh.h"
43 #include "bgeot_geotrans_inv.h"
44 #include "getfem_context.h"
45 #include "getfem_mesh_region.h"
46 #include <memory>
47 
48 
49 namespace getfem {
50 
51  /* Version counter for convexes. */
52  gmm::uint64_type APIDECL act_counter();
53 
54  class integration_method;
55  typedef std::shared_ptr<const integration_method> pintegration_method;
56 
57  /**@addtogroup mesh*/
58  /**@{*/
59 
60  /** Describe a mesh (collection of convexes (elements) and points).
61  Note that mesh object have no copy constructor, use
62  mesh::copy_from instead. This class inherits from
63  bgeot::mesh<base_node> and bgeot::mesh_structure. Compared to
64  bgeot::mesh, It provides some additional methods for convenience
65  (add_simplex_by_points etc...), and it is able to send
66  "messages" to objects that depend on it (for example
67  getfem::mesh_fem) when something in the mesh has been
68  modified. For example, when a convex is removed from the mesh,
69  the mesh_fem that uses this mesh automagically removes any finite
70  element on the convex, and renumbers the dof.
71 
72  Points and convex numbering:
73 
74  The numbering of points and convexes is not dynamical, which
75  means that when a point or a convex has been removed from the mesh,
76  there might be holes in the numbering. To loop over the set of
77  valid points in the mesh, one should use
78 
79  @code
80  for (dal::bv_visitor ip(mesh.points_index()); !ip.finished(); ++ip) {
81  ...
82  }
83  @endcode
84 
85  instead of
86 
87  @code
88  for (size_type ip = 0; ip < mesh.nb_points(); ++ip) {
89  }
90  @endcode
91  (same thing for the convexes, always use convex_index()).
92  */
93 
94  class APIDECL mesh : public bgeot::basic_mesh,
95  public context_dependencies,
96  virtual public dal::static_stored_object,
97  public std::enable_shared_from_this<mesh>
98  {
99  public :
100 
102  typedef bgeot::mesh_structure::ind_cv_ct ind_cv_ct;
103  typedef bgeot::mesh_structure::ind_set ind_set;
109 
110 
111  protected :
112  /*
113  * When a new field is added, do NOT forget to add it in copy_from method!
114  */
115 
116  mutable std::map<size_type, mesh_region> cvf_sets;
117  mutable dal::bit_vector valid_cvf_sets;
118  void handle_region_refinement(size_type, const std::vector<size_type> &,
119  bool);
120 
121  mutable bool cuthill_mckee_uptodate;
123  mutable std::vector<size_type> cmk_order; // cuthill-mckee
124  void init();
125 
126 #if GETFEM_PARA_LEVEL > 1
127  mutable bool modified;
128  mutable mesh_region mpi_region;
129  mutable std::map<size_type, mesh_region> mpi_sub_region;
130  // mutable dal::dynamic_array<mesh_region> mpi_sub_region;
131  mutable dal::bit_vector valid_sub_regions;
132 
133  void touch() {
134  modified = true; cuthill_mckee_uptodate = false;
135  context_dependencies::touch();
136  }
137  void compute_mpi_region() const ;
138  void compute_mpi_sub_region(size_type) const;
139 
140  public :
141 
142  const mesh_region& get_mpi_region() const
143  { if (modified) compute_mpi_region(); return mpi_region; }
144  const mesh_region& get_mpi_sub_region(size_type n) const {
145  if (modified) compute_mpi_region();
146  if (n == size_type(-1)) return mpi_region;
147  if (!(valid_sub_regions.is_in(n))) compute_mpi_sub_region(n);
148  return mpi_sub_region[n];
149  }
150  void intersect_with_mpi_region(mesh_region &rg) const;
151 #else
152  void touch()
153  { cuthill_mckee_uptodate = false; context_dependencies::touch(); }
154  public :
155  const mesh_region get_mpi_region() const
156  { return mesh_region::all_convexes(); }
157  const mesh_region get_mpi_sub_region(size_type n) const {
158  if (n == size_type(-1)) return get_mpi_region();
159  return cvf_sets[n];
160  }
161  void intersect_with_mpi_region(mesh_region &) const {}
162 #endif
163 
164  /// Constructor.
165  explicit mesh(const std::string name = "");
166  mesh(const bgeot::basic_mesh &m, const std::string name = "");
167  mesh(const mesh &m);
168  mesh &operator=(const mesh &m);
169 
170  inline std::string get_name() const {return name_;}
171  void update_from_context() const {}
172  /// Mesh dimension.
173  using basic_mesh::dim;
174  /// Return the array of PT.
175  using basic_mesh::points;
176  PT_TAB &points() { return pts; } // non-const version
177 
178  /// Return a (pseudo)container of the points of a given convex
179  using basic_mesh::points_of_convex;
180  /// Return a (pseudo)container of points of face of a given convex
182  short_type f) const {
183  ind_pt_face_ct rct = ind_points_of_face_of_convex(ic, f);
184  return ref_mesh_face_pt_ct(pts.begin(), rct.begin(), rct.end());
185  }
186 
187  /// return a bgeot::convex object for the convex number ic.
189  { return ref_convex(structure_of_convex(ic), points_of_convex(ic)); }
190 
191  using basic_mesh::add_point;
192  /// Give the number of geometrical nodes in the mesh.
193  size_type nb_points() const { return pts.card(); }
194  /// Return the points index
195  const dal::bit_vector &points_index() const { return pts.index(); }
196  /** Delete the point of index i from the mesh if it is not linked to a
197  convex.
198  */
199  void sup_point(size_type i) { if (!is_point_valid(i)) pts.sup_node(i); }
200  /// Swap the indexes of points of index i and j in the whole structure.
202  { if (i != j) { pts.swap_points(i,j); mesh_structure::swap_points(i,j); } }
203  /** Search a point given its coordinates.
204  @param pt the point that is searched.
205  @return the point index if pt was found in (or approximatively in)
206  the mesh nodes, size_type(-1) if not found.
207  */
208  size_type search_point(const base_node &pt, const scalar_type radius=0) const
209  { return pts.search_node(pt,radius); }
210 
211  using basic_mesh::trans_of_convex;
212 
213  /// return the version number of the convex ic.
214  gmm::uint64_type convex_version_number(size_type ic) const
215  { return cvs_v_num[ic]; }
216 
217  /** Add a convex to the mesh.
218  This methods assume that the convex nodes have already been
219  added to the mesh.
220  @param pgt the geometric transformation of the convex.
221  @param ipts an iterator to a set of point index.
222  @return the number of the new convex.
223  */
224  template<class ITER>
226  bool present;
227  size_type i = bgeot::mesh_structure::add_convex(pgt->structure(),
228  ipts, &present);
229  gtab[i] = pgt; trans_exists[i] = true;
230  if (!present) { cvs_v_num[i] = act_counter(); touch(); }
231  return i;
232  }
233 
234  /** Add a convex to the mesh, given a geometric transformation and a
235  list of point coordinates.
236 
237  As a side-effect, the points are also added to the mesh (if
238  they were not already in the mesh).
239 
240  @param pgt the geometric transformation of the convex.
241  @param ipts an iterator on a set of getfem::base_node.
242  @return the number of the new convex.
243  */
244  template<class ITER>
245  size_type add_convex_by_points(bgeot::pgeometric_trans pgt, ITER ipts,
246  const scalar_type tol=scalar_type(0));
247 
248  /** Add a simplex to the mesh, given its dimension and point numbers.
249  */
250  template<class ITER>
251  size_type add_simplex(dim_type di, ITER ipts)
252  { return add_convex(bgeot::simplex_geotrans(di, 1), ipts); }
253  /** Add a simplex to the mesh, given its dimension and point coordinates.
254  @see add_convex_by_points.
255  */
256  template<class ITER>
257  size_type add_simplex_by_points(dim_type dim, ITER ipts);
258  /** Add a segment to the mesh, given the point id of its vertices. */
259  size_type add_segment(size_type a, size_type b);
260  /** Add a segment to the mesh, given the coordinates of its vertices. */
261  size_type add_segment_by_points(const base_node &pt1,
262  const base_node &pt2) {
263  size_type ipt1 = add_point(pt1);
264  size_type ipt2 = add_point(pt2);
265  return add_segment(ipt1, ipt2);
266  }
267  /** Add a triangle to the mesh, given the point id of its vertices. */
268  size_type add_triangle(size_type a,size_type b, size_type c);
269  /** Add a triangle to the mesh, given the coordinates of its vertices. */
270  size_type add_triangle_by_points(const base_node &p1,
271  const base_node &p2,
272  const base_node &p3);
273  /** Add a tetrahedron to the mesh, given the point id of its vertices. */
274  size_type add_tetrahedron(size_type a,
275  size_type b, size_type c, size_type d);
276  /** Add a tetrahedron to the mesh, given the coordinates of its vertices. */
277  size_type add_tetrahedron_by_points(const base_node &p1,
278  const base_node &p2,
279  const base_node &p3,
280  const base_node &p4);
281  /** Add a pyramid to the mesh, given the point id of its vertices. */
282  size_type add_pyramid(size_type a,
284  /** Add a parallelepiped to the mesh.
285  @param di dimension of the parallelepiped
286  @param ipts iterator on the list of point id.
287  @return the number of the new convex.
288  */
289  template<class ITER>
290  size_type add_parallelepiped(dim_type di, const ITER &ipts);
291  /** Add a parallelepiped to the mesh.
292  @param di dimension of the parallelepiped
293  @param ps iterator on the list of point coordinates.
294  @return the number of the new convex.
295  */
296  template<class ITER>
297  size_type add_parallelepiped_by_points(dim_type di, const ITER &ps);
298  /* Add a parallelepiped of dimension dim to the
299  * mesh. org is the point of type base_node representing
300  * the origine and "it" is an iterator on a list of
301  * vectors of type base_vector.
302  * Return the index of the convex in the mesh.
303 
304  template<class ITER>
305  size_type add_parallelepiped_by_vectors(dim_type di,
306  const base_node &org, const ITER &vects);
307  */
308  /** Add a prism to the mesh.
309  @param di dimension of the prism
310  @param ipts iterator on the list of point id.
311  @return the number of the new convex.
312  */
313  template<class ITER>
314  size_type add_prism(dim_type di, const ITER &ipts);
315 
316  /** Add a prism to the mesh.
317  @param di dimension of the prism
318  @param ps iterator on the list of point coordinates.
319  @return the number of the new convex.
320  */
321  template<class ITER>
322  size_type add_prism_by_points(dim_type di, const ITER &ps);
323 
324  size_type add_face_of_convex(size_type, short_type)
325  { GMM_ASSERT1(false, "Sorry, to be done"); }
326 
327  void add_faces_of_convex(size_type)
328  { GMM_ASSERT1(false, "Sorry, to be done"); }
329 
330 
331  /** Merge all convexes from a another mesh, possibly restricted to a
332  mesh region.
333  */
334  void merge_convexes_from_mesh(const mesh &m, size_type rg=size_type(-1),
335  scalar_type tol=scalar_type(0));
336 
337  /// Delete the convex of index ic from the mesh.
338  void sup_convex(size_type ic, bool sup_points = false);
339  /** Swap the indexes of the convex of indexes i and j
340  * in the whole structure.
341  */
342  void swap_convex(size_type i, size_type j);
343 
344  /** Return the normal of the given convex face, evaluated at the point pt.
345  @param ic the convex number.
346  @param f the face number.
347  @param pt the point at which the normal is taken, in the
348  reference convex. This point should of course be on the
349  correspounding face of the reference convex, except if the
350  geometric transformation is linear: in that case, the normal
351  constant.
352  @return the face normal.
353  */
354  base_small_vector normal_of_face_of_convex(size_type ic, short_type f,
355  const base_node &pt) const;
356  /* Return the normal of the given convex face.
357 
358  @param ic the convex number.
359  @param f the face number.
360  @param n local point index in the reference convex, at which the normal will be evaluated -- should be faster than the function above since it uses a bgeot::geotrans_precomp .
361 
362  @return the face normal.
363  */
364  base_small_vector normal_of_face_of_convex(size_type ic, short_type f,
365  size_type n=0) const;
366 
367 
368  /* Return the outward unit normal vector of the given element
369  face computed on the mean of the geometrical nodes of the face.
370 
371  @param ic the convex number.
372  @param f the face number.
373  @return the face normal.
374  */
375  base_small_vector mean_normal_of_face_of_convex(size_type ic,
376  short_type f) const;
377  /** Return a local basis for the convex face.
378  @param ic the convex number.
379  @param f the face number.
380  @param pt the point coordinates.
381  @return a matrix with the normal vector in its first column, and the
382  tangent vectors in the other columns.
383  */
384  base_matrix local_basis_of_face_of_convex(size_type ic, short_type f,
385  const base_node &pt) const;
386  /** Return a local basis for the convex face.
387  @param ic the convex number.
388  @param f the face number.
389  @param n local point index in the reference convex, at which the basis will be evaluated -- should be faster than the function above since it uses a bgeot::geotrans_precomp .
390  @return a matrix with the normal vector in its first column, and the
391  tangent vectors in the other columns.
392  */
393  base_matrix local_basis_of_face_of_convex(size_type ic, short_type f,
394  size_type n) const;
395  /** Return an estimate of the convex quality (0 <= Q <= 1). */
396  scalar_type convex_quality_estimate(size_type ic) const;
397  /** Return an estimate of the convex area.
398  @param ic the convex number
399  @param degree the degree of the approximation integration
400  method used to compute the area.
401  */
402  scalar_type convex_area_estimate(size_type ic, size_type degree=2) const;
403  /** Return an estimate of the convex largest dimension. @see getfem::convex_quality_estimate */
404  virtual scalar_type convex_radius_estimate(size_type ic) const;
405  /** Return an estimate of the convex smallest dimension. @see getfem::convex_radius_estimate */
406  scalar_type minimal_convex_radius_estimate() const;
407  /** Return an estimate of the convex largest dimension. @see getfem::convex_radius_estimate */
408  scalar_type maximal_convex_radius_estimate() const;
409  /** Apply the given translation to each mesh node. */
410  void translation(const base_small_vector &);
411  /** apply the given matrix transformation to each mesh node. */
412  void transformation(const base_matrix &);
413  /** Return the bounding box [Pmin - Pmax] of the mesh. */
414  void bounding_box(base_node& Pmin, base_node &Pmax) const;
415  /** Return the region of index 'id'. Regions stored in mesh are
416  updated when the mesh is modified (i.e. when a convex is
417  removed from the mesh, it is also removed from all regions of
418  the mesh.
419  */
420  const mesh_region region(size_type id) const {
421  if (id == mesh_region::all_convexes().id())
422  return mesh_region::all_convexes();
423  else if (!has_region(id)) {
424  valid_cvf_sets.add(id);
425  cvf_sets[id] = mesh_region(const_cast<mesh&>(*this),id);
426  }
427  return cvf_sets[id];
428  }
429  /* Return a reference such that operator= works as expected */
430  mesh_region &region(size_type id) {
431  if (!has_region(id)) {
432  valid_cvf_sets.add(id);
433  cvf_sets[id] = mesh_region(*this,id);
434  }
435  return cvf_sets[id];
436  }
437  /** Return true if the region of index 's' exists in the mesh */
438  bool has_region(size_type s) const { return valid_cvf_sets[s]; }
439  void add_face_to_set(size_type s, size_type c, short_type f) IS_DEPRECATED;
440  const dal::bit_vector &regions_index() const { return valid_cvf_sets; }
441 
442  /*
443  void sup_face_from_set(size_type b, size_type c, short_type f)
444  { if (valid_cvf_sets[b]) { cvf_sets[b].sup(c,f); touch(); } }*/
445  /** Remove the region of index b. */
447  if (valid_cvf_sets[b])
448  { cvf_sets[b].clear(); valid_cvf_sets.sup(b); touch(); }
449  }
450  /** Remove all references to a convex from all regions stored in the mesh.
451  @param cv the convex number.*/
452  void sup_convex_from_regions(size_type cv);
453  /** Pack the mesh : renumber convexes and nodes such that there
454  is no holes in their numbering. Do NOT do the Cuthill-McKee. */
455  void optimize_structure(bool with_renumbering = true);
456  /// Return the list of convex IDs for a Cuthill-McKee ordering
457  const std::vector<size_type> &cuthill_mckee_ordering() const;
458  /// Erase the mesh.
459  void clear();
460  /** Write the mesh to a file. The format is getfem-specific.
461  @param name the file name.
462  */
463  void write_to_file(const std::string &name) const;
464  /** Write the mesh to a text stream.
465  @param ost the stream.
466  */
467  void write_to_file(std::ostream &ost) const;
468  /** Load the mesh from a file.
469  @param name the file name.
470  @see getfem::import_mesh.
471  */
472  void read_from_file(const std::string &name);
473  /** Load the mesh from a text stream.
474  @param ist the text stream.
475  @see getfem::import_mesh.
476  */
477  void read_from_file(std::istream &ist);
478  /** Clone a mesh */
479  void copy_from(const mesh& m); /* might be the copy constructor */
480  size_type memsize() const;
481 
482  friend class mesh_region;
483  private:
484  void swap_convex_in_regions(size_type c1, size_type c2);
485  void touch_from_region(size_type /*id*/) { touch(); }
486  void to_edges() {} /* to be done, the to_edges of mesh_structure does */
487  /* not handle geotrans */
488 
489  //
490  // Bank refinement
491  //
492 
493  // TODO : - mise a jour des mesh_regions ;
494  // - dialogue avec mesh_im et mesh_fem
495 
496  struct green_simplex {
498  std::vector<size_type> sub_simplices;
500  std::vector<size_type> ipt_loc;
501  };
502  struct edge {
503  size_type i0, i1, i2;
504  bool operator <(const edge &e) const;
505  edge(size_type a, size_type b) : i0(0), i1(a), i2(b) {}
506  edge(size_type a, size_type b, size_type c) : i0(a), i1(b), i2(c) {}
507  };
508  typedef std::set<edge> edge_set;
509 
510  struct Bank_info_struct {
511  dal::bit_vector is_green_simplex; // indices of green simplices.
512  std::map<size_type, size_type> num_green_simplex;
513  dal::dynamic_tas<green_simplex> green_simplices;
514  edge_set edges;
515  };
516 
517  std::unique_ptr<Bank_info_struct> Bank_info;
518 
519  std::string name_; //optional name of the mesh
520  void set_name(const std::string&);
521 
522  void Bank_convex_with_edge(size_type, size_type,
523  std::vector<size_type> &);
524  bool Bank_is_convex_having_points(size_type,
525  const std::vector<size_type> &);
526  void Bank_sup_convex_from_green(size_type);
527  void Bank_swap_convex(size_type, size_type);
528  void Bank_build_first_mesh(mesh &, size_type);
529  void Bank_basic_refine_convex(size_type);
530  void Bank_refine_normal_convex(size_type);
531  size_type Bank_test_and_refine_convex(size_type, dal::bit_vector &,
532  bool = true);
533  void Bank_build_green_simplexes(size_type, std::vector<size_type> &);
534 
535  public :
536 
537  /** Use the Bank strategy to refine some convexes. */
538  void Bank_refine(dal::bit_vector);
539 
540  };
541 
542  /* Dummy mesh for default parameter of functions. */
543  const mesh &dummy_mesh();
544 
545 
546  inline void APIDECL mesh::add_face_to_set(size_type s, size_type c, short_type f) {
547  region(s).add(c,f);
548  }
549 
550  /**
551  * build a N+1 dimensions mesh from a N-dimensions mesh by extrusion.
552  */
553  void APIDECL extrude(const mesh& in, mesh& out, size_type nb_layers,
554  short_type degree=short_type(1));
555 
556  template<class ITER>
558  ITER ipts, const scalar_type tol)
559  {
560  short_type nb = short_type(pgt->nb_points());
561  std::vector<size_type> ind(nb);
562  for (short_type i = 0; i < nb; ++ipts, ++i) ind[i] = add_point(*ipts, tol);
563  return add_convex(pgt, ind.begin());
564  }
565 
566  template<class ITER>
567  size_type mesh::add_simplex_by_points(dim_type di, ITER ipts)
568  { return add_convex_by_points(bgeot::simplex_geotrans(di, 1), ipts); }
569 
570  template<class ITER>
571  size_type mesh::add_parallelepiped(dim_type di, const ITER &ipts)
572  { return add_convex(bgeot::parallelepiped_geotrans(di, 1), ipts); }
573 
574  template<class ITER>
576  (dim_type di, const ITER &ps)
577  { return add_convex_by_points(bgeot::parallelepiped_geotrans(di, 1), ps); }
578 
579  template<class ITER>
580  size_type mesh::add_prism(dim_type di, const ITER &ipts)
581  { return add_convex(bgeot::prism_geotrans(di, 1), ipts); }
582 
583  template<class ITER>
585  (dim_type di, const ITER &ps)
586  { return add_convex_by_points(bgeot::prism_geotrans(di, 1), ps); }
587 
588  /** rough estimate of the convex area.
589  @param pgt the geometric transformation.
590  @param pts the convex nodes.
591  @param pai the approximate integration used for the computation
592  of the convex area.
593  */
594  scalar_type APIDECL convex_area_estimate(bgeot::pgeometric_trans pgt,
595  const base_matrix& pts,
596  pintegration_method pim);
597 
598  /** rough estimate of the maximum value of the condition
599  * number of the jacobian of the geometric transformation */
600  scalar_type APIDECL convex_quality_estimate(bgeot::pgeometric_trans pgt,
601  const base_matrix& pts);
602 
603  /** rough estimate of the radius of the convex using the largest eigenvalue
604  * of the jacobian of the geometric transformation */
605  scalar_type APIDECL convex_radius_estimate(bgeot::pgeometric_trans pgt,
606  const base_matrix& pts);
607 
608  /* stores a convex face. if f == -1, it is the whole convex. */
609  struct convex_face;
610  typedef std::vector<convex_face> convex_face_ct;
611  struct APIDECL convex_face {
612  size_type cv;
613  short_type f;
614  inline bool operator < (const convex_face &e) const {
615  if (cv < e.cv) return true;
616  if (cv > e.cv) return false;
617  if (f < e.f) return true; else if (f > e.f) return false;
618  return false;
619  }
620  bool is_face() const { return f != short_type(-1); }
621  convex_face(size_type cv_, short_type f_=short_type(-1)) : cv(cv_), f(f_) {}
622  convex_face() : cv(size_type(-1)), f(short_type(-1)) {}
623  }; // IS_DEPRECATED;
624 
625  /** returns a list of "exterior" faces of a mesh
626  * (i.e. faces which are not shared by two convexes)
627  * + convexes whose dimension is smaller that m.dim()
628  */
629  void APIDECL outer_faces_of_mesh(const mesh &m, const dal::bit_vector &cvlst,
630  convex_face_ct &flist);
631 
632  inline void outer_faces_of_mesh(const mesh &m, convex_face_ct &flist)
633  { outer_faces_of_mesh(m, m.convex_index(), flist); }
634 
635  void APIDECL outer_faces_of_mesh(const mesh &m, const mesh_region &cvlst,
636  mesh_region &flist);
637 
638  inline void APIDECL outer_faces_of_mesh(const mesh &m, mesh_region &flist)
639  { outer_faces_of_mesh(m, m.convex_index(), flist); }
640 
641 
642  inline mesh_region APIDECL outer_faces_of_mesh(const mesh &m) {
643  mesh_region fl;
644  outer_faces_of_mesh(m, m.convex_index(), fl);
645  return fl;
646  }
647 
648  /** Select all the faces sharing at least two element of the given mesh
649  region. Each face is represented only once and is arbitrarily chosen
650  between the two neighbor elements.
651  */
652  mesh_region APIDECL
653  inner_faces_of_mesh(const mesh &m,
654  const mesh_region &mr = mesh_region::all_convexes());
655 
656  /** Select all the faces of the given mesh region. The faces are represented*
657  twice if they are shared by two neighbor elements.
658  */
659  mesh_region APIDECL
660  all_faces_of_mesh(const mesh &m,
661  const mesh_region &mr = mesh_region::all_convexes());
662 
663  /** Select in the region mr the faces of the mesh m with their unit
664  outward vector having a maximal angle "angle" with the vector V.
665  */
666  mesh_region APIDECL
667  select_faces_of_normal(const mesh &m,
668  const mesh_region &mr,
669  const base_small_vector &V,
670  scalar_type angle);
671 
672  /** Select in the region mr the faces of the mesh m lying entirely in the
673  box delimated by pt1 and pt2.
674  */
675  mesh_region APIDECL
676  select_faces_in_box(const mesh &m, const mesh_region &mr,
677  const base_node &pt1,
678  const base_node &pt2);
679 
680  /** Select in the region mr the faces of the mesh m lying entirely in the
681  ball delimated by pt1 and radius.
682  */
683  mesh_region APIDECL
684  select_faces_in_ball(const mesh &m, const mesh_region &mr,
685  const base_node &center,
686  scalar_type radius);
687 
688  mesh_region APIDECL
689  select_convexes_in_box(const mesh &m, const mesh_region &mr,
690  const base_node &pt1,
691  const base_node &pt2);
692 
693  inline mesh_region APIDECL
694  select_convexes_in_box(const mesh &m,
695  const base_node &pt1,
696  const base_node &pt2)
697  { return select_convexes_in_box(m, m.region(-1), pt1, pt2); }
698 
699  ///@}
700 } /* end of namespace getfem. */
701 
702 
703 #endif /* GETFEM_MESH_H__ */
"File Tools"
Inversion of geometric transformations.
Basic mesh definition.
generic definition of a convex ( bgeot::convex_structure + vertices coordinates )
Definition: bgeot_convex.h:49
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).
Definition: getfem_mesh.h:98
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.
Definition: getfem_mesh.h:195
size_type nb_points() const
Give the number of geometrical nodes in the mesh.
Definition: getfem_mesh.h:193
ref_convex convex(size_type ic) const
return a bgeot::convex object for the convex number ic.
Definition: getfem_mesh.h:188
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.
Definition: getfem_mesh.h:181
void update_from_context() const
this function has to be defined and should update the object when the context is modified.
Definition: getfem_mesh.h:171
size_type search_point(const base_node &pt, const scalar_type radius=0) const
Search a point given its coordinates.
Definition: getfem_mesh.h:208
void sup_region(size_type b)
Remove the region of index b.
Definition: getfem_mesh.h:446
gmm::uint64_type convex_version_number(size_type ic) const
return the version number of the convex ic.
Definition: getfem_mesh.h:214
const mesh_region region(size_type id) const
Return the region of index 'id'.
Definition: getfem_mesh.h:420
size_type add_convex(bgeot::pgeometric_trans pgt, ITER ipts)
Add a convex to the mesh.
Definition: getfem_mesh.h:225
void sup_point(size_type i)
Delete the point of index i from the mesh if it is not linked to a convex.
Definition: getfem_mesh.h:199
bool has_region(size_type s) const
Return true if the region of index 's' exists in the mesh.
Definition: getfem_mesh.h:438
void swap_points(size_type i, size_type j)
Swap the indexes of points of index i and j in the whole structure.
Definition: getfem_mesh.h:201
size_type add_simplex(dim_type di, ITER ipts)
Add a simplex to the mesh, given its dimension and point numbers.
Definition: getfem_mesh.h:251
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.
Definition: getfem_mesh.h:261
indexed array reference (given a container X, and a set of indexes I, this class provides a pseudo-co...
Definition: gmm_ref.h:303
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...
Definition: getfem_mesh.cc:797
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...
Definition: getfem_mesh.cc:771
size_type add_simplex_by_points(dim_type dim, ITER ipts)
Add a simplex to the mesh, given its dimension and point coordinates.
Definition: getfem_mesh.h:567
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...
Definition: getfem_mesh.cc:922
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.
Definition: getfem_mesh.h:557
mesh_region APIDECL select_faces_in_ball(const mesh &m, const mesh_region &mr, const base_node &center, scalar_type radius)
Select in the region mr the faces of the mesh m lying entirely in the ball delimated by pt1 and radiu...
Definition: getfem_mesh.cc:959
scalar_type APIDECL convex_area_estimate(bgeot::pgeometric_trans pgt, const base_matrix &pts, pintegration_method pim)
rough estimate of the convex area.
Definition: getfem_mesh.cc:744
size_type add_prism(dim_type di, const ITER &ipts)
Add a prism to the mesh.
Definition: getfem_mesh.h:580
size_type add_parallelepiped(dim_type di, const ITER &ipts)
Add a parallelepiped to the mesh.
Definition: getfem_mesh.h:571
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.
Definition: getfem_mesh.cc:876
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.
Definition: getfem_mesh.cc:821
size_type add_parallelepiped_by_points(dim_type di, const ITER &ps)
Add a parallelepiped to the mesh.
Definition: getfem_mesh.h:576
size_type add_prism_by_points(dim_type di, const ITER &ps)
Add a prism to the mesh.
Definition: getfem_mesh.h:585
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.
Definition: getfem_mesh.cc:936
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.
Definition: getfem_mesh.cc:857
gmm::uint16_type short_type
used as the common short type integer in the library
Definition: bgeot_config.h:72
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:48
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
GEneric Tool for Finite Element Methods.