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