28 #if defined(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 defined(GETFEM_HAVE_METIS)
37 gmm::uint64_type act_counter(
void) {
38 static gmm::uint64_type c = gmm::uint64_type(1);
43 for (dal::bv_visitor i(valid_cvf_sets); !i.finished(); ++i)
44 cvf_sets[i].sup_all(c);
49 for (dal::bv_visitor i(valid_cvf_sets); !i.finished(); ++i)
54 void mesh::handle_region_refinement(
size_type ic,
55 const std::vector<size_type> &icv,
61 for (dal::bv_visitor ir(valid_cvf_sets); !ir.finished(); ++ir) {
62 mesh_region &r = cvf_sets[ir];
65 if (refine && r[ic].any()) {
67 for (
size_type jc = 0; jc < icv.size(); ++jc)
71 for (
short_type f = 0; f < pgt->structure()->nb_faces(); ++f) {
74 for (
size_type jc = 0; jc < icv.size(); ++jc) {
76 for (
short_type fsub = 0; fsub < pgtsub->structure()->nb_faces();
79 ind_set::const_iterator it = s.begin(), ite = s.end();
81 for (; it != ite; ++it)
82 if (std::find(icv.begin(), icv.end(), *it) != icv.end())
83 { found =
true;
break; }
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]));
90 giv.init(points_of_convex(ic), pgt);
91 giv.
invert(pt, barycentre);
94 if (gmm::abs(pgt->convex_ref()->is_in_face(f,barycentre)) < 0.001)
102 for (
size_type jc = 0; jc < icv.size(); ++jc)
103 if (!refine && r[icv[jc]].any()) {
108 for (
short_type fsub = 0; fsub < pgtsub->structure()->nb_faces();
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]));
115 giv.init(points_of_convex(ic), pgt);
116 giv.
invert(pt, barycentre);
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; }
126 void mesh::init(
void) {
127 #if GETFEM_PARA_LEVEL > 1
130 cuthill_mckee_uptodate =
false;
135 mesh::mesh(
const bgeot::basic_mesh &m,
const std::string name)
136 :
bgeot::basic_mesh(m), name_(name) { init(); }
138 mesh::mesh(
const mesh &m) : context_dependencies(),
139 std::enable_shared_from_this<
getfem::mesh>()
142 mesh &mesh::operator=(
const mesh &m) {
143 clear_dependencies();
149 #if GETFEM_PARA_LEVEL > 1
151 void mesh::compute_mpi_region()
const {
153 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
154 MPI_Comm_size(MPI_COMM_WORLD, &size);
158 mpi_region.from_mesh(*
this);
161 std::vector<int> xadj(ne+1), adjncy, numelt(ne), npart(ne);
164 double t_ref = MPI_Wtime();
168 for (dal::bv_visitor ic(
convex_index()); !ic.finished(); ++ic, ++j) {
173 for (dal::bv_visitor ic(
convex_index()); !ic.finished(); ++ic, ++j) {
176 for (ind_set::iterator it = s.begin();
177 it != s.end(); ++it) { adjncy.push_back(indelt[*it]); ++k; }
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]));
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]));
195 if (npart[i] == rank) mpi_region.add(numelt[i]);
198 cout <<
"Partition time "<< MPI_Wtime()-t_ref << endl;
201 valid_sub_regions.clear();
204 void mesh::compute_mpi_sub_region(
size_type n)
const {
205 if (valid_cvf_sets.is_in(n)) {
209 mpi_sub_region[n] = mesh_region();
210 valid_sub_regions.add(n);
213 void mesh::intersect_with_mpi_region(mesh_region &rg)
const {
215 rg = get_mpi_region();
216 }
else if (
int(rg.id()) >= 0) {
217 rg = get_mpi_sub_region(rg.id());
226 for (i = 0; i < j; i++)
227 if (!convex_tab.index_valid(i))
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;
236 if (with_renumbering) {
237 std::vector<size_type> cmk, iord(nbc), iordinv(nbc);
238 for (i = 0; i < nbc; ++i) iord[i] = iordinv[i] = i;
241 for (i = 0; i < nbc; ++i) {
245 std::swap(iord[i], iord[j]);
246 std::swap(iordinv[iord[i]], iordinv[iord[j]]);
253 { pts.translation(V); touch(); }
256 pts.transformation(M);
257 Bank_info = std::unique_ptr<Bank_info_struct>();
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]);
274 mesh_structure::clear();
276 gtab.clear(); trans_exists.clear();
277 cvf_sets.clear(); valid_cvf_sets.clear();
284 size_type ipt[2]; ipt[0] = a; ipt[1] = b;
285 return add_convex(bgeot::simplex_geotrans(1, 1), &(ipt[0]));
290 size_type ipt[3]; ipt[0] = a; ipt[1] = b; ipt[2] = c;
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)); }
300 size_type ipt[4]; ipt[0] = a; ipt[1] = b; ipt[2] = c; ipt[3] = d;
307 return add_convex(bgeot::pyramid_QK_geotrans(1), &(ipt[0]));
311 (
const base_node &p1,
const base_node &p2,
312 const base_node &p3,
const base_node &p4) {
314 add_point(p3), add_point(p4));
322 "Please call the optimize_structure() function before "
323 "merging elements from another mesh");
325 "The provided mesh region should only contain convexes");
328 const dal::bit_vector &convexes = (rg ==
size_type(-1))
331 for (dal::bv_visitor cv(convexes); !cv.finished(); ++cv) {
336 GMM_ASSERT1(nb == rct.size(),
"Internal error");
337 std::vector<size_type> ind(nb);
343 base_node pt = msource.points()[old_pid];
345 if (new_pid < next_pid && new_pid >= nbpt) {
347 new_pid = pts.add_node(pt, -1.);
348 GMM_ASSERT1(new_pid == next_pid,
"Internal error");
350 old2new[old_pid] = new_pid;
359 static std::vector<size_type> ipt;
362 ipt.assign(ct.begin(), ct.end());
367 trans_exists.sup(ic);
369 if (Bank_info.get()) Bank_sup_convex_from_green(ic);
376 trans_exists.swap(i, j);
378 swap_convex_in_regions(i, j);
379 if (Bank_info.get()) Bank_swap_convex(i,j);
380 cvs_v_num[i] = cvs_v_num[j] = act_counter(); touch();
385 const base_node &pt)
const {
387 base_matrix G(dim(),pgt->nb_points());
388 vectors_to_base_matrix(G, points_of_convex(ic));
397 bgeot::pgeotrans_precomp pgp
398 = bgeot::geotrans_precomp(pgt, pgt->pgeometric_nodes(), 0);
400 vectors_to_base_matrix(G, points_of_convex(ic));
402 c(pgp,pgt->structure()->ind_points_of_face(f)[n], G);
406 base_small_vector mesh::mean_normal_of_face_of_convex(
size_type ic,
409 base_matrix G; vectors_to_base_matrix(G,points_of_convex(ic));
410 base_small_vector ptmean(dim());
411 size_type nbpt = pgt->structure()->nb_points_of_face(f);
413 gmm::add(pgt->geometric_nodes()[pgt->structure()->ind_points_of_face(f)[i]], ptmean);
414 ptmean /= scalar_type(nbpt);
417 n /= gmm::vect_norm2(n);
422 const base_node &pt)
const {
424 base_matrix G(dim(),pgt->nb_points());
425 vectors_to_base_matrix(G,points_of_convex(ic));
433 bgeot::pgeotrans_precomp pgp
434 = bgeot::geotrans_precomp(pgt, pgt->pgeometric_nodes(), 0);
435 base_matrix G(dim(),pgt->nb_points());
436 vectors_to_base_matrix(G,points_of_convex(ic));
438 c(pgp,pgt->structure()->ind_points_of_face(f)[n], G);
444 bgeot::vectors_to_base_matrix(G, points_of_convex(ic));
452 bgeot::vectors_to_base_matrix(G, points_of_convex(ic));
453 auto pgt = trans_of_convex(ic);
455 pgt = pgt_torus->get_original_transformation();
456 G.resize(pgt->dim(), G.ncols());
463 bgeot::vectors_to_base_matrix(G, points_of_convex(ic));
470 for (dal::bv_visitor cv(
convex_index()); !cv.finished(); ++cv) {
479 for (dal::bv_visitor cv(
convex_index()); !cv.finished(); ++cv) {
485 void mesh::set_name(
const std::string& name){name_=name;}
490 bgeot::basic_mesh::operator=(m);
491 for (
const auto &kv : m.cvf_sets) {
492 if (kv.second.get_parent_mesh() != 0)
493 cvf_sets[kv.first].set_parent_mesh(
this);
494 cvf_sets[kv.first] = kv.second;
496 valid_cvf_sets = m.valid_cvf_sets;
498 gmm::uint64_type d = act_counter();
499 for (dal::bv_visitor i(
convex_index()); !i.finished(); ++i)
501 Bank_info = std::unique_ptr<Bank_info_struct>();
502 if (m.Bank_info.get())
503 Bank_info = std::make_unique<Bank_info_struct>(*(m.Bank_info));
506 struct mesh_convex_structure_loc {
508 std::vector<size_type> pts;
512 gmm::stream_standard_locale sl(ist);
516 bool te =
false, please_get =
true;
520 ist.seekg(0);ist.clear();
521 bgeot::read_until(ist,
"BEGIN POINTS LIST");
526 if (!bgeot::casecmp(tmp,
"END"))
528 else if (!bgeot::casecmp(tmp,
"POINT")) {
530 if (!bgeot::casecmp(tmp,
"COUNT")) {
535 GMM_ASSERT1(!npt.is_in(ip),
536 "Two points with the same index. loading aborted.");
539 while (isdigit(tmp[0]) || tmp[0] ==
'-' || tmp[0] ==
'+'
544 for (
size_type i = 0; i < d; i++) v[i] = tmpv[i];
547 GMM_ASSERT1(!npt.is_in(ipl),
"Two points [#" << ip <<
" and #"
548 << ipl <<
"] with the same coords "<< v
549 <<
". loading aborted.");
553 }
else if (tmp.size()) {
554 GMM_ASSERT1(
false,
"Syntax error in file, at token '" << tmp
555 <<
"', pos=" << std::streamoff(ist.tellg()));
556 }
else if (ist.eof()) {
557 GMM_ASSERT1(
false,
"Unexpected end of stream while reading mesh");
566 if (!bgeot::read_until(ist,
"BEGIN MESH STRUCTURE DESCRIPTION"))
567 GMM_ASSERT1(
false,
"This seems not to be a mesh file");
571 if (!bgeot::casecmp(tmp,
"END"))
573 else if (!bgeot::casecmp(tmp,
"CONVEX")) {
576 if (!bgeot::casecmp(tmp,
"COUNT")) {
579 ic = gmm::abs(atoi(tmp.c_str()));
580 GMM_ASSERT1(!ncv.is_in(ic),
581 "Negative or repeated index, loading aborted.");
587 while (!isspace(c)) { tmp.push_back(c); ist.get(c); }
593 cv[ic].cstruct = pgt;
594 cv[ic].pts.resize(nb);
597 cv[ic].pts[i] = gmm::abs(atoi(tmp.c_str()));
601 else if (tmp.size()) {
602 GMM_ASSERT1(
false,
"Syntax error reading a mesh file "
603 " at pos " << std::streamoff(ist.tellg())
604 <<
"(expecting 'CONVEX' or 'END', found '"<< tmp <<
"')");
605 }
else if (ist.eof()) {
606 GMM_ASSERT1(
false,
"Unexpected end of stream "
607 <<
"(missing BEGIN MESH/END MESH ?)");
610 ist >> bgeot::skip(
"MESH STRUCTURE DESCRIPTION");
612 for (dal::bv_visitor ic(ncv); !ic.finished(); ++ic) {
621 if (bgeot::casecmp(tmp,
"BEGIN")==0) {
623 if (bgeot::casecmp(tmp,
"BOUNDARY")==0 ||
624 bgeot::casecmp(tmp,
"REGION")==0) {
629 if (bgeot::casecmp(tmp,
"END")!=0) {
634 if (!bgeot::casecmp(tmp,
"END"))
break;
635 int f = atoi(tmp.c_str());
641 if (!bgeot::casecmp(tmp,
"END"))
break;
656 std::ifstream o(name.c_str());
657 GMM_ASSERT1(o,
"Mesh file '" << name <<
"' does not exist");
663 void write_tab_to_file_(std::ostream &ost,
const ITER& b_,
const ITER& e)
664 {
for (ITER b(b_) ; b != e; ++b) ost <<
" " << *b; }
667 static void write_convex_to_file_(
const mesh &ms,
670 for ( ; b != e; ++b) {
672 ost <<
" CONVEX " << i <<
" \'"
675 write_tab_to_file_(ost, ms.ind_points_of_convex(i).begin(),
676 ms.ind_points_of_convex(i).end() );
681 template<
class ITER>
static void write_point_to_file_(std::ostream &ost,
683 {
for ( ; b != e; ++b) ost <<
" " << *b; ost <<
'\n'; }
687 gmm::stream_standard_locale sl(ost);
688 ost <<
'\n' <<
"BEGIN POINTS LIST" <<
'\n' <<
'\n';
689 ost <<
" POINT COUNT " << points().index().last_true()+1 <<
'\n';
692 ost <<
" POINT " << i;
693 write_point_to_file_(ost, pts[i].begin(), pts[i].end());
695 ost <<
'\n' <<
"END POINTS LIST" <<
'\n' <<
'\n' <<
'\n';
697 ost <<
'\n' <<
"BEGIN MESH STRUCTURE DESCRIPTION" <<
'\n' <<
'\n';
698 ost <<
" CONVEX COUNT " <<
convex_index().last_true()+1 <<
'\n';
699 write_convex_to_file_(*
this, ost, convex_tab.tas_begin(),
700 convex_tab.tas_end());
701 ost <<
'\n' <<
"END MESH STRUCTURE DESCRIPTION" <<
'\n';
703 for (dal::bv_visitor bnum(valid_cvf_sets); !bnum.finished(); ++bnum) {
704 ost <<
"BEGIN REGION " << bnum <<
"\n" <<
region(bnum) <<
"\n"
705 <<
"END REGION " << bnum <<
"\n";
710 std::ofstream o(name.c_str());
711 GMM_ASSERT1(o,
"impossible to write to file '" << name <<
"'");
712 o <<
"% GETFEM MESH FILE " <<
'\n';
713 o <<
"% GETFEM VERSION " << GETFEM_VERSION <<
'\n' <<
'\n' <<
'\n';
720 + pts.memsize() + (pts.index().last_true()+1)*dim()*
sizeof(scalar_type)
721 +
sizeof(
mesh) + trans_exists.memsize() + gtab.memsize()
722 + valid_cvf_sets.card()*
sizeof(
mesh_region) + valid_cvf_sets.memsize();
725 struct equilateral_to_GT_PK_grad_aux :
public std::vector<base_matrix> {};
726 static const base_matrix &equilateral_to_GT_PK_grad(dim_type N) {
727 std::vector<base_matrix>
729 if (N > pbm.size()) pbm.resize(N);
730 if (pbm[N-1].empty()) {
733 base_matrix G(N,N+1);
734 vectors_to_base_matrix
737 (pgt, pgt->pgeometric_nodes(), 0)->grad(0), Gr);
746 const base_matrix& G,
747 pintegration_method pi) {
750 static bgeot::pgeotrans_precomp pgp = 0;
751 static pintegration_method pim_old = 0;
752 papprox_integration pai = get_approx_im_or_fail(pi);
753 if (pgt_old != pgt || pim_old != pi) {
756 pgp = bgeot::geotrans_precomp
757 (pgt, pai->pintegration_points(), pi);
760 for (
size_type i = 0; i < pai->nb_points_on_convex(); ++i) {
762 area += pai->coeff(i) * gic.
J();
773 const base_matrix& G) {
775 static bgeot::pgeotrans_precomp pgp = 0;
776 if (pgt_old != pgt) {
778 pgp=bgeot::geotrans_precomp(pgt, pgt->pgeometric_nodes(), 0);
781 size_type n = (pgt->is_linear()) ? 1 : pgt->nb_points();
783 dim_type N = dim_type(G.nrows()), P = pgt->structure()->dim();
786 gmm::mult(G, pgp->grad(ip), K);
792 gmm::mult(base_matrix(K),equilateral_to_GT_PK_grad(P),K);
793 q = std::max(q, gmm::condition_number(K));
799 const base_matrix& G) {
801 static bgeot::pgeotrans_precomp pgp = 0;
802 if (pgt_old != pgt) {
804 pgp=bgeot::geotrans_precomp(pgt, pgt->pgeometric_nodes(), 0);
807 size_type n = (pgt->is_linear()) ? 1 : pgt->nb_points();
809 base_matrix K(pgp->grad(0).ncols(),G.nrows());
811 gmm::mult(gmm::transposed(pgp->grad(ip)), gmm::transposed(G), K);
812 scalar_type emax, emin; gmm::condition_number(K,emax,emin);
813 q = std::max(q, emax);
815 return q * sqrt(scalar_type(N)) / scalar_type(2);
823 convex_face_ct& flist) {
824 for (dal::bv_visitor ic(cvlst); !ic.finished(); ++ic) {
828 flist.push_back(convex_face(ic,f));
832 flist.push_back(convex_face(ic,
short_type(-1)));
838 const mesh_region &cvlst,
839 mesh_region &flist) {
840 cvlst.error_if_not_convexes();
841 for (mr_visitor i(cvlst); !i.finished(); ++i) {
842 if (m.structure_of_convex(i.cv())->dim() == m.dim()) {
843 for (
short_type f = 0; f < m.structure_of_convex(i.cv())->nb_faces();
845 size_type cv2 = m.neighbor_of_convex(i.cv(), f);
846 if (cv2 ==
size_type(-1) || !cvlst.is_in(cv2)) {
851 else flist.add(i.cv());
861 mr.error_if_not_convexes();
880 mr.error_if_not_convexes();
881 dal::bit_vector visited;
882 bgeot::mesh_structure::ind_set neighbors;
887 bool neighbor_visited =
false;
890 for (
size_type j = 0; j < neighbors.size(); ++j)
891 if (visited.is_in(neighbors[j]))
892 { neighbor_visited =
true;
break; }
894 if (!neighbor_visited) {
897 if (cv2 !=
size_type(-1) && mr.is_in(cv2)) mrr.add(cv1,f);
905 if (!(visited.is_in(cv1))) {
910 for (
size_type j = 0; j < neighbors.size(); ++j) {
911 if (visited.is_in(neighbors[j])) { ok =
false;
break; }
912 if (mr.is_in(neighbors[j])) { ok =
true; }
914 if (ok) { mrr.add(cv1,f); }
924 const base_small_vector &V,
927 scalar_type threshold = gmm::vect_norm2(V)*cos(angle);
930 base_node un = m.mean_normal_of_face_of_convex(i.cv(), i.f());
931 if (gmm::vect_sp(V, un) - threshold >= -1E-8)
932 mrr.add(i.cv(), i.f());
938 const base_node &pt1,
const base_node &pt2) {
941 GMM_ASSERT1(pt1.size() == N && pt2.size() == N,
"Wrong dimensions");
949 it != pt.end(); ++it) {
951 if (m.points()[*it][j] < pt1[j] || m.points()[*it][j] > pt2[j])
952 { is_in =
false;
break; }
955 if (is_in) mrr.add(i.cv(), i.f());
961 const base_node ¢er, scalar_type radius) {
964 GMM_ASSERT1(center.size() == N,
"Wrong dimensions");
972 it != pt.end(); ++it) {
973 scalar_type checked_radius = scalar_type(0.0);
975 checked_radius += pow(m.points()[*it][j] - center[j], 2);
976 checked_radius = std::sqrt(checked_radius);
977 if (checked_radius > radius) { is_in =
false;
break; }
979 if (is_in) mrr.add(i.cv(), i.f());
984 mesh_region select_convexes_in_box(
const mesh &m,
const mesh_region &mr,
985 const base_node &pt1,
const base_node &pt2) {
988 GMM_ASSERT1(pt1.size() == N && pt2.size() == N,
"Wrong dimensions");
991 bgeot::mesh_structure::ind_cv_ct pt = m.ind_points_of_convex(i.cv());
994 for (bgeot::mesh_structure::ind_cv_ct::iterator it = pt.begin();
995 it != pt.end(); ++it) {
997 if (m.points()[*it][j] < pt1[j] || m.points()[*it][j] > pt2[j])
998 { is_in =
false;
break; }
1001 if (is_in) mrr.add(i.cv());
1007 dim_type dim = in.dim();
1008 base_node pt(dim+1);
1010 size_type nbpt = in.points().index().last()+1;
1011 GMM_ASSERT1(nbpt == in.points().index().card(),
1012 "please call the optimize_structure() method before using "
1013 "the mesh as a base for prismatic mesh");
1014 size_type nb_layers_total = nb_layers * degree;
1015 for (
const base_node &inpt : in.points()) {
1016 std::copy(inpt.begin(), inpt.end(), pt.begin());
1018 for (
size_type j = 0; j <= nb_layers_total;
1019 ++j, pt[dim] += scalar_type(1) / scalar_type(nb_layers_total))
1023 std::vector<size_type> tab;
1024 for (dal::bv_visitor cv(in.
convex_index()); !cv.finished(); ++cv) {
1026 tab.resize((degree+1)*nbp);
1027 for (
size_type j = 0; j < nb_layers; ++j) {
1032 bgeot::product_geotrans(in.trans_of_convex(cv),
1033 bgeot::simplex_geotrans(1,degree));
1052 bool mesh::edge::operator <(
const edge &e)
const {
1053 if (i0 < e.i0)
return true;
1054 if (i0 > e.i0)
return false;
1055 if (i1 < e.i1)
return true;
1056 if (i1 > e.i1)
return false;
1057 if (i2 < e.i2)
return true;
1061 void mesh::Bank_sup_convex_from_green(
size_type i) {
1062 if (Bank_info.get() && Bank_info->is_green_simplex.is_in(i)) {
1063 size_type igs = Bank_info->num_green_simplex[i];
1064 green_simplex &gs = Bank_info->green_simplices[igs];
1065 for (
size_type j = 0; j < gs.sub_simplices.size(); ++j) {
1066 Bank_info->num_green_simplex.erase(gs.sub_simplices[j]);
1067 Bank_info->is_green_simplex.sup(gs.sub_simplices[j]);
1069 Bank_info->green_simplices.sup(igs);
1074 if (Bank_info.get()) {
1075 Bank_info->is_green_simplex.swap(i, j);
1076 std::map<size_type, size_type>::iterator
1077 iti = Bank_info->num_green_simplex.find(i);
1078 std::map<size_type, size_type>::iterator
1079 ite = Bank_info->num_green_simplex.end();
1080 std::map<size_type, size_type>::iterator
1081 itj = Bank_info->num_green_simplex.find(j);
1084 { numi = iti->second; Bank_info->num_green_simplex.erase(i); }
1086 { numj = itj->second; Bank_info->num_green_simplex.erase(j); }
1088 Bank_info->num_green_simplex[j] = numi;
1089 green_simplex &gs = Bank_info->green_simplices[numi];
1090 for (
size_type k = 0; k < gs.sub_simplices.size(); ++k)
1091 if (gs.sub_simplices[k] == i) gs.sub_simplices[k] = j;
1092 else if (gs.sub_simplices[k] == j) gs.sub_simplices[k] = i;
1095 Bank_info->num_green_simplex[i] = numj;
1096 if (iti == ite || numi != numj) {
1097 green_simplex &gs = Bank_info->green_simplices[numj];
1098 for (
size_type k = 0; k < gs.sub_simplices.size(); ++k)
1099 if (gs.sub_simplices[k] == i) gs.sub_simplices[k] = j;
1100 else if (gs.sub_simplices[k] == j) gs.sub_simplices[k] = i;
1106 void mesh::Bank_build_first_mesh(mesh &m,
size_type n) {
1109 for (
size_type ip = 0; ip < pcr->nb_points(); ++ip)
1110 m.add_point(pcr->points()[ip]);
1112 size_type nbc = bgeot::refinement_simplexe_tab(n, &tab);
1113 for (
size_type ic = 0; ic < nbc; ++ic, tab+=(n+1))
1114 m.add_simplex(dim_type(n), tab);
1117 struct mesh_cache_for_Bank_basic_refine_convex :
public mesh {};
1119 void mesh::Bank_basic_refine_convex(
size_type i) {
1121 size_type n = pgt->basic_structure()->dim();
1127 static bgeot::pstored_point_tab pspt = 0;
1128 static bgeot::pgeotrans_precomp pgp = 0;
1129 static std::vector<size_type> ipt, ipt2, icl;
1134 Bank_build_first_mesh(mesh1, n);
1137 ipt.resize(pgt->nb_points());
1138 for (
size_type ic = 0; ic < mesh1.nb_convex(); ++ic) {
1139 assert(mesh1.convex_index().is_in(ic));
1141 for (
size_type ip = 0; ip < pgt->nb_points(); ++ip)
1143 mesh2.add_point(pgt2->transform(pgt->geometric_nodes()[ip],
1144 mesh1.points_of_convex(ic)));
1145 mesh2.add_convex(pgt, &ipt[0]);
1149 pspt = bgeot::store_point_tab(mesh2.points());
1150 pgp = bgeot::geotrans_precomp(pgt, pspt, 0);
1153 base_node pt(dim());
1154 ipt.resize(pspt->size());
1155 for (
size_type ip = 0; ip < pspt->size(); ++ip) {
1156 pgp->transform(points_of_convex(i), ip, pt);
1157 ipt[ip] = add_point(pt);
1160 ipt2.resize(pgt->nb_points()); icl.resize(0);
1161 for (
size_type ic = 0; ic < mesh2.nb_convex(); ++ic) {
1162 for (
size_type j = 0; j < pgt->nb_points(); ++j)
1163 ipt2[j] = ipt[mesh2.ind_points_of_convex(ic)[j]];
1164 icl.push_back(
add_convex(pgt, ipt2.begin()));
1166 handle_region_refinement(i, icl,
true);
1171 std::vector<size_type> &ipt) {
1174 size_type cv = points_tab[i1][k], found = 0;
1175 const std::vector<size_type> &loc_ind = trans_of_convex(cv)->vertices();
1176 for (
size_type i = 0; i < loc_ind.size(); ++i) {
1177 size_type ipp = convex_tab[cv].pts[loc_ind[i]];
1178 if (ipp == i1) ++found;
1179 if (ipp == i2) ++found;
1181 GMM_ASSERT1(found <= 2,
"Invalid convex with repeated nodes ");
1182 if (found == 2) ipt.push_back(cv);
1186 bool mesh::Bank_is_convex_having_points(
size_type cv,
1187 const std::vector<size_type> &ipt) {
1189 const std::vector<size_type> &loc_ind = trans_of_convex(cv)->vertices();
1190 for (
size_type i = 0; i < loc_ind.size(); ++i)
1191 if (std::find(ipt.begin(), ipt.end(),
1192 convex_tab[cv].pts[loc_ind[i]]) != ipt.end()) ++found;
1193 return (found >= ipt.size());
1196 void mesh::Bank_refine_normal_convex(
size_type i) {
1199 "Sorry, refinement is only working with simplices.");
1201 const std::vector<size_type> &loc_ind = pgt->vertices();
1202 for (
size_type ip1 = 0; ip1 < loc_ind.size(); ++ip1)
1203 for (
size_type ip2 = ip1+1; ip2 < loc_ind.size(); ++ip2)
1206 Bank_basic_refine_convex(i);
1210 dal::bit_vector &b,
bool ref) {
1211 if (Bank_info->is_green_simplex[i]) {
1212 size_type igs = Bank_info->num_green_simplex[i];
1213 green_simplex &gs = Bank_info->green_simplices[igs];
1216 handle_region_refinement(icc, gs.sub_simplices,
false);
1217 for (
size_type ic = 0; ic < gs.sub_simplices.size(); ++ic) {
1219 b.sup(gs.sub_simplices[ic]);
1222 for (
size_type ip = 0; ip < gs.ipt_loc.size(); ++ip)
1223 for (
size_type jp = ip+1; jp < gs.ipt_loc.size(); ++jp)
1224 Bank_info->edges.insert
1227 Bank_sup_convex_from_green(i);
1228 if (ref) { Bank_refine_normal_convex(icc);
return size_type(-1); }
1231 else if (ref) Bank_refine_normal_convex(i);
1236 struct mesh_cache_for_Bank_build_green_simplexes :
public mesh {};
1238 void mesh::Bank_build_green_simplexes(
size_type ic,
1239 std::vector<size_type> &ipt) {
1240 GMM_ASSERT1(
convex_index().is_in(ic),
"Internal error");
1241 size_type igs = Bank_info->green_simplices.add(green_simplex());
1242 green_simplex &gs(Bank_info->green_simplices[igs]);
1244 ref_mesh_pt_ct ptab = points_of_convex(ic);
1245 pt_tab.assign(ptab.begin(), ptab.end());
1252 static bgeot::pstored_point_tab pspt1 = 0;
1257 Bank_build_first_mesh(mesh1, d);
1258 pspt1 = bgeot::store_point_tab(mesh1.points());
1261 const std::vector<size_type> &loc_ind = pgt->vertices();
1264 gs.ipt_loc.resize(ipt.size());
1265 std::vector<size_type> ipt_other;
1267 for (
size_type ip = 0; ip < loc_ind.size(); ++ip) {
1269 for (
size_type i = 0; i < ipt.size(); ++i)
1270 if (ct[loc_ind[ip]] == ipt[i])
1271 { gs.ipt_loc[i] = ip; found =
true; ++nb_found;
break; }
1272 if (!found) ipt_other.push_back(ip);
1274 assert(nb_found == ipt.size());
1277 for (
size_type ip = 0; ip < loc_ind.size(); ++ip)
1278 mesh2.add_point(pgt->geometric_nodes()[loc_ind[ip]]);
1279 size_type ic1 = mesh2.add_simplex(dim_type(d), gs.ipt_loc.begin());
1281 for (
size_type i = 0; i < ipt.size(); ++i)
1282 gs.ipt_loc[i] = loc_ind[gs.ipt_loc[i]];
1284 bgeot::pgeotrans_precomp pgp = bgeot::geotrans_precomp(pgt1, pspt1, 0);
1286 std::vector<size_type> ipt1(pspt1->size());
1287 base_node pt(dim());
1288 for (
size_type i = 0; i < pspt1->size(); ++i) {
1289 pgp->transform(mesh2.points_of_convex(ic1), i, pt);
1290 ipt1[i] = mesh2.add_point(pt);
1292 mesh2.sup_convex(ic1);
1294 std::vector<size_type> ipt2(n+1);
1295 for (
size_type i = 0; i < mesh1.nb_convex(); ++i) {
1297 ipt2[j] = ipt1[mesh1.ind_points_of_convex(i)[j]];
1299 ipt2[j] = ipt_other[j-d-1];
1300 mesh2.add_simplex(dim_type(n), ipt2.begin());
1304 ipt1.resize(pgt->nb_points());
1305 for (dal::bv_visitor i(mesh2.convex_index()); !i.finished(); ++i) {
1307 for (
size_type ip = 0; ip < pgt->nb_points(); ++ip)
1309 mesh3.add_point(pgt2->transform(pgt->geometric_nodes()[ip],
1310 mesh2.points_of_convex(i)));
1311 mesh3.add_convex(pgt, ipt1.begin());
1315 bgeot::pstored_point_tab pspt3 = bgeot::store_point_tab(mesh3.points());
1316 pgp = bgeot::geotrans_precomp(pgt, pspt3, 0);
1318 ipt1.resize(pspt3->size());
1319 for (
size_type ip = 0; ip < pspt3->size(); ++ip) {
1320 pgp->transform(points_of_convex(ic), ip, pt);
1321 ipt1[ip] = add_point(pt);
1325 ipt2.resize(pgt->nb_points());
1326 for (
size_type icc = 0; icc < mesh3.nb_convex(); ++icc) {
1327 for (
size_type j = 0; j < pgt->nb_points(); ++j)
1328 ipt2[j] = ipt1[mesh3.ind_points_of_convex(icc)[j]];
1330 gs.sub_simplices.push_back(i);
1331 Bank_info->is_green_simplex.add(i);
1332 Bank_info->num_green_simplex[i] = igs;
1335 for (
size_type ip1 = 0; ip1 < ipt.size(); ++ip1)
1336 for (
size_type ip2 = ip1+1; ip2 < ipt.size(); ++ip2)
1337 Bank_info->edges.insert(edge(ipt[ip1], ipt[ip2]));
1339 handle_region_refinement(ic, gs.sub_simplices,
true);
1344 if (!(Bank_info.get())) Bank_info = std::make_unique<Bank_info_struct>();
1347 if (b.card() == 0)
return;
1349 Bank_info->edges.clear();
1350 while (b.card() > 0)
1351 Bank_test_and_refine_convex(b.take_first(), b);
1353 std::vector<size_type> ipt;
1354 edge_set marked_convexes;
1355 while (Bank_info->edges.size()) {
1356 marked_convexes.clear();
1358 edge_set::const_iterator it = Bank_info->edges.begin();
1359 edge_set::const_iterator ite = Bank_info->edges.end(), it2=it;
1360 assert(!Bank_info->edges.empty());
1361 for (; it != ite; it = it2) {
1362 if (it2 != ite) ++it2;
1363 Bank_convex_with_edge(it->i1, it->i2, ipt);
1364 if (ipt.size() == 0) ;
1365 else for (
size_type ic = 0; ic < ipt.size(); ++ic)
1366 marked_convexes.insert(edge(ipt[ic], it->i1, it->i2));
1369 it = marked_convexes.begin(); ite = marked_convexes.end();
1370 if (it == ite)
break;
1375 while (it != ite && it->i0 == ic) {
1376 bool found1 =
false, found2 =
false;
1377 for (
size_type j = 0; j < ipt.size(); ++j) {
1378 if (ipt[j] == it->i1) found1 =
true;
1379 if (ipt[j] == it->i2) found2 =
true;
1381 if (!found1) ipt.push_back(it->i1);
1382 if (!found2) ipt.push_back(it->i2);
1387 Bank_test_and_refine_convex(ic, b);
1388 else if (Bank_info->is_green_simplex[ic]) {
1389 size_type icc = Bank_test_and_refine_convex(ic, b,
false);
1390 if (!Bank_is_convex_having_points(icc, ipt)) {
1391 Bank_test_and_refine_convex(icc, b);
1394 else Bank_build_green_simplexes(ic, ipt);
1398 Bank_info->edges.clear();
1401 struct dummy_mesh_ {
1403 dummy_mesh_() : m() {}
1406 const mesh &dummy_mesh()
the geotrans_interpolation_context structure is passed as the argument of geometric transformation in...
void set_ii(size_type ii__)
change the current point (assuming a geotrans_precomp_ is used)
scalar_type J() const
get the Jacobian of the geometric trans (taken at point xref() )
does the inversion of the geometric transformation for a given convex
bool invert(const base_node &n, base_node &n_ref, scalar_type IN_EPS=1e-12, bool project_into_element=false)
given the node on the real element, returns the node on the reference element (even if it is outside ...
Mesh structure definition.
short_type nb_points_of_convex(size_type ic) const
Return the number of points of convex ic.
const dal::bit_vector & convex_index() const
Return the list of valid convex IDs.
size_type neighbor_of_convex(size_type ic, short_type f) const
Return a neighbor convex of a given convex face.
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.
size_type nb_convex() const
The total number of convexes in the mesh.
void neighbors_of_convex(size_type ic, short_type f, ind_set &s) const
Return in s a list of neighbors of a given convex face.
void swap_convex(size_type cv1, size_type cv2)
Exchange two convex IDs.
pconvex_structure structure_of_convex(size_type ic) const
Return the pconvex_structure of the convex ic.
void optimize_structure()
Reorder the convex IDs and point IDs, such that there is no hole in their numbering.
const ind_set & ind_points_of_convex(size_type ic) const
Return a container to the list of points attached to convex ic.
bool is_point_valid(size_type i) const
Return true if the point i is used by at least one convex.
size_type nb_allocated_convex() const
The number of convex indexes from 0 to the index of the last convex.
void sup_convex(size_type ic)
Remove the convex ic.
size_type add_node(const base_node &pt, const scalar_type radius=0, bool remove_duplicated_nodes=true)
Add a point to the array or use an existing point, located within a distance smaller than radius.
iterator begin(void)
Iterator on the first element.
size_type size(void) const
Number of allocated elements.
void clear(void)
Clear and desallocate all the elements.
static T & instance()
Instance from the current thread.
"iterator" class for regions.
structure used to hold a set of convexes and/or convex faces.
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...
static mesh_region all_convexes()
provide a default value for the mesh_region parameters of assembly procedures etc.
static mesh_region intersection(const mesh_region &a, const mesh_region &b)
return the intersection of two mesh regions
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.
Describe a mesh (collection of convexes (elements) and points).
mesh(const std::string name="")
Constructor.
virtual scalar_type convex_radius_estimate(size_type ic) const
Return an estimate of the convex largest dimension.
void sup_convex_from_regions(size_type cv)
Remove all references to a convex from all regions stored in the mesh.
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.
void sup_convex(size_type ic, bool sup_points=false)
Delete the convex of index ic from the mesh.
const dal::bit_vector & points_index() const
Return the points index.
size_type nb_points() const
Give the number of geometrical nodes in the mesh.
void transformation(const base_matrix &)
apply the given matrix transformation to each mesh node.
scalar_type convex_quality_estimate(size_type ic) const
Return an estimate of the convex quality (0 <= Q <= 1).
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.
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.
scalar_type convex_area_estimate(size_type ic, size_type degree=2) const
Return an estimate of the convex area.
void write_to_file(const std::string &name) const
Write the mesh to a file.
void bounding_box(base_node &Pmin, base_node &Pmax) const
Return the bounding box [Pmin - Pmax] of the mesh.
void swap_convex(size_type i, size_type j)
Swap the indexes of the convex of indexes i and j in the whole structure.
void Bank_refine(dal::bit_vector)
Use the Bank strategy to refine some convexes.
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.
size_type add_segment(size_type a, size_type b)
Add a segment to the mesh, given the point id of its vertices.
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.
scalar_type minimal_convex_radius_estimate() const
Return an estimate of the convex smallest dimension.
const mesh_region region(size_type id) const
Return the region of index 'id'.
size_type add_convex(bgeot::pgeometric_trans pgt, ITER ipts)
Add a convex to the mesh.
void read_from_file(const std::string &name)
Load the mesh from a file.
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.
void sup_point(size_type i)
Delete the point of index i from the mesh if it is not linked to a convex.
void swap_points(size_type i, size_type j)
Swap the indexes of points of index i and j in the whole structure.
size_type add_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.
void copy_from(const mesh &m)
Clone a mesh.
size_type add_simplex(dim_type di, ITER ipts)
Add a simplex to the mesh, given its dimension and point numbers.
scalar_type maximal_convex_radius_estimate() const
Return an estimate of the convex largest dimension.
void merge_convexes_from_mesh(const mesh &m, size_type rg=size_type(-1), scalar_type tol=scalar_type(0))
Merge all convexes from a another mesh, possibly restricted to a mesh region.
void clear()
Erase the mesh.
void translation(const base_small_vector &)
Apply the given translation to each mesh node.
indexed array reference (given a container X, and a set of indexes I, this class provides a pseudo-co...
A simple singleton implementation.
Integration methods (exact and approximated) on convexes.
Define a getfem::getfem_mesh object.
void mult(const L1 &l1, const L2 &l2, L3 &l3)
*/
computation of the condition number of dense matrices.
void lu_inverse(const DenseMatrixLU &LU, const Pvector &pvector, const DenseMatrix &AInv_)
Given an LU factored matrix, build the inverse of the matrix.
scalar_type APIDECL convex_radius_estimate(bgeot::pgeometric_trans pgt, const base_matrix &pts)
rough estimate of the radius of the convex using the largest eigenvalue of the jacobian of the geomet...
scalar_type APIDECL convex_quality_estimate(bgeot::pgeometric_trans pgt, const base_matrix &pts)
rough estimate of the maximum value of the condition number of the jacobian of the geometric transfor...
mesh_region APIDECL select_faces_of_normal(const mesh &m, const mesh_region &mr, const base_small_vector &V, scalar_type angle)
Select in the region mr the faces of the mesh m with their unit outward vector having a maximal angle...
size_type add_convex_by_points(bgeot::pgeometric_trans pgt, ITER ipts, const scalar_type tol=scalar_type(0))
Add a convex to the mesh, given a geometric transformation and a list of point coordinates.
mesh_region APIDECL select_faces_in_ball(const mesh &m, const mesh_region &mr, const base_node ¢er, scalar_type radius)
Select in the region mr the faces of the mesh m lying entirely in the ball delimated by pt1 and radiu...
scalar_type APIDECL convex_area_estimate(bgeot::pgeometric_trans pgt, const base_matrix &pts, pintegration_method pim)
rough estimate of the convex area.
mesh_region APIDECL inner_faces_of_mesh(const mesh &m, const mesh_region &mr=mesh_region::all_convexes())
Select all the faces sharing at least two element of the given mesh region.
void APIDECL outer_faces_of_mesh(const mesh &m, const dal::bit_vector &cvlst, convex_face_ct &flist)
returns a list of "exterior" faces of a mesh (i.e.
mesh_region APIDECL select_faces_in_box(const mesh &m, const mesh_region &mr, const base_node &pt1, const base_node &pt2)
Select in the region mr the faces of the mesh m lying entirely in the box delimated by pt1 and pt2.
void APIDECL extrude(const mesh &in, mesh &out, size_type nb_layers, short_type degree=short_type(1))
build a N+1 dimensions mesh from a N-dimensions mesh by extrusion.
mesh_region APIDECL all_faces_of_mesh(const mesh &m, const mesh_region &mr=mesh_region::all_convexes())
Select all the faces of the given mesh region.
gmm::uint16_type short_type
used as the common short type integer in the library
std::string name_of_geometric_trans(pgeometric_trans p)
Get the string name of a geometric transformation.
pconvex_ref equilateral_simplex_of_reference(dim_type nc)
equilateral simplex (degree 1).
pconvex_ref simplex_of_reference(dim_type nc, short_type K)
returns a simplex of reference of dimension nc and degree k
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 languages with a "MATLAB like" syn...
void cuthill_mckee_on_convexes(const bgeot::mesh_structure &ms, std::vector< size_type > &cmk)
Return the cuthill_mc_kee ordering on the convexes.
base_matrix compute_local_basis(const geotrans_interpolation_context &c, size_type face)
return the local basis (i.e.
pconvex_structure simplex_structure(dim_type nc)
Give a pointer on the structures of a simplex of dimension d.
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...
pgeometric_trans geometric_trans_descriptor(std::string name)
Get the geometric transformation from its string name.
size_t size_type
used as the common size type in the library
pconvex_structure basic_structure(pconvex_structure cv)
Original structure (if concerned)
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
GEneric Tool for Finite Element Methods.
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...
An adaptor that adapts a two dimensional geometric_trans to include radial dimension.
iterator over a gmm::tab_ref_index_ref<ITER,ITER_INDEX>