GetFEM  5.5
getfem_mesh_region.cc
1 /*===========================================================================
2 
3  Copyright (C) 2005-2026 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, see https://www.gnu.org/licenses/.
18 
19 ===========================================================================*/
20 
22 #include "getfem/getfem_mesh.h"
23 #include "getfem/getfem_omp.h"
24 
25 namespace getfem {
26 
27  using face_bitset = mesh_region::face_bitset;
28 
29  mesh_region::mesh_region(const mesh_region &other)
30  : p(std::make_shared<impl>()), id_(size_type(-2)), parent_mesh(0) {
31  this->operator=(other);
32  }
33 
34  mesh_region::mesh_region()
35  : p(std::make_shared<impl>()), id_(size_type(-2)), type_(size_type(-1)),
36  partitioning_allowed{true}, parent_mesh(nullptr){
37  if (me_is_multithreaded_now()) prohibit_partitioning();
38  mark_region_changed();
39  }
40 
41  mesh_region::mesh_region(size_type id__) : id_(id__), type_(size_type(-1)),
42  partitioning_allowed{true}, parent_mesh(nullptr){
43  mark_region_changed();
44  }
45 
46  mesh_region::mesh_region(mesh& m, size_type id__, size_type type) :
47  p(std::make_shared<impl>()), id_(id__), type_(type),
48  partitioning_allowed{true}, parent_mesh(&m){
50  mark_region_changed();
51  }
52 
53  mesh_region::mesh_region(const dal::bit_vector &bv)
54  : p(std::make_shared<impl>()), id_(size_type(-2)), type_(size_type(-1)),
55  partitioning_allowed{true}, parent_mesh(nullptr){
57  add(bv);
58  mark_region_changed();
59  }
60 
61  void mesh_region::mark_region_changed() const{
62  index_updated.all_threads() = false;
63  partitions_updated.all_threads() = false;
64  serial_index_updated = false;
65  }
66 
67  void mesh_region::touch_parent_mesh(){
68  if (parent_mesh) parent_mesh->touch_from_region(id_);
69  }
70 
71  const mesh_region& mesh_region::from_mesh(const mesh &m) const{
72  if (!p)
73  {
74  auto r = const_cast<mesh_region*>(this);
75  if (id_ == size_type(-1))
76  {
77  r->p = std::make_shared<impl>();
78  r->add(m.convex_index());
79  }
80  else if (id_ != size_type(-2))
81  {
82  *r = m.region(id_);
83  }
84  }
85  mark_region_changed();
86  return *this;
87  }
88 
89  mesh_region& mesh_region::operator=(const mesh_region &from){
90  if (!parent_mesh && !from.parent_mesh){
91  id_ = from.id_;
92  type_ = from.type_;
93  partitioning_allowed.store(from.partitioning_allowed.load());
94  if (from.p) {
95  if (!p) p = std::make_shared<impl>();
96  wp() = from.rp();
97  }
98  else p = nullptr;
99  }
100  else if (!parent_mesh) {
101  p = from.p;
102  id_ = from.id_;
103  type_ = from.type_;
104  parent_mesh = from.parent_mesh;
105  partitioning_allowed.store(from.partitioning_allowed.load());
106  }
107  else {
108  if (from.p){
109  wp() = from.rp();
110  type_= from.get_type();
111  partitioning_allowed.store(from.partitioning_allowed.load());
112  }
113  else if (from.id_ == size_type(-1)) {
114  clear();
115  add(parent_mesh->convex_index());
116  type_ = size_type(-1);
117  partitioning_allowed = true;
118  }
119  touch_parent_mesh();
120  }
121  mark_region_changed();
122  return *this;
123  }
124 
125  bool mesh_region::compare(const mesh &m1,
126  const mesh_region &mr,
127  const mesh &m2) const {
128  if (&m1 != &m2) return false;
129  if (!p && !mr.p) return (id_ == mr.id_);
130  from_mesh(m1);
131  mr.from_mesh(m2);
132  if (p && !(mr.p)) return false;
133  if (!p && mr.p) return false;
134  if (p)
135  if (p->m != mr.p->m) return false;
136  return true;
137  }
138 
139  face_bitset mesh_region::operator[](size_t cv) const{
140  auto it = rp().m.find(cv);
141  if (it != rp().m.end()) return it->second;
142  else return {};
143  }
144 
145  void mesh_region::update_partition_iterators() const{
146  if ((partitions_updated.thrd_cast() == true)) return;
147  itbegin = partition_begin();
148  itend = partition_end ();
149  partitions_updated = true;
150  }
151  mesh_region::const_iterator
152  mesh_region::partition_begin( ) const{
153  auto region_size = rp().m.size();
154  if (region_size < partitions_updated.num_threads()){
155  //for small regions: put the whole region into zero thread
156  if (partitions_updated.this_thread() == 0) return rp().m.begin();
157  else return rp().m.end();
158  }
159  auto partition_size = static_cast<size_type>
160  (std::ceil(static_cast<scalar_type>(region_size)/
161  static_cast<scalar_type >(partitions_updated.num_threads())));
162  auto index_begin = partition_size * partitions_updated.this_thread();
163  if (index_begin >= region_size ) return rp().m.end();
164 
165  auto it = rp().m.begin();
166  for (size_type i = 0; i != index_begin; ++i) ++it;
167  return it;
168  }
169 
170  mesh_region::const_iterator
171  mesh_region::partition_end( ) const{
172  auto region_size = rp().m.size();
173  if (region_size< partitions_updated.num_threads()) return rp().m.end();
174 
175  auto partition_size = static_cast<size_type>
176  (std::ceil(static_cast<scalar_type>(region_size)/
177  static_cast<scalar_type >(partitions_updated.num_threads())));
178  auto index_end = partition_size * (partitions_updated.this_thread() + 1);
179  if (index_end >= region_size ) return rp().m.end();
180 
181  auto it = rp().m.begin();
182  for (size_type i = 0; i < index_end && it != rp().m.end(); ++i) ++it;
183  return it;
184  }
185 
186  mesh_region::const_iterator mesh_region::begin() const{
187  GMM_ASSERT1(p != 0, "Internal error");
188  if (me_is_multithreaded_now() && partitioning_allowed){
189  update_partition_iterators();
190  return itbegin;
191  }
192  else return rp().m.begin();
193  }
194 
195  mesh_region::const_iterator mesh_region::end() const{
196  if (me_is_multithreaded_now() && partitioning_allowed){
197  update_partition_iterators();
198  return itend;
199  }
200  else return rp().m.end();
201  }
202 
204  partitioning_allowed = true;
205  }
206 
207  void mesh_region::bounding_box(base_node& Pmin, base_node& Pmax) const{
208  auto &mesh = *get_parent_mesh();
209  for (auto &&cv : dal::bv_iterable_c{index()}) {
210  for (const auto &pt : mesh.points_of_convex(cv)) {
211  for (auto j = 0; j != Pmin.size(); ++j){
212  Pmin[j] = std::min(Pmin[j], pt[j]);
213  Pmax[j] = std::max(Pmax[j], pt[j]);
214  }
215  }
216  }
217  }
218 
220  partitioning_allowed = false;
221  }
222 
223  bool mesh_region::is_partitioning_allowed() const{
224  return partitioning_allowed;
225  }
226 
227  void mesh_region::update_index() const{
228  auto& convex_index = me_is_multithreaded_now() ?
229  rp().index_.thrd_cast() : rp().serial_index_;
230  if (convex_index.card() != 0) convex_index.clear();
231  for (auto &&pair : *this){
232  if (pair.second.any()) convex_index.add(pair.first);
233  }
234  }
235 
236  const dal::bit_vector& mesh_region::index() const{
237  GMM_ASSERT1(p, "Use from_mesh on that region before");
238  if (me_is_multithreaded_now()) {
239  if (!(index_updated.thrd_cast() == true)){
240  update_index();
241  index_updated = true;
242  }
243  return rp().index_;
244  }
245  else {
246  if (!serial_index_updated){
247  update_index();
248  serial_index_updated = true;
249  }
250  return rp().serial_index_;
251  }
252  }
253 
254  void mesh_region::add(const dal::bit_vector &bv){
255  for (dal::bv_visitor i(bv); !i.finished(); ++i){
256  wp().m[i].set(0, 1);
257  }
258  touch_parent_mesh();
259  mark_region_changed();
260  }
261 
262  void mesh_region::add(size_type cv, short_type f){
263  wp().m[cv].set(short_type(f + 1), 1);
264  touch_parent_mesh();
265  mark_region_changed();
266  }
267 
268  void mesh_region::sup_all(size_type cv){
269  auto it = wp().m.find(cv);
270  if (it != wp().m.end()){
271  wp().m.erase(it);
272  touch_parent_mesh();
273  mark_region_changed();
274  }
275  }
276 
277  void mesh_region::sup(size_type cv, short_type f){
278  auto it = wp().m.find(cv);
279  if (it != wp().m.end()) {
280  it->second.set(short_type(f + 1), 0);
281  if (it->second.none()) wp().m.erase(it);
282  touch_parent_mesh();
283  mark_region_changed();
284  }
285  }
286 
287  void mesh_region::clear(){
288  wp().m.clear();
289  touch_parent_mesh();
290  mark_region_changed();
291  }
292 
293  void mesh_region::clean(){
294  for (map_t::iterator it = wp().m.begin(), itn;
295  it != wp().m.end(); it = itn) {
296  itn = it;
297  ++itn;
298  if (!(*it).second.any()) wp().m.erase(it);
299  }
300  touch_parent_mesh();
301  mark_region_changed();
302  }
303 
304  void mesh_region::swap_convex(size_type cv1, size_type cv2){
305  auto it1 = wp().m.find(cv1);
306  auto it2 = wp().m.find(cv2);
307  auto ite = wp().m.end();
308  face_bitset f1, f2;
309 
310  if (it1 != ite) f1 = it1->second;
311  if (it2 != ite) f2 = it2->second;
312  if (!f1.none()) wp().m[cv2] = f1;
313  else if (it2 != ite) wp().m.erase(it2);
314  if (!f2.none()) wp().m[cv1] = f2;
315  else if (it1 != ite) wp().m.erase(it1);
316  touch_parent_mesh();
317  mark_region_changed();
318  }
319 
320  bool mesh_region::is_in(size_type cv, short_type f) const{
321  GMM_ASSERT1(p, "Use from mesh on that region before");
322  auto it = rp().m.find(cv);
323  if (it == rp().m.end() || short_type(f+1) >= MAX_FACES_PER_CV) return false;
324  return ((*it).second)[short_type(f+1)];
325  }
326 
327  bool mesh_region::is_in(size_type cv, short_type f, const mesh &m) const{
328  if (p) {
329  auto it = rp().m.find(cv);
330  if (it == rp().m.end() || short_type(f+1) >= MAX_FACES_PER_CV)
331  return false;
332  return ((*it).second)[short_type(f+1)];
333  }
334  else{
335  if (id() == size_type(-1)) return true;
336  else return m.region(id()).is_in(cv, f);
337  }
338  }
339 
340  bool mesh_region::is_empty() const{
341  return rp().m.empty();
342  }
343 
345  return is_empty() ||
346  (or_mask()[0] == true && or_mask().count() == 1);
347  }
348 
350  return (id() != size_type(-1)) && (is_empty() || (and_mask()[0] == false));
351  }
352 
353  face_bitset mesh_region::faces_of_convex(size_type cv) const{
354  auto it = rp().m.find(cv);
355  if (it != rp().m.end()) return ((*it).second) >> 1;
356  else return face_bitset();
357  }
358 
359  face_bitset mesh_region::and_mask() const{
360  face_bitset bs;
361  if (rp().m.empty()) return bs;
362  bs.set();
363  for (auto it = rp().m.begin(); it != rp().m.end(); ++it)
364  if ( (*it).second.any() ) bs &= (*it).second;
365  return bs;
366  }
367 
368  face_bitset mesh_region::or_mask() const{
369  face_bitset bs;
370  if (rp().m.empty()) return bs;
371  for (auto it = rp().m.begin(); it != rp().m.end(); ++it)
372  if ( (*it).second.any() ) bs |= (*it).second;
373  return bs;
374  }
375 
377  size_type sz = 0;
378  for (auto it = begin(); it != end(); ++it)
379  sz += (*it).second.count();
380  return sz;
381  }
382 
383  size_type mesh_region::unpartitioned_size() const{
384  size_type sz = 0;
385  for (auto it = rp().m.begin(); it != rp().m.end(); ++it)
386  sz += (*it).second.count();
387  return sz;
388  }
389 
391  const mesh_region &b){
392  GMM_TRACE4("intersection of "<<a.id()<<" and "<<b.id());
393  mesh_region r;
394  /* we do not allow the "all_convexes" kind of regions
395  for these operations as there are not intended to be manipulated
396  (they only exist to provide a default argument to the mesh_region
397  parameters of assembly procedures etc. */
398  GMM_ASSERT1(a.id() != size_type(-1)||
399  b.id() != size_type(-1), "the 'all_convexes' regions "
400  "are not supported for set operations");
401  if (a.id() == size_type(-1)){
402  for (const_iterator it = b.begin();it != b.end(); ++it) r.wp().m.insert(*it);
403  return r;
404  }
405  else if (b.id() == size_type(-1)){
406  for (const_iterator it = a.begin();it != a.end(); ++it) r.wp().m.insert(*it);
407  return r;
408  }
409 
410  auto ita = a.begin(), enda = a.end(),
411  itb = b.begin(), endb = b.end();
412 
413  while (ita != enda && itb != endb) {
414  if (ita->first < itb->first) ++ita;
415  else if (ita->first > itb->first) ++itb;
416  else {
417  face_bitset maska = ita->second, maskb = itb->second, bs;
418  if (maska[0] && !maskb[0]) bs = maskb;
419  else if (maskb[0] && !maska[0]) bs = maska;
420  else bs = maska & maskb;
421  if (bs.any()) r.wp().m.insert(r.wp().m.end(), std::make_pair(ita->first,bs));
422  ++ita; ++itb;
423  }
424  }
425  return r;
426  }
427 
429  const mesh_region &b){
430  GMM_TRACE4("Merger of " << a.id() << " and " << b.id());
431  mesh_region r;
432  GMM_ASSERT1(a.id() != size_type(-1) &&
433  b.id() != size_type(-1), "the 'all_convexes' regions "
434  "are not supported for set operations");
435  for (auto it = a.begin();it != a.end(); ++it){
436  r.wp().m.insert(*it);
437  }
438  for (auto it = b.begin();it != b.end(); ++it){
439  r.wp().m[it->first] |= it->second;
440  }
441  return r;
442  }
443 
444 
446  const mesh_region &b){
447  GMM_TRACE4("subtraction of "<<a.id()<<" and "<<b.id());
448  mesh_region r;
449  GMM_ASSERT1(a.id() != size_type(-1) &&
450  b.id() != size_type(-1), "the 'all_convexes' regions "
451  "are not supported for set operations");
452  for (auto ita = a.begin();ita != a.end(); ++ita)
453  r.wp().m.insert(*ita);
454 
455  for (auto itb = b.begin();itb != b.end(); ++itb){
456  auto cv = itb->first;
457  auto it = r.wp().m.find(cv);
458  if (it != r.wp().m.end()){
459  it->second &= ~(itb->second);
460  if (it->second.none()) r.wp().m.erase(it);
461  }
462  }
463  return r;
464  }
465 
467  const mesh_region &rg2,
468  const getfem::mesh& m2) const{
469  if (&m1 != &m2) return 0;
470  int r = 1, partially = 0;
471  for (mr_visitor cv(*this, m1); !cv.finished(); cv.next())
472  if (cv.is_face() && rg2.is_in(cv.cv(),short_type(-1), m2))
473  partially = -1;
474  else
475  r = 0;
476  return r == 1 ? 1 : partially;
477  }
478 
480  return m.regions_index().last_true()+1;
481  }
482 
483  void mesh_region::error_if_not_faces() const{
484  GMM_ASSERT1(is_only_faces(), "Expecting a set of faces, not convexes");
485  }
486 
487  void mesh_region::error_if_not_convexes() const{
488  GMM_ASSERT1(is_only_convexes(), "Expecting a set of convexes, not faces");
489  }
490 
491  void mesh_region::error_if_not_homogeneous() const{
492  GMM_ASSERT1(is_only_faces() || is_only_convexes(), "Expecting a set "
493  "of convexes or a set of faces, but not a mixed set");
494  }
495 
496 
497 #if GETFEM_PARA_LEVEL > 1
498 
499  mesh_region::visitor::visitor(const mesh_region &s, const mesh &m,
500  bool intersect_with_mpi) :
501  cv_(size_type(-1)), f_(short_type(-1)), finished_(false)
502  {
503  if ((me_is_multithreaded_now() && s.partitioning_allowed)) {
504  s.from_mesh(m);
505  init(s);
506  } else {
507  if (s.id() == size_type(-1)) {
508  if (intersect_with_mpi)
509  init(m.get_mpi_region());
510  else
511  init(m.convex_index());
512  } else if (s.id() == size_type(-2) && s.p.get()) {
513  if (intersect_with_mpi) {
514  mpi_rg = std::make_unique<mesh_region>(s);
515  mpi_rg->from_mesh(m);
516  m.intersect_with_mpi_region(*mpi_rg);
517  init(*mpi_rg);
518  } else
519  init(s);
520  } else {
521  GMM_ASSERT1(s.id() != size_type(-2), "Internal error");
522  if (intersect_with_mpi)
523  init(m.get_mpi_sub_region(s.id()));
524  else
525  init(m.region(s.id()));
526  }
527  }
528  }
529 
530 #else
531 
532  mesh_region::visitor::visitor(const mesh_region &s, const mesh &m, bool)
533  :cv_(size_type(-1)), f_(short_type(-1)), finished_(false){
534  if ((me_is_multithreaded_now() && s.partitioning_allowed)) {
535  s.from_mesh(m);
536  init(s);
537  }
538  else {
539  if (s.id() == size_type(-1)) {
540  init(m.convex_index());
541  } else if (s.id() == size_type(-2) && s.p.get()) {
542  init(s);
543  } else {
544  GMM_ASSERT1(s.id() != size_type(-2), "Internal error");
545  init(m.region(s.id()));
546  }
547  }
548  }
549 
550 #endif
551 
552  bool mesh_region::visitor::next(){
553  if (whole_mesh) {
554  if (itb == iteb) {
555  finished_ = true;
556  return false;
557  }
558  cv_ = itb.index();
559  c = 0;
560  f_ = 0;
561  ++itb;
562  while (itb != iteb && !(*itb)) ++itb;
563  return true;
564  }
565  while (c.none()){
566  if (it == ite) { finished_=true; return false; }
567  cv_ = it->first;
568  c = it->second;
569  f_ = short_type(-1);
570  ++it;
571  if (c.none()) continue;
572  }
573  next_face();
574  return true;
575  }
576 
577  mesh_region::visitor::visitor(const mesh_region &s) :
578  cv_(size_type(-1)), f_(short_type(-1)), finished_(false){
579  init(s);
580  }
581 
582  void mesh_region::visitor::init(const dal::bit_vector &bv){
583  whole_mesh = true;
584  itb = bv.begin(); iteb = bv.end();
585  while (itb != iteb && !(*itb)) ++itb;
586  next();
587  }
588 
589  void mesh_region::visitor::init(const mesh_region &s){
590  whole_mesh = false;
591  it = s.begin();
592  ite = s.end();
593  next();
594  }
595 
596  std::ostream & operator <<(std::ostream &os, const mesh_region &w){
597  if (w.id() == size_type(-1))
598  os << " ALL_CONVEXES";
599  else if (w.p){
600  for (mr_visitor cv(w); !cv.finished(); cv.next()){
601  os << cv.cv();
602  if (cv.is_face()) os << "/" << cv.f();
603  os << " ";
604  }
605  }
606  else{
607  os << " region " << w.id();
608  }
609  return os;
610  }
611 
612  struct dummy_mesh_region_{
613  mesh_region mr;
614  };
615 
618  }
619 
620 } /* end of namespace getfem. */
const dal::bit_vector & convex_index() const
Return the list of valid convex IDs.
static T & instance()
Instance from the current thread.
"iterator" class for regions.
structure used to hold a set of convexes and/or convex faces.
static mesh_region merge(const mesh_region &a, const mesh_region &b)
return the union of two mesh_regions
const mesh_region & from_mesh(const mesh &m) const
For regions which have been built with just a number 'id', from_mesh(m) sets the current region to 'm...
void allow_partitioning()
In multithreaded part of the program makes only a partition of the region visible in the index() and ...
size_type size() const
Region size, or the size of the region partition on the current thread if the region is partitioned.
void prohibit_partitioning()
Disregard partitioning, which means being able to see the whole region in multithreaded code.
bool is_only_faces() const
Return true if the region do contain only convex faces.
int region_is_faces_of(const getfem::mesh &m1, const mesh_region &rg2, const getfem::mesh &m2) const
Test if the region is a boundary of a list of faces of elements of region rg.
static mesh_region intersection(const mesh_region &a, const mesh_region &b)
return the intersection of two mesh regions
void bounding_box(base_node &Pmin, base_node &Pmax) const
Return the bounding box [Pmin - Pmax] of the mesh_region.
bool is_only_convexes() const
Return true if the region do not contain any convex face.
const dal::bit_vector & index() const
Index of the region convexes, or the convexes from the partition on the current thread.
static size_type free_region_id(const getfem::mesh &m)
Extract the next region number that does not yet exists in the mesh.
static mesh_region subtract(const mesh_region &a, const mesh_region &b)
subtract the second region from the first one
Describe a mesh (collection of convexes (elements) and points).
Definition: getfem_mesh.h:98
const mesh_region region(size_type id) const
Return the region of index 'id'.
Definition: getfem_mesh.h:420
Define a getfem::getfem_mesh object.
region objects (set of convexes and/or convex faces)
Tools for multithreaded, OpenMP and Boost based parallelization.
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
GEneric Tool for Finite Element Methods.
bool me_is_multithreaded_now()
is the program running in the parallel section
Definition: getfem_omp.cc:105
const mesh_region & dummy_mesh_region()
Dummy mesh_region for default parameter of functions.