GetFEM++  5.3
getfem_mesh.cc
1 /*===========================================================================
2 
3  Copyright (C) 1999-2017 Yves Renard
4 
5  This file is a part of GetFEM++
6 
7  GetFEM++ is free software; you can redistribute it and/or modify it
8  under the terms of the GNU Lesser General Public License as published
9  by the Free Software Foundation; either version 3 of the License, or
10  (at your option) any later version along with the GCC Runtime Library
11  Exception either version 3.1 or (at your option) any later version.
12  This program is distributed in the hope that it will be useful, but
13  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
14  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
15  License and GCC Runtime Library Exception for more details.
16  You should have received a copy of the GNU Lesser General Public License
17  along with this program; if not, write to the Free Software Foundation,
18  Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.
19 
20 ===========================================================================*/
21 
22 #include "getfem/bgeot_torus.h"
23 #include "getfem/dal_singleton.h"
25 #include "getfem/getfem_mesh.h"
27 
28 #if GETFEM_HAVE_METIS_OLD_API
29 extern "C" void METIS_PartGraphKway(int *, int *, int *, int *, int *, int *,
30  int *, int *, int *, int *, int *);
31 #elif GETFEM_HAVE_METIS
32 # include <metis.h>
33 #endif
34 
35 namespace getfem {
36 
37  gmm::uint64_type act_counter(void) {
38  static gmm::uint64_type c = gmm::uint64_type(1);
39  return ++c;
40  }
41 
43  for (dal::bv_visitor i(valid_cvf_sets); !i.finished(); ++i)
44  cvf_sets[i].sup_all(c);
45  touch();
46  }
47 
48  void mesh::swap_convex_in_regions(size_type c1, size_type c2) {
49  for (dal::bv_visitor i(valid_cvf_sets); !i.finished(); ++i)
50  cvf_sets[i].swap_convex(c1, c2);
51  touch();
52  }
53 
54  void mesh::handle_region_refinement(size_type ic,
55  const std::vector<size_type> &icv,
56  bool refine) {
57 
58  bgeot::pgeometric_trans pgt = trans_of_convex(ic);
59  ind_set s;
60 
61  for (dal::bv_visitor ir(valid_cvf_sets); !ir.finished(); ++ir) {
62  mesh_region &r = cvf_sets[ir];
63 
64 
65  if (refine && r[ic].any()) {
66  if (r[ic][0])
67  for (size_type jc = 0; jc < icv.size(); ++jc)
68  r.add(icv[jc]);
69 
71  for (short_type f = 0; f < pgt->structure()->nb_faces(); ++f) {
72  if (r[ic][f+1]) {
73 
74  for (size_type jc = 0; jc < icv.size(); ++jc) {
75  bgeot::pgeometric_trans pgtsub = trans_of_convex(icv[jc]);
76  for (short_type fsub = 0; fsub < pgtsub->structure()->nb_faces();
77  ++fsub) {
78  neighbours_of_convex(icv[jc], fsub, s);
79  ind_set::const_iterator it = s.begin(), ite = s.end();
80  bool found = false;
81  for (; it != ite; ++it)
82  if (std::find(icv.begin(), icv.end(), *it) != icv.end())
83  { found = true; break; }
84  if (found) continue;
85 
86  base_node pt, barycentre
87  =gmm::mean_value(pgtsub->convex_ref()->points_of_face(fsub));
88  pt = pgtsub->transform(barycentre, points_of_convex(icv[jc]));
89 
90  giv.init(points_of_convex(ic), pgt);
91  giv.invert(pt, barycentre);
92 
93 
94  if (gmm::abs(pgt->convex_ref()->is_in_face(f,barycentre)) < 0.001)
95  r.add(icv[jc], fsub);
96  }
97  }
98  }
99  }
100  }
101 
102  for (size_type jc = 0; jc < icv.size(); ++jc)
103  if (!refine && r[icv[jc]].any()) {
104  if (r[icv[jc]][0])
105  r.add(ic);
107  bgeot::pgeometric_trans pgtsub = trans_of_convex(icv[jc]);
108  for (short_type fsub = 0; fsub < pgtsub->structure()->nb_faces();
109  ++fsub)
110  if (r[icv[jc]][fsub+1]) {
111  base_node pt, barycentre
112  = gmm::mean_value(pgtsub->convex_ref()->points_of_face(fsub));
113  pt = pgtsub->transform(barycentre, points_of_convex(icv[jc]));
114 
115  giv.init(points_of_convex(ic), pgt);
116  giv.invert(pt, barycentre);
117 
118  for (short_type f = 0; f < pgt->structure()->nb_faces(); ++f)
119  if (gmm::abs(pgt->convex_ref()->is_in_face(f,barycentre)) < 0.001)
120  { r.add(ic, f); break; }
121  }
122  }
123  }
124  }
125 
126  void mesh::init(void) {
127 #if GETFEM_PARA_LEVEL > 1
128  modified = true;
129 #endif
130  cuthill_mckee_uptodate = false;
131  }
132 
133  mesh::mesh(const std::string name) : name_(name) { init(); }
134 
135  mesh::mesh(const bgeot::basic_mesh &m, const std::string name)
136  : bgeot::basic_mesh(m), name_(name) { init(); }
137 
138  mesh::mesh(const mesh &m) : context_dependencies(),
139  std::enable_shared_from_this<getfem::mesh>()
140  { copy_from(m); }
141 
142  mesh &mesh::operator=(const mesh &m) {
143  clear_dependencies();
144  clear();
145  copy_from(m);
146  return *this;
147  }
148 
149 #if GETFEM_PARA_LEVEL > 1
150 
151  void mesh::compute_mpi_region() const {
152  int size, rank;
153  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
154  MPI_Comm_size(MPI_COMM_WORLD, &size);
155 
156  if (size < 2) {
157  mpi_region = mesh_region::all_convexes();
158  mpi_region.from_mesh(*this);
159  } else {
160  int ne = int(nb_convex());
161  std::vector<int> xadj(ne+1), adjncy, numelt(ne), npart(ne);
162  std::vector<int> indelt(nb_allocated_convex());
163 
164  double t_ref = MPI_Wtime();
165 
166  int j = 0, k = 0;
167  ind_set s;
168  for (dal::bv_visitor ic(convex_index()); !ic.finished(); ++ic, ++j) {
169  numelt[j] = ic;
170  indelt[ic] = j;
171  }
172  j = 0;
173  for (dal::bv_visitor ic(convex_index()); !ic.finished(); ++ic, ++j) {
174  xadj[j] = k;
175  neighbours_of_convex(ic, s);
176  for (ind_set::iterator it = s.begin();
177  it != s.end(); ++it) { adjncy.push_back(indelt[*it]); ++k; }
178  }
179  xadj[j] = k;
180 
181 #ifdef GETFEM_HAVE_METIS_OLD_API
182  int wgtflag = 0, numflag = 0, edgecut;
183  int options[5] = {0,0,0,0,0};
184  METIS_PartGraphKway(&ne, &(xadj[0]), &(adjncy[0]), 0, 0, &wgtflag,
185  &numflag, &size, options, &edgecut, &(npart[0]));
186 #else
187  int ncon = 1, edgecut;
188  int options[METIS_NOPTIONS] = { 0 };
189  METIS_SetDefaultOptions(options);
190  METIS_PartGraphKway(&ne, &ncon, &(xadj[0]), &(adjncy[0]), 0, 0, 0,
191  &size, 0, 0, options, &edgecut, &(npart[0]));
192 #endif
193 
194  for (size_type i = 0; i < size_type(ne); ++i)
195  if (npart[i] == rank) mpi_region.add(numelt[i]);
196 
197  if (MPI_IS_MASTER())
198  cout << "Partition time "<< MPI_Wtime()-t_ref << endl;
199  }
200  modified = false;
201  valid_sub_regions.clear();
202  }
203 
204  void mesh::compute_mpi_sub_region(size_type n) const {
205  if (valid_cvf_sets.is_in(n)) {
206  mpi_sub_region[n] = mesh_region::intersection(cvf_sets[n], mpi_region);
207  }
208  else
209  mpi_sub_region[n] = mesh_region();
210  valid_sub_regions.add(n);
211  }
212 
213  void mesh::intersect_with_mpi_region(mesh_region &rg) const {
214  if (rg.id() == mesh_region::all_convexes().id()) {
215  rg = get_mpi_region();
216  } else if (int(rg.id()) >= 0) {
217  rg = get_mpi_sub_region(rg.id());
218  } else
219  rg = mesh_region::intersection(rg, get_mpi_region());
220  }
221 #endif
222 
223  void mesh::optimize_structure(bool with_renumbering) {
224  pts.resort();
225  size_type i, j = nb_convex(), nbc = j;
226  for (i = 0; i < j; i++)
227  if (!convex_tab.index_valid(i))
228  swap_convex(i, convex_tab.ind_last());
229  if (pts.size())
230  for (i = 0, j = pts.size()-1;
231  i < j && j != ST_NIL; ++i, --j) {
232  while (i < j && j != ST_NIL && pts.index()[i]) ++i;
233  while (i < j && j != ST_NIL && !(pts.index()[j])) --j;
234  if (i < j && j != ST_NIL ) swap_points(i, j);
235  }
236  if (with_renumbering) { // Could be optimized no using only swap_convex
237  std::vector<size_type> cmk, iord(nbc), iordinv(nbc);
238  for (i = 0; i < nbc; ++i) iord[i] = iordinv[i] = i;
239 
241  for (i = 0; i < nbc; ++i) {
242  j = iordinv[cmk[i]];
243  if (i != j) {
244  swap_convex(i, j);
245  std::swap(iord[i], iord[j]);
246  std::swap(iordinv[iord[i]], iordinv[iord[j]]);
247  }
248  }
249  }
250  }
251 
252  void mesh::translation(const base_small_vector &V)
253  { pts.translation(V); touch(); }
254 
255  void mesh::transformation(const base_matrix &M) {
256  pts.transformation(M);
257  Bank_info = std::unique_ptr<Bank_info_struct>();
258  touch();
259  }
260 
261  void mesh::bounding_box(base_node& Pmin, base_node& Pmax) const {
262  bool is_first = true;
263  Pmin.clear(); Pmax.clear();
264  for (dal::bv_visitor i(pts.index()); !i.finished(); ++i) {
265  if (is_first) { Pmin = Pmax = pts[i]; is_first = false; }
266  else for (dim_type j=0; j < dim(); ++j) {
267  Pmin[j] = std::min(Pmin[j], pts[i][j]);
268  Pmax[j] = std::max(Pmax[j], pts[i][j]);
269  }
270  }
271  }
272 
273  void mesh::clear(void) {
274  mesh_structure::clear();
275  pts.clear();
276  gtab.clear(); trans_exists.clear();
277  cvf_sets.clear(); valid_cvf_sets.clear();
278  cvs_v_num.clear();
279  Bank_info = nullptr;
280  touch();
281  }
282 
284  size_type ipt[2]; ipt[0] = a; ipt[1] = b;
285  return add_convex(bgeot::simplex_geotrans(1, 1), &(ipt[0]));
286  }
287 
289  size_type b, size_type c) {
290  size_type ipt[3]; ipt[0] = a; ipt[1] = b; ipt[2] = c;
291  return add_simplex(2, &(ipt[0]));
292  }
293 
295  (const base_node &p1, const base_node &p2, const base_node &p3)
296  { return add_triangle(add_point(p1), add_point(p2), add_point(p3)); }
297 
299  size_type c, size_type d) {
300  size_type ipt[4]; ipt[0] = a; ipt[1] = b; ipt[2] = c; ipt[3] = d;
301  return add_simplex(3, &(ipt[0]));
302  }
303 
305  size_type c, size_type d, size_type e) {
306  size_type ipt[5] = {a, b, c, d, e};
307  return add_convex(bgeot::pyramid_QK_geotrans(1), &(ipt[0]));
308  }
309 
311  (const base_node &p1, const base_node &p2,
312  const base_node &p3, const base_node &p4) {
313  return add_tetrahedron(add_point(p1), add_point(p2),
314  add_point(p3), add_point(p4));
315  }
316 
317  void mesh::sup_convex(size_type ic, bool sup_points) {
318  static std::vector<size_type> ipt;
319  if (sup_points) {
320  const ind_cv_ct &ct = ind_points_of_convex(ic);
321  ipt.assign(ct.begin(), ct.end());
322  }
324  if (sup_points)
325  for (size_type ip = 0; ip < ipt.size(); ++ip) sup_point(ipt[ip]);
326  trans_exists.sup(ic);
328  if (Bank_info.get()) Bank_sup_convex_from_green(ic);
329  touch();
330  }
331 
333  if (i != j) {
335  trans_exists.swap(i, j);
336  gtab.swap(i,j);
337  swap_convex_in_regions(i, j);
338  if (Bank_info.get()) Bank_swap_convex(i,j);
339  cvs_v_num[i] = cvs_v_num[j] = act_counter(); touch();
340  }
341  }
342 
344  const base_node &pt) const {
345  bgeot::pgeometric_trans pgt = trans_of_convex(ic);
346  base_matrix G(dim(),pgt->nb_points());
347  vectors_to_base_matrix(G, points_of_convex(ic));
348  bgeot::geotrans_interpolation_context c(trans_of_convex(ic), pt, G);
349  return bgeot::compute_normal(c, f);
350  }
351 
352 
353  base_small_vector mesh::normal_of_face_of_convex(size_type ic, short_type f,
354  size_type n) const {
355  bgeot::pgeometric_trans pgt = trans_of_convex(ic);
356  bgeot::pgeotrans_precomp pgp
357  = bgeot::geotrans_precomp(pgt, pgt->pgeometric_nodes(), 0);
358  base_matrix G;
359  vectors_to_base_matrix(G, points_of_convex(ic));
361  c(pgp,pgt->structure()->ind_points_of_face(f)[n], G);
362  return bgeot::compute_normal(c, f);
363  }
364 
365  base_small_vector mesh::mean_normal_of_face_of_convex(size_type ic,
366  short_type f) const {
367  bgeot::pgeometric_trans pgt = trans_of_convex(ic);
368  base_matrix G; vectors_to_base_matrix(G,points_of_convex(ic));
369  base_small_vector ptmean(dim());
370  size_type nbpt = pgt->structure()->nb_points_of_face(f);
371  for (size_type i = 0; i < nbpt; ++i)
372  gmm::add(pgt->geometric_nodes()[pgt->structure()->ind_points_of_face(f)[i]], ptmean);
373  ptmean /= scalar_type(nbpt);
374  bgeot::geotrans_interpolation_context c(pgt,ptmean, G);
375  base_small_vector n = bgeot::compute_normal(c, f);
376  n /= gmm::vect_norm2(n);
377  return n;
378  }
379 
381  const base_node &pt) const {
382  bgeot::pgeometric_trans pgt = trans_of_convex(ic);
383  base_matrix G(dim(),pgt->nb_points());
384  vectors_to_base_matrix(G,points_of_convex(ic));
385  bgeot::geotrans_interpolation_context c(trans_of_convex(ic), pt, G);
386  return bgeot::compute_local_basis(c, f);
387  }
388 
390  size_type n) const {
391  bgeot::pgeometric_trans pgt = trans_of_convex(ic);
392  bgeot::pgeotrans_precomp pgp
393  = bgeot::geotrans_precomp(pgt, pgt->pgeometric_nodes(), 0);
394  base_matrix G(dim(),pgt->nb_points());
395  vectors_to_base_matrix(G,points_of_convex(ic));
397  c(pgp,pgt->structure()->ind_points_of_face(f)[n], G);
398  return bgeot::compute_local_basis(c, f);
399  }
400 
401  scalar_type mesh::convex_area_estimate(size_type ic, size_type deg) const {
402  base_matrix G;
403  bgeot::vectors_to_base_matrix(G, points_of_convex(ic));
405  (trans_of_convex(ic), G, classical_approx_im(trans_of_convex(ic),
406  dim_type(deg)));
407  }
408 
409  scalar_type mesh::convex_quality_estimate(size_type ic) const {
410  base_matrix G;
411  bgeot::vectors_to_base_matrix(G, points_of_convex(ic));
412  auto pgt = trans_of_convex(ic);
413  if (auto pgt_torus = dynamic_cast<const bgeot::torus_geom_trans*>(pgt.get())) {
414  pgt = pgt_torus->get_original_transformation();
415  G.resize(pgt->dim(), G.ncols());
416  }
417  return getfem::convex_quality_estimate(pgt, G);
418  }
419 
420  scalar_type mesh::convex_radius_estimate(size_type ic) const {
421  base_matrix G;
422  bgeot::vectors_to_base_matrix(G, points_of_convex(ic));
423  return getfem::convex_radius_estimate(trans_of_convex(ic), G);
424  }
425 
427  if (convex_index().empty()) return 1;
428  scalar_type r = convex_radius_estimate(convex_index().first_true());
429  for (dal::bv_visitor cv(convex_index()); !cv.finished(); ++cv) {
430  r = std::min(r, convex_radius_estimate(cv));
431  }
432  return r;
433  }
434 
436  if (convex_index().empty()) return 1;
437  scalar_type r = convex_radius_estimate(convex_index().first_true());
438  for (dal::bv_visitor cv(convex_index()); !cv.finished(); ++cv) {
439  r = std::max(r, convex_radius_estimate(cv));
440  }
441  return r;
442  }
443 
444  void mesh::set_name(const std::string& name){name_=name;}
445 
446  void mesh::copy_from(const mesh& m) {
447  clear();
448  set_name(m.name_);
449  bgeot::basic_mesh::operator=(m);
450  cvf_sets = m.cvf_sets;
451  valid_cvf_sets = m.valid_cvf_sets;
452  for (std::map<size_type, mesh_region>::iterator it = cvf_sets.begin();
453  it != cvf_sets.end(); ++it)
454  if (it->second.get_parent_mesh() != 0) it->second.set_parent_mesh(this);
455  cvs_v_num.clear();
456  gmm::uint64_type d = act_counter();
457  for (dal::bv_visitor i(convex_index()); !i.finished(); ++i)
458  cvs_v_num[i] = d;
459  Bank_info = std::unique_ptr<Bank_info_struct>();
460  if (m.Bank_info.get())
461  Bank_info = std::make_unique<Bank_info_struct>(*(m.Bank_info));
462  }
463 
464  struct mesh_convex_structure_loc {
465  bgeot::pgeometric_trans cstruct;
466  std::vector<size_type> pts;
467  };
468 
469  void mesh::read_from_file(std::istream &ist) {
470  gmm::stream_standard_locale sl(ist);
471  dal::bit_vector npt;
473  std::string tmp;
474  bool te = false, please_get = true;
475 
476  ist.precision(16);
477  clear();
478  ist.seekg(0);ist.clear();
479  bgeot::read_until(ist, "BEGIN POINTS LIST");
480 
481  while (!te) {
482  if (please_get) bgeot::get_token(ist, tmp); else please_get = true;
483 
484  if (!bgeot::casecmp(tmp, "END"))
485  { te = true; }
486  else if (!bgeot::casecmp(tmp, "POINT")) {
487  bgeot::get_token(ist, tmp);
488  if (!bgeot::casecmp(tmp, "COUNT")) {
489  bgeot::get_token(ist, tmp); // Ignored. Used in some applications
490  } else {
491  size_type ip = atoi(tmp.c_str());
492  dim_type d = 0;
493  GMM_ASSERT1(!npt.is_in(ip),
494  "Two points with the same index. loading aborted.");
495  npt.add(ip);
496  bgeot::get_token(ist, tmp);
497  while (isdigit(tmp[0]) || tmp[0] == '-' || tmp[0] == '+'
498  || tmp[0] == '.')
499  { tmpv[d++] = atof(tmp.c_str()); bgeot::get_token(ist, tmp); }
500  please_get = false;
501  base_node v(d);
502  for (size_type i = 0; i < d; i++) v[i] = tmpv[i];
503  size_type ipl = add_point(v);
504  if (ip != ipl) {
505  GMM_ASSERT1(!npt.is_in(ipl), "Two points [#" << ip << " and #"
506  << ipl << "] with the same coords "<< v
507  << ". loading aborted.");
508  swap_points(ip, ipl);
509  }
510  }
511  } else if (tmp.size()) {
512  GMM_ASSERT1(false, "Syntax error in file, at token '" << tmp
513  << "', pos=" << std::streamoff(ist.tellg()));
514  } else if (ist.eof()) {
515  GMM_ASSERT1(false, "Unexpected end of stream while reading mesh");
516  }
517  }
518 
519  bool tend = false;
521  dal::bit_vector ncv;
522 
523  ist.seekg(0);
524  if (!bgeot::read_until(ist, "BEGIN MESH STRUCTURE DESCRIPTION"))
525  GMM_ASSERT1(false, "This seems not to be a mesh file");
526 
527  while (!tend) {
528  tend = !bgeot::get_token(ist, tmp);
529  if (!bgeot::casecmp(tmp, "END"))
530  { tend = true; }
531  else if (!bgeot::casecmp(tmp, "CONVEX")) {
532  size_type ic;
533  bgeot::get_token(ist, tmp);
534  if (!bgeot::casecmp(tmp, "COUNT")) {
535  bgeot::get_token(ist, tmp); // Ignored. Used in some applications
536  } else {
537  ic = gmm::abs(atoi(tmp.c_str()));
538  GMM_ASSERT1(!ncv.is_in(ic),
539  "Negative or repeated index, loading aborted.");
540  ncv.add(ic);
541 
542  int rgt = bgeot::get_token(ist, tmp);
543  if (rgt != 3) { // for backward compatibility with version 1.7
544  char c; ist.get(c);
545  while (!isspace(c)) { tmp.push_back(c); ist.get(c); }
546  }
547 
549  size_type nb = pgt->nb_points();
550 
551  cv[ic].cstruct = pgt;
552  cv[ic].pts.resize(nb);
553  for (size_type i = 0; i < nb; i++) {
554  bgeot::get_token(ist, tmp);
555  cv[ic].pts[i] = gmm::abs(atoi(tmp.c_str()));
556  }
557  }
558  }
559  else if (tmp.size()) {
560  GMM_ASSERT1(false, "Syntax error reading a mesh file "
561  " at pos " << std::streamoff(ist.tellg())
562  << "(expecting 'CONVEX' or 'END', found '"<< tmp << "')");
563  } else if (ist.eof()) {
564  GMM_ASSERT1(false, "Unexpected end of stream "
565  << "(missing BEGIN MESH/END MESH ?)");
566  }
567  }
568  ist >> bgeot::skip("MESH STRUCTURE DESCRIPTION");
569 
570  for (dal::bv_visitor ic(ncv); !ic.finished(); ++ic) {
571  size_type i = add_convex(cv[ic].cstruct, cv[ic].pts.begin());
572  if (i != ic) swap_convex(i, ic);
573  }
574 
575  tend = false;
576  while (!tend) {
577  tend = !bgeot::get_token(ist, tmp);
578  // bool error = false;
579  if (bgeot::casecmp(tmp, "BEGIN")==0) {
580  bgeot::get_token(ist, tmp);
581  if (bgeot::casecmp(tmp, "BOUNDARY")==0 ||
582  bgeot::casecmp(tmp, "REGION")==0) {
583  bgeot::get_token(ist, tmp);
584  size_type bnum = atoi(tmp.c_str());
585  bgeot::get_token(ist, tmp);
586  while (true) {
587  if (bgeot::casecmp(tmp, "END")!=0) {
588  size_type ic = atoi(tmp.c_str());
589  bgeot::get_token(ist, tmp);
590  if (tmp[0] == '/') {
591  bgeot::get_token(ist, tmp);
592  if (!bgeot::casecmp(tmp, "END")) break;
593  int f = atoi(tmp.c_str());
594  region(bnum).add(ic, short_type(f));
595  bgeot::get_token(ist, tmp);
596  }
597  else {
598  region(bnum).add(ic);
599  if (!bgeot::casecmp(tmp, "END")) break;
600  }
601  } else break;
602  }
603  bgeot::get_token(ist, tmp);
604  bgeot::get_token(ist, tmp);
605  } else tend = true;
606  /*else GMM_ASSERT1(false, "Syntax error in file at token '"
607  << tmp << "' [pos=" << std::streamoff(ist.tellg())
608  << "]");*/
609  } else tend=true;
610  }
611  }
612 
613  void mesh::read_from_file(const std::string &name) {
614  std::ifstream o(name.c_str());
615  GMM_ASSERT1(o, "Mesh file '" << name << "' does not exist");
616  read_from_file(o);
617  o.close();
618  }
619 
620  template<class ITER>
621  void write_tab_to_file_(std::ostream &ost, const ITER& b_, const ITER& e)
622  { for (ITER b(b_) ; b != e; ++b) ost << " " << *b; }
623 
624  template<class ITER>
625  static void write_convex_to_file_(const mesh &ms,
626  std::ostream &ost,
627  ITER b, ITER e) {
628  for ( ; b != e; ++b) {
629  size_type i = b.index();
630  ost << " CONVEX " << i << " \'"
631  << bgeot::name_of_geometric_trans(ms.trans_of_convex(i)).c_str()
632  << "\' ";
633  write_tab_to_file_(ost, ms.ind_points_of_convex(i).begin(),
634  ms.ind_points_of_convex(i).end() );
635  ost << '\n';
636  }
637  }
638 
639  template<class ITER> static void write_point_to_file_(std::ostream &ost,
640  ITER b, ITER e)
641  { for ( ; b != e; ++b) ost << " " << *b; ost << '\n'; }
642 
643  void mesh::write_to_file(std::ostream &ost) const {
644  ost.precision(16);
645  gmm::stream_standard_locale sl(ost);
646  ost << '\n' << "BEGIN POINTS LIST" << '\n' << '\n';
647  ost << " POINT COUNT " << points().index().last_true()+1 << '\n';
648  for (size_type i = 0; i < points_tab.size(); ++i)
649  if (is_point_valid(i) ) {
650  ost << " POINT " << i;
651  write_point_to_file_(ost, pts[i].begin(), pts[i].end());
652  }
653  ost << '\n' << "END POINTS LIST" << '\n' << '\n' << '\n';
654 
655  ost << '\n' << "BEGIN MESH STRUCTURE DESCRIPTION" << '\n' << '\n';
656  ost << " CONVEX COUNT " << convex_index().last_true()+1 << '\n';
657  write_convex_to_file_(*this, ost, convex_tab.tas_begin(),
658  convex_tab.tas_end());
659  ost << '\n' << "END MESH STRUCTURE DESCRIPTION" << '\n';
660 
661  for (dal::bv_visitor bnum(valid_cvf_sets); !bnum.finished(); ++bnum) {
662  ost << "BEGIN REGION " << bnum << "\n" << region(bnum) << "\n"
663  << "END REGION " << bnum << "\n";
664  }
665  }
666 
667  void mesh::write_to_file(const std::string &name) const {
668  std::ofstream o(name.c_str());
669  GMM_ASSERT1(o, "impossible to write to file '" << name << "'");
670  o << "% GETFEM MESH FILE " << '\n';
671  o << "% GETFEM VERSION " << GETFEM_VERSION << '\n' << '\n' << '\n';
672  write_to_file(o);
673  o.close();
674  }
675 
676  size_type mesh::memsize(void) const {
677  return bgeot::mesh_structure::memsize() - sizeof(bgeot::mesh_structure)
678  + pts.memsize() + (pts.index().last_true()+1)*dim()*sizeof(scalar_type)
679  + sizeof(mesh) + trans_exists.memsize() + gtab.memsize()
680  + valid_cvf_sets.card()*sizeof(mesh_region) + valid_cvf_sets.memsize();
681  }
682 
683  struct equilateral_to_GT_PK_grad_aux : public std::vector<base_matrix> {};
684  static const base_matrix &equilateral_to_GT_PK_grad(dim_type N) {
685  std::vector<base_matrix>
687  if (N > pbm.size()) pbm.resize(N);
688  if (pbm[N-1].empty()) {
689  bgeot::pgeometric_trans pgt = bgeot::simplex_geotrans(N,1);
690  base_matrix Gr(N,N);
691  base_matrix G(N,N+1);
692  vectors_to_base_matrix
693  (G, bgeot::equilateral_simplex_of_reference(N)->points());
694  gmm::mult(G, bgeot::geotrans_precomp
695  (pgt, pgt->pgeometric_nodes(), 0)->grad(0), Gr);
696  gmm::lu_inverse(Gr);
697  pbm[N-1].swap(Gr);
698  }
699  return pbm[N-1];
700  }
701 
702 
704  const base_matrix& G,
705  pintegration_method pi) {
706  double area(0);
707  static bgeot::pgeometric_trans pgt_old = 0;
708  static bgeot::pgeotrans_precomp pgp = 0;
709  static pintegration_method pim_old = 0;
710  papprox_integration pai = get_approx_im_or_fail(pi);
711  if (pgt_old != pgt || pim_old != pi) {
712  pgt_old = pgt;
713  pim_old = pi;
714  pgp = bgeot::geotrans_precomp
715  (pgt, pai->pintegration_points(), pi);
716  }
718  for (size_type i = 0; i < pai->nb_points_on_convex(); ++i) {
719  gic.set_ii(i);
720  area += pai->coeff(i) * gic.J();
721  }
722  return area;
723  }
724 
725  /* TODO : use the geotrans from an "equilateral" reference element to
726  the real element
727  check if the sign of the determinants does change
728  (=> very very bad quality of the convex)
729  */
731  const base_matrix& G) {
732  static bgeot::pgeometric_trans pgt_old = 0;
733  static bgeot::pgeotrans_precomp pgp = 0;
734  if (pgt_old != pgt) {
735  pgt_old=pgt;
736  pgp=bgeot::geotrans_precomp(pgt, pgt->pgeometric_nodes(), 0);
737  }
738 
739  size_type n = (pgt->is_linear()) ? 1 : pgt->nb_points();
740  scalar_type q = 1;
741  dim_type N = dim_type(G.nrows()), P = pgt->structure()->dim();
742  base_matrix K(N,P);
743  for (size_type ip=0; ip < n; ++ip) {
744  gmm::mult(G, pgp->grad(ip), K);
745  /* TODO : this is an ugly fix for simplexes only.. there should be
746  a transformation of any pgt to the equivalent equilateral pgt
747  (for prisms etc) */
748  if (bgeot::basic_structure(pgt->structure())
750  gmm::mult(base_matrix(K),equilateral_to_GT_PK_grad(P),K);
751  q = std::max(q, gmm::condition_number(K));
752  }
753  return 1./q;
754  }
755 
757  const base_matrix& G) {
758  static bgeot::pgeometric_trans pgt_old = 0;
759  static bgeot::pgeotrans_precomp pgp = 0;
760  if (pgt_old != pgt) {
761  pgt_old=pgt;
762  pgp=bgeot::geotrans_precomp(pgt, pgt->pgeometric_nodes(), 0);
763  }
764  size_type N = G.nrows();
765  size_type n = (pgt->is_linear()) ? 1 : pgt->nb_points();
766  scalar_type q = 0;
767  base_matrix K(pgp->grad(0).ncols(),G.nrows());
768  for (size_type ip=0; ip < n; ++ip) {
769  gmm::mult(gmm::transposed(pgp->grad(ip)), gmm::transposed(G), K);
770  scalar_type emax, emin; gmm::condition_number(K,emax,emin);
771  q = std::max(q, emax);
772  }
773  return q * sqrt(scalar_type(N)) / scalar_type(2);
774  }
775 
776  /* extract faces of convexes which are not shared
777  + convexes whose dimension is smaller that m.dim()
778  */
779  void
780  outer_faces_of_mesh(const getfem::mesh &m, const dal::bit_vector& cvlst,
781  convex_face_ct& flist) {
782  for (dal::bv_visitor ic(cvlst); !ic.finished(); ++ic) {
783  if (m.structure_of_convex(ic)->dim() == m.dim()) {
784  for (short_type f = 0; f < m.structure_of_convex(ic)->nb_faces(); f++) {
785  if (m.neighbour_of_convex(ic,f) == size_type(-1)) {
786  flist.push_back(convex_face(ic,f));
787  }
788  }
789  } else {
790  flist.push_back(convex_face(ic, short_type(-1)));
791  }
792  }
793  }
794 
795  void outer_faces_of_mesh(const mesh &m,
796  const mesh_region &cvlst,
797  mesh_region &flist) {
798  cvlst.error_if_not_convexes();
799  for (mr_visitor i(cvlst); !i.finished(); ++i) {
800  if (m.structure_of_convex(i.cv())->dim() == m.dim()) {
801  for (short_type f = 0; f < m.structure_of_convex(i.cv())->nb_faces();
802  f++) {
803  size_type cv2 = m.neighbour_of_convex(i.cv(), f);
804  if (cv2 == size_type(-1) || !cvlst.is_in(cv2)) {
805  flist.add(i.cv(),f);
806  }
807  }
808  }
809  else flist.add(i.cv());
810  }
811  }
812 
813  /* Select all the faces sharing at least two element of the given mesh
814  region. Each face is represented only once and is arbitrarily chosen
815  between the two neighbour elements. Try to minimize the number of
816  elements.
817  */
819  mesh_region mrr;
820  mr.from_mesh(m);
821  mr.error_if_not_convexes();
822  dal::bit_vector visited;
823  bgeot::mesh_structure::ind_set neighbours;
824 
825  for (mr_visitor i(mr); !i.finished(); ++i) {
826  size_type cv1 = i.cv();
827  short_type nbf = m.structure_of_convex(i.cv())->nb_faces();
828  bool neighbour_visited = false;
829  for (short_type f = 0; f < nbf; ++f) {
830  neighbours.resize(0); m.neighbours_of_convex(cv1, f, neighbours);
831  for (size_type j = 0; j < neighbours.size(); ++j)
832  if (visited.is_in(neighbours[j]))
833  { neighbour_visited = true; break; }
834  }
835  if (!neighbour_visited) {
836  for (short_type f = 0; f < nbf; ++f) {
837  size_type cv2 = m.neighbour_of_convex(cv1, f);
838  if (cv2 != size_type(-1) && mr.is_in(cv2)) mrr.add(cv1,f);
839  }
840  visited.add(cv1);
841  }
842  }
843 
844  for (mr_visitor i(mr); !i.finished(); ++i) {
845  size_type cv1 = i.cv();
846  if (!(visited.is_in(cv1))) {
847  short_type nbf = m.structure_of_convex(i.cv())->nb_faces();
848  for (short_type f = 0; f < nbf; ++f) {
849  neighbours.resize(0); m.neighbours_of_convex(cv1, f, neighbours);
850  bool ok = false;
851  for (size_type j = 0; j < neighbours.size(); ++j) {
852  if (visited.is_in(neighbours[j])) { ok = false; break; }
853  if (mr.is_in(neighbours[j])) { ok = true; }
854  }
855  if (ok) { mrr.add(cv1,f); }
856  }
857  visited.add(cv1);
858  }
859  }
860  return mrr;
861  }
862 
863 
865  const base_small_vector &V,
866  scalar_type angle) {
867  mesh_region mrr;
868  scalar_type threshold = gmm::vect_norm2(V)*cos(angle);
869  for (getfem::mr_visitor i(mr); !i.finished(); ++i)
870  if (i.f() != short_type(-1)) {
871  base_node un = m.mean_normal_of_face_of_convex(i.cv(), i.f());
872  if (gmm::vect_sp(V, un) - threshold >= -1E-8)
873  mrr.add(i.cv(), i.f());
874  }
875  return mrr;
876  }
877 
879  const base_node &pt1, const base_node &pt2) {
880  mesh_region mrr;
881  size_type N = m.dim();
882  GMM_ASSERT1(pt1.size() == N && pt2.size() == N, "Wrong dimensions");
883  for (getfem::mr_visitor i(mr, m); !i.finished(); ++i)
884  if (i.f() != short_type(-1)) {
886  = m.ind_points_of_face_of_convex(i.cv(), i.f());
887 
888  bool is_in = true;
890  it != pt.end(); ++it) {
891  for (size_type j = 0; j < N; ++j)
892  if (m.points()[*it][j] < pt1[j] || m.points()[*it][j] > pt2[j])
893  { is_in = false; break; }
894  if (!is_in) break;
895  }
896  if (is_in) mrr.add(i.cv(), i.f());
897  }
898  return mrr;
899  }
900 
901  mesh_region select_convexes_in_box(const mesh &m, const mesh_region &mr,
902  const base_node &pt1, const base_node &pt2) {
903  mesh_region mrr;
904  size_type N = m.dim();
905  GMM_ASSERT1(pt1.size() == N && pt2.size() == N, "Wrong dimensions");
906  for (getfem::mr_visitor i(mr, m); !i.finished(); ++i)
907  if (i.f() == short_type(-1)) {
908  bgeot::mesh_structure::ind_cv_ct pt = m.ind_points_of_convex(i.cv());
909 
910  bool is_in = true;
911  for (bgeot::mesh_structure::ind_cv_ct::iterator it = pt.begin();
912  it != pt.end(); ++it) {
913  for (size_type j = 0; j < N; ++j)
914  if (m.points()[*it][j] < pt1[j] || m.points()[*it][j] > pt2[j])
915  { is_in = false; break; }
916  if (!is_in) break;
917  }
918  if (is_in) mrr.add(i.cv());
919  }
920  return mrr;
921  }
922 
923  void extrude(const mesh& in, mesh& out, size_type nb_layers, short_type degree) {
924  dim_type dim = in.dim();
925  base_node pt(dim+1);
926  out.clear();
927  size_type nbpt = in.points().index().last()+1;
928  GMM_ASSERT1(nbpt == in.points().index().card(),
929  "please optimize the mesh before using "
930  "it as a base for prismatic mesh");
931  size_type nb_layers_total = nb_layers * degree;
932  for (const base_node &inpt : in.points()) {
933  std::copy(inpt.begin(), inpt.end(), pt.begin());
934  pt[dim] = 0.0;
935  for (size_type j = 0; j <= nb_layers_total;
936  ++j, pt[dim] += scalar_type(1) / scalar_type(nb_layers_total))
937  out.add_point(pt);
938  }
939 
940  std::vector<size_type> tab;
941  for (dal::bv_visitor cv(in.convex_index()); !cv.finished(); ++cv) {
942  size_type nbp = in.nb_points_of_convex(cv);
943  tab.resize((degree+1)*nbp);
944  for (size_type j = 0; j < nb_layers; ++j) {
945  for (short_type d = 0; d < degree+1; ++d)
946  for (size_type k = 0; k < nbp; ++k)
947  tab[k+nbp*d] = (nb_layers_total+1)*in.ind_points_of_convex(cv)[k] + j*degree + d;
949  bgeot::product_geotrans(in.trans_of_convex(cv),
950  bgeot::simplex_geotrans(1,degree));
951  out.add_convex(pgt, tab.begin());
952  }
953  }
954  }
955 }
956 
957 
958 //
959 // Bank refinement
960 //
961 
962 
963 namespace bgeot {
964  size_type refinement_simplexe_tab(size_type n, size_type **tab);
965 }
966 
967 namespace getfem {
968 
969  bool mesh::edge::operator <(const edge &e) const {
970  if (i0 < e.i0) return true;
971  if (i0 > e.i0) return false;
972  if (i1 < e.i1) return true;
973  if (i1 > e.i1) return false;
974  if (i2 < e.i2) return true;
975  return false;
976  }
977 
978  void mesh::Bank_sup_convex_from_green(size_type i) {
979  if (Bank_info.get() && Bank_info->is_green_simplex.is_in(i)) {
980  size_type igs = Bank_info->num_green_simplex[i];
981  green_simplex &gs = Bank_info->green_simplices[igs];
982  for (size_type j = 0; j < gs.sub_simplices.size(); ++j) {
983  Bank_info->num_green_simplex.erase(gs.sub_simplices[j]);
984  Bank_info->is_green_simplex.sup(gs.sub_simplices[j]);
985  }
986  Bank_info->green_simplices.sup(igs);
987  }
988  }
989 
990  void mesh::Bank_swap_convex(size_type i, size_type j) {
991  if (Bank_info.get()) {
992  Bank_info->is_green_simplex.swap(i, j);
993  std::map<size_type, size_type>::iterator
994  iti = Bank_info->num_green_simplex.find(i);
995  std::map<size_type, size_type>::iterator
996  ite = Bank_info->num_green_simplex.end();
997  std::map<size_type, size_type>::iterator
998  itj = Bank_info->num_green_simplex.find(j);
999  size_type numi(0), numj(0);
1000  if (iti != ite)
1001  { numi = iti->second; Bank_info->num_green_simplex.erase(i); }
1002  if (itj != ite)
1003  { numj = itj->second; Bank_info->num_green_simplex.erase(j); }
1004  if (iti != ite) {
1005  Bank_info->num_green_simplex[j] = numi;
1006  green_simplex &gs = Bank_info->green_simplices[numi];
1007  for (size_type k = 0; k < gs.sub_simplices.size(); ++k)
1008  if (gs.sub_simplices[k] == i) gs.sub_simplices[k] = j;
1009  else if (gs.sub_simplices[k] == j) gs.sub_simplices[k] = i;
1010  }
1011  if (itj != ite) {
1012  Bank_info->num_green_simplex[i] = numj;
1013  if (iti == ite || numi != numj) {
1014  green_simplex &gs = Bank_info->green_simplices[numj];
1015  for (size_type k = 0; k < gs.sub_simplices.size(); ++k)
1016  if (gs.sub_simplices[k] == i) gs.sub_simplices[k] = j;
1017  else if (gs.sub_simplices[k] == j) gs.sub_simplices[k] = i;
1018  }
1019  }
1020  }
1021  }
1022 
1023  void mesh::Bank_build_first_mesh(mesh &m, size_type n) {
1024  bgeot::pconvex_ref pcr = bgeot::simplex_of_reference(dim_type(n), 2);
1025  m.clear();
1026  for (size_type ip = 0; ip < pcr->nb_points(); ++ip)
1027  m.add_point(pcr->points()[ip]);
1028  size_type *tab;
1029  size_type nbc = bgeot::refinement_simplexe_tab(n, &tab);
1030  for (size_type ic = 0; ic < nbc; ++ic, tab+=(n+1))
1031  m.add_simplex(dim_type(n), tab);
1032  }
1033 
1034  struct mesh_cache_for_Bank_basic_refine_convex : public mesh {};
1035 
1036  void mesh::Bank_basic_refine_convex(size_type i) {
1037  bgeot::pgeometric_trans pgt = trans_of_convex(i);
1038  size_type n = pgt->basic_structure()->dim();
1039 
1040  static bgeot::pgeometric_trans pgt1 = 0;
1041 
1043 
1044  static bgeot::pstored_point_tab pspt = 0;
1045  static bgeot::pgeotrans_precomp pgp = 0;
1046  static std::vector<size_type> ipt, ipt2, icl;
1047 
1048  if (pgt != pgt1) {
1049  pgt1 = pgt;
1050  mesh mesh1;
1051  Bank_build_first_mesh(mesh1, n);
1052 
1053  mesh2.clear();
1054  ipt.resize(pgt->nb_points());
1055  for (size_type ic = 0; ic < mesh1.nb_convex(); ++ic) {
1056  assert(mesh1.convex_index().is_in(ic));
1057  bgeot::pgeometric_trans pgt2 = mesh1.trans_of_convex(ic);
1058  for (size_type ip = 0; ip < pgt->nb_points(); ++ip)
1059  ipt[ip] =
1060  mesh2.add_point(pgt2->transform(pgt->geometric_nodes()[ip],
1061  mesh1.points_of_convex(ic)));
1062  mesh2.add_convex(pgt, &ipt[0]);
1063  }
1064 
1065  // if (pspt) dal::del_stored_object(pspt);
1066  pspt = bgeot::store_point_tab(mesh2.points());
1067  pgp = bgeot::geotrans_precomp(pgt, pspt, 0);
1068  }
1069 
1070  base_node pt(n);
1071  ipt.resize(pspt->size());
1072  for (size_type ip = 0; ip < pspt->size(); ++ip) {
1073  pgp->transform(points_of_convex(i), ip, pt);
1074  ipt[ip] = add_point(pt);
1075  }
1076 
1077  ipt2.resize(pgt->nb_points()); icl.resize(0);
1078  for (size_type ic = 0; ic < mesh2.nb_convex(); ++ic) {
1079  for (size_type j = 0; j < pgt->nb_points(); ++j)
1080  ipt2[j] = ipt[mesh2.ind_points_of_convex(ic)[j]];
1081  icl.push_back(add_convex(pgt, ipt2.begin()));
1082  }
1083  handle_region_refinement(i, icl, true);
1084  sup_convex(i, true);
1085  }
1086 
1087  void mesh::Bank_convex_with_edge(size_type i1, size_type i2,
1088  std::vector<size_type> &ipt) {
1089  ipt.resize(0);
1090  for (size_type k = 0; k < points_tab[i1].size(); ++k) {
1091  size_type cv = points_tab[i1][k], found = 0;
1092  const std::vector<size_type> &loc_ind = trans_of_convex(cv)->vertices();
1093  for (size_type i = 0; i < loc_ind.size(); ++i) {
1094  size_type ipp = convex_tab[cv].pts[loc_ind[i]];
1095  if (ipp == i1) ++found;
1096  if (ipp == i2) ++found;
1097  }
1098  GMM_ASSERT1(found <= 2, "Invalid convex with repeated nodes ");
1099  if (found == 2) ipt.push_back(cv);
1100  }
1101  }
1102 
1103  bool mesh::Bank_is_convex_having_points(size_type cv,
1104  const std::vector<size_type> &ipt) {
1105  size_type found = 0;
1106  const std::vector<size_type> &loc_ind = trans_of_convex(cv)->vertices();
1107  for (size_type i = 0; i < loc_ind.size(); ++i)
1108  if (std::find(ipt.begin(), ipt.end(),
1109  convex_tab[cv].pts[loc_ind[i]]) != ipt.end()) ++found;
1110  return (found >= ipt.size());
1111  }
1112 
1113  void mesh::Bank_refine_normal_convex(size_type i) {
1114  bgeot::pgeometric_trans pgt = trans_of_convex(i);
1115  GMM_ASSERT1(pgt->basic_structure() == bgeot::simplex_structure(pgt->dim()),
1116  "Sorry, refinement is only working with simplices.");
1117 
1118  const std::vector<size_type> &loc_ind = pgt->vertices();
1119  for (size_type ip1 = 0; ip1 < loc_ind.size(); ++ip1)
1120  for (size_type ip2 = ip1+1; ip2 < loc_ind.size(); ++ip2)
1121  Bank_info->edges.insert(edge(ind_points_of_convex(i)[loc_ind[ip1]],
1122  ind_points_of_convex(i)[loc_ind[ip2]]));
1123  Bank_basic_refine_convex(i);
1124  }
1125 
1126  size_type mesh::Bank_test_and_refine_convex(size_type i,
1127  dal::bit_vector &b, bool ref) {
1128  if (Bank_info->is_green_simplex[i]) {
1129  size_type igs = Bank_info->num_green_simplex[i];
1130  green_simplex &gs = Bank_info->green_simplices[igs];
1131 
1132  size_type icc = add_convex_by_points(gs.pgt, gs.cv.points().begin());
1133  handle_region_refinement(icc, gs.sub_simplices, false);
1134  for (size_type ic = 0; ic < gs.sub_simplices.size(); ++ic) {
1135  sup_convex(gs.sub_simplices[ic], true);
1136  b.sup(gs.sub_simplices[ic]);
1137  }
1138  if (!ref)
1139  for (size_type ip = 0; ip < gs.ipt_loc.size(); ++ip)
1140  for (size_type jp = ip+1; jp < gs.ipt_loc.size(); ++jp)
1141  Bank_info->edges.insert
1142  (edge(ind_points_of_convex(icc)[gs.ipt_loc[ip]],
1143  ind_points_of_convex(icc)[gs.ipt_loc[jp]]));
1144  Bank_sup_convex_from_green(i);
1145  if (ref) { Bank_refine_normal_convex(icc); return size_type(-1); }
1146  else return icc;
1147  }
1148  else if (ref) Bank_refine_normal_convex(i);
1149  return size_type(-1);
1150  }
1151 
1152  struct mesh_cache_for_Bank_build_green_simplexes : public mesh {};
1153 
1154  void mesh::Bank_build_green_simplexes(size_type ic,
1155  std::vector<size_type> &ipt) {
1156  GMM_ASSERT1(convex_index().is_in(ic), "Internal error");
1157  size_type igs = Bank_info->green_simplices.add(green_simplex());
1158  green_simplex &gs(Bank_info->green_simplices[igs]);
1159  std::vector<base_node> pt_tab(nb_points_of_convex(ic));
1160  ref_mesh_pt_ct ptab = points_of_convex(ic);
1161  pt_tab.assign(ptab.begin(), ptab.end());
1162  gs.cv = bgeot::convex<base_node>(structure_of_convex(ic), pt_tab);
1163 
1164  bgeot::pgeometric_trans pgt = gs.pgt = trans_of_convex(ic);
1165 
1166  size_type d = ipt.size() - 1, n = structure_of_convex(ic)->dim();
1167  static size_type d0 = 0;
1168  static bgeot::pstored_point_tab pspt1 = 0;
1169  mesh &mesh1
1171  if (d0 != d) {
1172  d0 = d;
1173  Bank_build_first_mesh(mesh1, d);
1174  pspt1 = bgeot::store_point_tab(mesh1.points());
1175  }
1176 
1177  const std::vector<size_type> &loc_ind = pgt->vertices();
1178  const bgeot::mesh_structure::ind_cv_ct &ct = ind_points_of_convex(ic);
1179 
1180  gs.ipt_loc.resize(ipt.size());
1181  std::vector<size_type> ipt_other;
1182  size_type nb_found = 0;
1183  for (size_type ip = 0; ip < loc_ind.size(); ++ip) {
1184  bool found = false;
1185  for (size_type i = 0; i < ipt.size(); ++i)
1186  if (ct[loc_ind[ip]] == ipt[i])
1187  { gs.ipt_loc[i] = ip; found = true; ++nb_found; break; }
1188  if (!found) ipt_other.push_back(ip);
1189  }
1190  assert(nb_found == ipt.size());
1191 
1192  mesh mesh2;
1193  for (size_type ip = 0; ip < loc_ind.size(); ++ip)
1194  mesh2.add_point(pgt->geometric_nodes()[loc_ind[ip]]);
1195  size_type ic1 = mesh2.add_simplex(dim_type(d), gs.ipt_loc.begin());
1196  bgeot::pgeometric_trans pgt1 = mesh2.trans_of_convex(ic1);
1197  for (size_type i = 0; i < ipt.size(); ++i)
1198  gs.ipt_loc[i] = loc_ind[gs.ipt_loc[i]];
1199 
1200  bgeot::pgeotrans_precomp pgp = bgeot::geotrans_precomp(pgt1, pspt1, 0);
1201 
1202  std::vector<size_type> ipt1(pspt1->size());
1203  base_node pt(n);
1204  for (size_type i = 0; i < pspt1->size(); ++i) {
1205  pgp->transform(mesh2.points_of_convex(ic1), i, pt);
1206  ipt1[i] = mesh2.add_point(pt);
1207  }
1208  mesh2.sup_convex(ic1);
1209 
1210  std::vector<size_type> ipt2(n+1);
1211  for (size_type i = 0; i < mesh1.nb_convex(); ++i) {
1212  for (size_type j = 0; j <= d; ++j)
1213  ipt2[j] = ipt1[mesh1.ind_points_of_convex(i)[j]];
1214  for (size_type j = d+1; j <= n; ++j)
1215  ipt2[j] = ipt_other[j-d-1];
1216  mesh2.add_simplex(dim_type(n), ipt2.begin());
1217  }
1218 
1219  mesh mesh3;
1220  ipt1.resize(pgt->nb_points());
1221  for (dal::bv_visitor i(mesh2.convex_index()); !i.finished(); ++i) {
1222  bgeot::pgeometric_trans pgt2 = mesh2.trans_of_convex(i);
1223  for (size_type ip = 0; ip < pgt->nb_points(); ++ip)
1224  ipt1[ip] =
1225  mesh3.add_point(pgt2->transform(pgt->geometric_nodes()[ip],
1226  mesh2.points_of_convex(i)));
1227  mesh3.add_convex(pgt, ipt1.begin());
1228  }
1229 
1230 
1231  bgeot::pstored_point_tab pspt3 = bgeot::store_point_tab(mesh3.points());
1232  pgp = bgeot::geotrans_precomp(pgt, pspt3, 0);
1233 
1234  ipt1.resize(pspt3->size());
1235  for (size_type ip = 0; ip < pspt3->size(); ++ip) {
1236  pgp->transform(points_of_convex(ic), ip, pt);
1237  ipt1[ip] = add_point(pt);
1238  }
1239  // dal::del_stored_object(pspt3);
1240 
1241  ipt2.resize(pgt->nb_points());
1242  for (size_type icc = 0; icc < mesh3.nb_convex(); ++icc) {
1243  for (size_type j = 0; j < pgt->nb_points(); ++j)
1244  ipt2[j] = ipt1[mesh3.ind_points_of_convex(icc)[j]];
1245  size_type i = add_convex(pgt, ipt2.begin());
1246  gs.sub_simplices.push_back(i);
1247  Bank_info->is_green_simplex.add(i);
1248  Bank_info->num_green_simplex[i] = igs;
1249  }
1250 
1251  for (size_type ip1 = 0; ip1 < ipt.size(); ++ip1)
1252  for (size_type ip2 = ip1+1; ip2 < ipt.size(); ++ip2)
1253  Bank_info->edges.insert(edge(ipt[ip1], ipt[ip2]));
1254 
1255  handle_region_refinement(ic, gs.sub_simplices, true);
1256  sup_convex(ic, true);
1257  }
1258 
1259  void mesh::Bank_refine(dal::bit_vector b) {
1260  if (!(Bank_info.get())) Bank_info = std::make_unique<Bank_info_struct>();
1261 
1262  b &= convex_index();
1263  if (b.card() == 0) return;
1264 
1265  Bank_info->edges.clear();
1266  while (b.card() > 0)
1267  Bank_test_and_refine_convex(b.take_first(), b);
1268 
1269  std::vector<size_type> ipt;
1270  edge_set marked_convexes;
1271  while (Bank_info->edges.size()) {
1272  marked_convexes.clear();
1273  b = convex_index();
1274  edge_set::const_iterator it = Bank_info->edges.begin();
1275  edge_set::const_iterator ite = Bank_info->edges.end(), it2=it;
1276  assert(!Bank_info->edges.empty());
1277  for (; it != ite; it = it2) {
1278  if (it2 != ite) ++it2;
1279  Bank_convex_with_edge(it->i1, it->i2, ipt);
1280  if (ipt.size() == 0) ; // Bank_info->edges.erase(it);
1281  else for (size_type ic = 0; ic < ipt.size(); ++ic)
1282  marked_convexes.insert(edge(ipt[ic], it->i1, it->i2));
1283  }
1284 
1285  it = marked_convexes.begin(); ite = marked_convexes.end();
1286  if (it == ite) break;
1287 
1288  while (it != ite) {
1289  size_type ic = it->i0;
1290  ipt.resize(0);
1291  while (it != ite && it->i0 == ic) {
1292  bool found1 = false, found2 = false;
1293  for (size_type j = 0; j < ipt.size(); ++j) {
1294  if (ipt[j] == it->i1) found1 = true;
1295  if (ipt[j] == it->i2) found2 = true;
1296  }
1297  if (!found1) ipt.push_back(it->i1);
1298  if (!found2) ipt.push_back(it->i2);
1299  ++it;
1300  }
1301  if (b.is_in(ic)) {
1302  if (ipt.size() > structure_of_convex(ic)->dim())
1303  Bank_test_and_refine_convex(ic, b);
1304  else if (Bank_info->is_green_simplex[ic]) {
1305  size_type icc = Bank_test_and_refine_convex(ic, b, false);
1306  if (!Bank_is_convex_having_points(icc, ipt)) {
1307  Bank_test_and_refine_convex(icc, b);
1308  }
1309  }
1310  else Bank_build_green_simplexes(ic, ipt);
1311  }
1312  }
1313  }
1314  Bank_info->edges.clear();
1315  }
1316 
1317  struct dummy_mesh_ {
1318  mesh m;
1319  dummy_mesh_() : m() {}
1320  };
1321 
1322  const mesh &dummy_mesh()
1324 
1325 } /* end of namespace getfem. */
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...
Definition: getfem_mesh.cc:878
base_matrix compute_local_basis(const geotrans_interpolation_context &c, size_type face)
return the local basis (i.e.
pconvex_ref equilateral_simplex_of_reference(dim_type nc)
equilateral simplex (degree 1).
void swap_convex(size_type cv1, size_type cv2)
Exchange two convex IDs.
indexed array reference (given a container X, and a set of indexes I, this class provides a pseudo-co...
Definition: gmm_ref.h:289
void write_to_file(const std::string &name) const
Write the mesh to a file.
Definition: getfem_mesh.cc:667
Define a getfem::getfem_mesh object.
short_type nb_points_of_convex(size_type ic) const
Return the number of points of convex ic.
static mesh_region all_convexes()
provide a default value for the mesh_region parameters of assembly procedures etc.
iterator over a gmm::tab_ref_index_ref<ITER,ITER_INDEX>
Definition: gmm_ref.h:226
ind_pt_face_ct ind_points_of_face_of_convex(size_type ic, short_type f) const
Return a container of the (global) point number for face f or convex ic.
scalar_type minimal_convex_radius_estimate() const
Return an estimate of the convex smallest dimension.
Definition: getfem_mesh.cc:426
void read_from_file(const std::string &name)
Load the mesh from a file.
Definition: getfem_mesh.cc:613
scalar_type convex_area_estimate(size_type ic, size_type degree=2) const
Return an estimate of the convex area.
Definition: getfem_mesh.cc:401
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:703
void optimize_structure()
Reorder the convex IDs and point IDs, such that there is no hole in their numbering.
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
computation of the condition number of dense matrices.
Describe a mesh (collection of convexes (elements) and points).
Definition: getfem_mesh.h:95
base_matrix local_basis_of_face_of_convex(size_type ic, short_type f, const base_node &pt) const
Return a local basis for the convex face.
Definition: getfem_mesh.cc:380
pconvex_structure structure_of_convex(size_type ic) const
Return the pconvex_structure of the convex ic.
size_type add_convex(bgeot::pgeometric_trans pgt, ITER ipts)
Add a convex to the mesh.
Definition: getfem_mesh.h:226
the geotrans_interpolation_context structure is passed as the argument of geometric transformation in...
does the inversion of the geometric transformation for a given convex
size_type add_triangle_by_points(const base_node &p1, const base_node &p2, const base_node &p3)
Add a triangle to the mesh, given the coordinates of its vertices.
Definition: getfem_mesh.cc:295
void clear()
Erase the mesh.
Definition: getfem_mesh.cc:273
Integration methods (exact and approximated) on convexes.
pgeometric_trans geometric_trans_descriptor(std::string name)
Get the geometric transformation from its string name.
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
void sup_convex_from_regions(size_type cv)
Remove all references to a convex from all regions stored in the mesh.
Definition: getfem_mesh.cc:42
static T & instance()
Instance from the current thread.
const mesh_region & from_mesh(const mesh &m) const
For regions which have been built with just a number &#39;id&#39;, from_mesh(m) sets the current region to &#39;m...
int get_token(std::istream &ist, std::string &st, bool ignore_cr, bool to_up, bool read_un_pm, int *linenb)
Very simple lexical analysis of general interest for reading small langages with a "MATLAB like" synt...
Definition: bgeot_ftool.cc:49
size_type nb_convex() const
The total number of convexes in the mesh.
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:49
void copy_from(const mesh &m)
Clone a mesh.
Definition: getfem_mesh.cc:446
scalar_type convex_quality_estimate(size_type ic) const
Return an estimate of the convex quality (0 <= Q <= 1).
Definition: getfem_mesh.cc:409
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.
Definition: getfem_mesh.cc:818
virtual scalar_type convex_radius_estimate(size_type ic) const
Return an estimate of the convex largest dimension.
Definition: getfem_mesh.cc:420
Dynamic Array.
Definition: dal_basic.h:47
A simple singleton implementation.
bool invert(const base_node &n, base_node &n_ref, scalar_type IN_EPS=1e-12)
given the node on the real element, returns the node on the reference element (even if it is outside ...
static mesh_region intersection(const mesh_region &a, const mesh_region &b)
return the intersection of two mesh regions
bool is_point_valid(size_type i) const
Return true if the point i is used by at least one convex.
void Bank_refine(dal::bit_vector)
Use the Bank strategy to refine some convexes.
"iterator" class for regions.
mesh(const std::string name="")
Constructor.
Definition: getfem_mesh.cc:133
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.
Definition: getfem_mesh.cc:923
size_type neighbour_of_convex(size_type ic, short_type f) const
Return a neighbour convex of a given convex face.
Deal with interdependencies of objects.
GEneric Tool for Finite Element Methods.
size_type add_segment(size_type a, size_type b)
Add a segment to the mesh, given the point id of its vertices.
Definition: getfem_mesh.cc:283
base_small_vector normal_of_face_of_convex(size_type ic, short_type f, const base_node &pt) const
Return the normal of the given convex face, evaluated at the point pt.
Definition: getfem_mesh.cc:343
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:756
gmm::uint16_type short_type
used as the common short type integer in the library
Definition: bgeot_config.h:79
std::string name_of_geometric_trans(pgeometric_trans p)
Get the string name of a geometric transformation.
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:730
size_type nb_allocated_convex() const
The number of convex indexes from 0 to the index of the last convex.
void swap_convex(size_type i, size_type j)
Swap the indexes of the convex of indexes i and j in the whole structure.
Definition: getfem_mesh.cc:332
base_small_vector compute_normal(const geotrans_interpolation_context &c, size_type face)
norm of returned vector is the ratio between the face surface on the real element and the face surfac...
pconvex_structure simplex_structure(dim_type nc)
Give a pointer on the structures of a simplex of dimension d.
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:864
void translation(const base_small_vector &)
Apply the given translation to each mesh node.
Definition: getfem_mesh.cc:252
void bounding_box(base_node &Pmin, base_node &Pmax) const
Return the bounding box [Pmin - Pmax] of the mesh.
Definition: getfem_mesh.cc:261
pintegration_method classical_approx_im(bgeot::pgeometric_trans pgt, dim_type degree)
try to find an approximate integration method for the geometric transformation pgt which is able to i...
pconvex_structure basic_structure(pconvex_structure cv)
Original structure (if concerned)
Mesh structure definition.
size_type add_triangle(size_type a, size_type b, size_type c)
Add a triangle to the mesh, given the point id of its vertices.
Definition: getfem_mesh.cc:288
Basic Geometric Tools.
void sup_convex(size_type ic, bool sup_points=false)
Delete the convex of index ic from the mesh.
Definition: getfem_mesh.cc:317
void cuthill_mckee_on_convexes(const bgeot::mesh_structure &ms, std::vector< size_type > &cmk)
Return the cuthill_mc_kee ordering on the convexes.
const ind_cv_ct & ind_points_of_convex(size_type ic) const
Return a container to the list of points attached to convex ic.
scalar_type maximal_convex_radius_estimate() const
Return an estimate of the convex largest dimension.
Definition: getfem_mesh.cc:435
void transformation(const base_matrix &)
apply the given matrix transformation to each mesh node.
Definition: getfem_mesh.cc:255
size_type add_pyramid(size_type a, size_type b, size_type c, size_type d, size_type e)
Add a pyramid to the mesh, given the point id of its vertices.
Definition: getfem_mesh.cc:304
void sup_convex(size_type ic)
Remove the convex ic.
iterator begin(void)
Iterator on the first element.
Definition: dal_basic.h:225
size_type add_tetrahedron_by_points(const base_node &p1, const base_node &p2, const base_node &p3, const base_node &p4)
Add a tetrahedron to the mesh, given the coordinates of its vertices.
Definition: getfem_mesh.cc:311
void clear(void)
Clear and desallocate all the elements.
Definition: dal_basic.h:298
const mesh_region region(size_type id) const
Return the region of index &#39;id&#39;.
Definition: getfem_mesh.h:414
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:551
void neighbours_of_convex(size_type ic, short_type f, ind_set &s) const
Return in s a list of neighbours of a given convex face.
pconvex_ref simplex_of_reference(dim_type nc, short_type K)
returns a simplex of reference of dimension nc and degree k
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
Provides mesh of torus.
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_tetrahedron(size_type a, size_type b, size_type c, size_type d)
Add a tetrahedron to the mesh, given the point id of its vertices.
Definition: getfem_mesh.cc:298
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:780