16 #include <deal.II/base/quadrature_lib.h> 17 #include <deal.II/base/thread_management.h> 18 #include <deal.II/base/mpi.h> 19 #include <deal.II/base/mpi.templates.h> 21 #include <deal.II/dofs/dof_handler.h> 22 #include <deal.II/dofs/dof_accessor.h> 23 #include <deal.II/dofs/dof_tools.h> 25 #include <deal.II/distributed/tria.h> 26 #include <deal.II/distributed/shared_tria.h> 28 #include <deal.II/fe/fe_nothing.h> 29 #include <deal.II/fe/fe_values.h> 30 #include <deal.II/fe/mapping_q1.h> 31 #include <deal.II/fe/mapping_q.h> 32 #include <deal.II/fe/mapping_q_generic.h> 33 #include <deal.II/fe/fe_q.h> 35 #include <deal.II/grid/filtered_iterator.h> 36 #include <deal.II/grid/tria.h> 37 #include <deal.II/grid/tria_accessor.h> 38 #include <deal.II/grid/tria_iterator.h> 39 #include <deal.II/grid/grid_tools.h> 40 #include <deal.II/grid/grid_tools_cache.h> 41 #include <deal.II/grid/grid_reordering.h> 42 #include <deal.II/grid/manifold.h> 44 #include <deal.II/lac/dynamic_sparsity_pattern.h> 45 #include <deal.II/lac/filtered_matrix.h> 46 #include <deal.II/lac/precondition.h> 47 #include <deal.II/lac/solver_cg.h> 48 #include <deal.II/lac/sparse_matrix.h> 49 #include <deal.II/lac/sparsity_pattern.h> 50 #include <deal.II/lac/sparsity_tools.h> 51 #include <deal.II/lac/vector.h> 52 #include <deal.II/lac/vector_memory.h> 54 #include <deal.II/numerics/matrix_tools.h> 56 #include <boost/random/uniform_real_distribution.hpp> 57 #include <boost/random/mersenne_twister.hpp> 65 #include <unordered_map> 68 DEAL_II_NAMESPACE_OPEN
74 template <
int dim,
int spacedim>
84 #if defined(DEAL_II_WITH_P4EST) && defined(DEBUG) 97 const std::vector<Point<spacedim> > &vertices = tria.
get_vertices ();
98 std::vector<bool> boundary_vertices (vertices.size(),
false);
104 for (; cell!=endc; ++cell)
105 for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
106 if (cell->face(face)->at_boundary ())
107 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_face; ++i)
108 boundary_vertices[cell->face(face)->vertex_index(i)] =
true;
112 double max_distance_sqr = 0;
113 std::vector<bool>::const_iterator pi = boundary_vertices.begin();
114 const unsigned int N = boundary_vertices.size();
115 for (
unsigned int i=0; i<N; ++i, ++pi)
117 std::vector<bool>::const_iterator pj = pi+1;
118 for (
unsigned int j=i+1; j<N; ++j, ++pj)
119 if ((*pi==
true) && (*pj==
true) &&
120 ((vertices[i]-vertices[j]).norm_square() > max_distance_sqr))
121 max_distance_sqr = (vertices[i]-vertices[j]).norm_square();
124 return std::sqrt(max_distance_sqr);
129 template <
int dim,
int spacedim>
135 unsigned int mapping_degree = 1;
137 mapping_degree = p->get_degree();
139 mapping_degree = p->get_degree();
142 const QGauss<dim> quadrature_formula (mapping_degree + 1);
143 const unsigned int n_q_points = quadrature_formula.
size();
153 cell = triangulation.begin_active(),
154 endc = triangulation.end();
156 double local_volume = 0;
159 for (; cell!=endc; ++cell)
160 if (cell->is_locally_owned())
163 for (
unsigned int q=0; q<n_q_points; ++q)
164 local_volume += fe_values.
JxW(q);
167 double global_volume = 0;
169 #ifdef DEAL_II_WITH_MPI 175 global_volume = local_volume;
177 return global_volume;
185 (
const std::vector<Point<1> > &all_vertices,
188 return all_vertices[vertex_indices[1]][0]
189 - all_vertices[vertex_indices[0]][0];
197 (
const std::vector<Point<3> > &all_vertices,
206 const double x[8] = { all_vertices[vertex_indices[0]](0),
207 all_vertices[vertex_indices[1]](0),
208 all_vertices[vertex_indices[2]](0),
209 all_vertices[vertex_indices[3]](0),
210 all_vertices[vertex_indices[4]](0),
211 all_vertices[vertex_indices[5]](0),
212 all_vertices[vertex_indices[6]](0),
213 all_vertices[vertex_indices[7]](0)
215 const double y[8] = { all_vertices[vertex_indices[0]](1),
216 all_vertices[vertex_indices[1]](1),
217 all_vertices[vertex_indices[2]](1),
218 all_vertices[vertex_indices[3]](1),
219 all_vertices[vertex_indices[4]](1),
220 all_vertices[vertex_indices[5]](1),
221 all_vertices[vertex_indices[6]](1),
222 all_vertices[vertex_indices[7]](1)
224 const double z[8] = { all_vertices[vertex_indices[0]](2),
225 all_vertices[vertex_indices[1]](2),
226 all_vertices[vertex_indices[2]](2),
227 all_vertices[vertex_indices[3]](2),
228 all_vertices[vertex_indices[4]](2),
229 all_vertices[vertex_indices[5]](2),
230 all_vertices[vertex_indices[6]](2),
231 all_vertices[vertex_indices[7]](2)
270 const double t3 = y[3]*x[2];
271 const double t5 = z[1]*x[5];
272 const double t9 = z[3]*x[2];
273 const double t11 = x[1]*y[0];
274 const double t14 = x[4]*y[0];
275 const double t18 = x[5]*y[7];
276 const double t20 = y[1]*x[3];
277 const double t22 = y[5]*x[4];
278 const double t26 = z[7]*x[6];
279 const double t28 = x[0]*y[4];
280 const double t34 = z[3]*x[1]*y[2]+t3*z[1]-t5*y[7]+y[7]*x[4]*z[6]+t9*y[6]-t11*z[4]-t5*y[3]-t14*z[2]+z[1]*x[4]*y[0]-t18*z[3]+t20*z[0]-t22*z[0]-y[0]*x[5]*z[4]-t26*y[3]+t28*z[2]-t9*y[1]-y[1]*x[4]*z[0]-t11*z[5];
281 const double t37 = y[1]*x[0];
282 const double t44 = x[1]*y[5];
283 const double t46 = z[1]*x[0];
284 const double t49 = x[0]*y[2];
285 const double t52 = y[5]*x[7];
286 const double t54 = x[3]*y[7];
287 const double t56 = x[2]*z[0];
288 const double t58 = x[3]*y[2];
289 const double t64 = -x[6]*y[4]*z[2]-t37*z[2]+t18*z[6]-x[3]*y[6]*z[2]+t11*z[2]+t5*y[0]+t44*z[4]-t46*y[4]-t20*z[7]-t49*z[6]-t22*z[1]+t52*z[3]-t54*z[2]-t56*y[4]-t58*z[0]+y[1]*x[2]*z[0]+t9*y[7]+t37*z[4];
290 const double t66 = x[1]*y[7];
291 const double t68 = y[0]*x[6];
292 const double t70 = x[7]*y[6];
293 const double t73 = z[5]*x[4];
294 const double t76 = x[6]*y[7];
295 const double t90 = x[4]*z[0];
296 const double t92 = x[1]*y[3];
297 const double t95 = -t66*z[3]-t68*z[2]-t70*z[2]+t26*y[5]-t73*y[6]-t14*z[6]+t76*z[2]-t3*z[6]+x[6]*y[2]*z[4]-z[3]*x[6]*y[2]+t26*y[4]-t44*z[3]-x[1]*y[2]*z[0]+x[5]*y[6]*z[4]+t54*z[5]+t90*y[2]-t92*z[2]+t46*y[2];
298 const double t102 = x[2]*y[0];
299 const double t107 = y[3]*x[7];
300 const double t114 = x[0]*y[6];
301 const double t125 = y[0]*x[3]*z[2]-z[7]*x[5]*y[6]-x[2]*y[6]*z[4]+t102*z[6]-t52*z[6]+x[2]*y[4]*z[6]-t107*z[5]-t54*z[6]+t58*z[6]-x[7]*y[4]*z[6]+t37*z[5]-t114*z[4]+t102*z[4]-z[1]*x[2]*y[0]+t28*z[6]-y[5]*x[6]*z[4]-z[5]*x[1]*y[4]-t73*y[7];
302 const double t129 = z[0]*x[6];
303 const double t133 = y[1]*x[7];
304 const double t145 = y[1]*x[5];
305 const double t156 = t90*y[6]-t129*y[4]+z[7]*x[2]*y[6]-t133*z[5]+x[5]*y[3]*z[7]-t26*y[2]-t70*z[3]+t46*y[3]+z[5]*x[7]*y[4]+z[7]*x[3]*y[6]-t49*z[4]+t145*z[7]-x[2]*y[7]*z[6]+t70*z[5]+t66*z[5]-z[7]*x[4]*y[6]+t18*z[4]+x[1]*y[4]*z[0];
306 const double t160 = x[5]*y[4];
307 const double t165 = z[1]*x[7];
308 const double t178 = z[1]*x[3];
309 const double t181 = t107*z[6]+t22*z[7]+t76*z[3]+t160*z[1]-x[4]*y[2]*z[6]+t70*z[4]+t165*y[5]+x[7]*y[2]*z[6]-t76*z[5]-t76*z[4]+t133*z[3]-t58*z[1]+y[5]*x[0]*z[4]+t114*z[2]-t3*z[7]+t20*z[2]+t178*y[7]+t129*y[2];
310 const double t207 = t92*z[7]+t22*z[6]+z[3]*x[0]*y[2]-x[0]*y[3]*z[2]-z[3]*x[7]*y[2]-t165*y[3]-t9*y[0]+t58*z[7]+y[3]*x[6]*z[2]+t107*z[2]+t73*y[0]-x[3]*y[5]*z[7]+t3*z[0]-t56*y[6]-z[5]*x[0]*y[4]+t73*y[1]-t160*z[6]+t160*z[0];
311 const double t228 = -t44*z[7]+z[5]*x[6]*y[4]-t52*z[4]-t145*z[4]+t68*z[4]+t92*z[5]-t92*z[0]+t11*z[3]+t44*z[0]+t178*y[5]-t46*y[5]-t178*y[0]-t145*z[0]-t20*z[5]-t37*z[3]-t160*z[7]+t145*z[3]+x[4]*y[6]*z[2];
313 return (t34+t64+t95+t125+t156+t181+t207+t228)/12.;
321 (
const std::vector<Point<2> > &all_vertices,
360 const double x[4] = { all_vertices[vertex_indices[0]](0),
361 all_vertices[vertex_indices[1]](0),
362 all_vertices[vertex_indices[2]](0),
363 all_vertices[vertex_indices[3]](0)
366 const double y[4] = { all_vertices[vertex_indices[0]](1),
367 all_vertices[vertex_indices[1]](1),
368 all_vertices[vertex_indices[2]](1),
369 all_vertices[vertex_indices[3]](1)
372 return (-x[1]*y[0]+x[1]*y[3]+y[0]*x[2]+x[0]*y[1]-x[0]*y[2]-y[1]*x[3]-x[2]*y[3]+x[3]*y[2])/2;
378 template <
int dim,
int spacedim>
383 const auto predicate = [](
const iterator &)
393 template <
int dim,
int spacedim>
400 ExcMessage(
"Invalid SubCellData supplied according to ::check_consistency(). " 401 "This is caused by data containing objects for the wrong dimension."));
404 std::vector<bool> vertex_used (vertices.size(),
false);
405 for (
unsigned int c=0; c<cells.size(); ++c)
406 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
408 Assert(cells[c].vertices[v] < vertices.size(),
409 ExcMessage(
"Invalid vertex index encountered! cells[" 415 +
" is invalid, because only " 417 +
" vertices were supplied."));
418 vertex_used[cells[c].vertices[v]] =
true;
425 std::vector<unsigned int> new_vertex_numbers (vertices.size(), invalid_vertex);
426 unsigned int next_free_number = 0;
427 for (
unsigned int i=0; i<vertices.size(); ++i)
428 if (vertex_used[i] ==
true)
430 new_vertex_numbers[i] = next_free_number;
435 for (
unsigned int c=0; c<cells.size(); ++c)
436 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
437 cells[c].vertices[v] = new_vertex_numbers[cells[c].vertices[v]];
441 for (
unsigned int v=0; v<GeometryInfo<1>::vertices_per_cell; ++v)
444 ExcMessage(
"Invalid vertex index in subcelldata.boundary_lines. " 445 "subcelldata.boundary_lines[" 451 +
" is invalid, because only " 453 +
" vertices were supplied."));
459 for (
unsigned int v=0; v<GeometryInfo<2>::vertices_per_cell; ++v)
462 ExcMessage(
"Invalid vertex index in subcelldata.boundary_quads. " 463 "subcelldata.boundary_quads[" 469 +
" is invalid, because only " 471 +
" vertices were supplied."));
479 std::vector<Point<spacedim> > tmp;
480 tmp.reserve (std::count(vertex_used.begin(), vertex_used.end(),
true));
481 for (
unsigned int v=0; v<vertices.size(); ++v)
482 if (vertex_used[v] ==
true)
483 tmp.push_back (vertices[v]);
484 swap (vertices, tmp);
489 template <
int dim,
int spacedim>
494 std::vector<unsigned int> &considered_vertices,
500 std::vector<unsigned int> new_vertex_numbers(vertices.size());
501 for (
unsigned int i=0; i<vertices.size(); ++i)
502 new_vertex_numbers[i] = i;
506 if (considered_vertices.size()==0)
507 considered_vertices = new_vertex_numbers;
509 Assert(considered_vertices.size() <= vertices.size(),
516 for (
unsigned int i=0; i<considered_vertices.size(); ++i)
518 Assert(considered_vertices[i]<vertices.size(),
520 if (new_vertex_numbers[considered_vertices[i]]!=considered_vertices[i])
531 for (
unsigned int j=i+1; j<considered_vertices.size(); ++j)
534 for (
unsigned int d=0; d<spacedim; ++d)
535 equal &= (fabs(vertices[considered_vertices[j]](d)-vertices[considered_vertices[i]](d))<tol);
538 new_vertex_numbers[considered_vertices[j]]=considered_vertices[i];
553 for (
unsigned int c=0; c<cells.size(); ++c)
554 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
555 cells[c].vertices[v]=new_vertex_numbers[cells[c].vertices[v]];
565 template <
int spacedim>
588 explicit Rotate2d (
const double angle)
594 return Point<2> (std::cos(angle)*p(0) - std::sin(angle) * p(1),
595 std::sin(angle)*p(0) + std::cos(angle) * p(1));
605 Rotate3d (
const double angle,
606 const unsigned int axis)
616 std::cos(angle)*p(1) - std::sin(angle) * p(2),
617 std::sin(angle)*p(1) + std::cos(angle) * p(2));
619 return Point<3> (std::cos(angle)*p(0) + std::sin(angle) * p(2),
621 -std::sin(angle)*p(0) + std::cos(angle) * p(2));
623 return Point<3> (std::cos(angle)*p(0) - std::sin(angle) * p(1),
624 std::sin(angle)*p(0) + std::cos(angle) * p(1),
629 const unsigned int axis;
632 template <
int spacedim>
636 explicit Scale (
const double factor)
650 template <
int dim,
int spacedim>
655 transform (Shift<spacedim>(shift_vector), triangulation);
664 transform (Rotate2d(angle), triangulation);
670 const unsigned int axis,
675 transform (Rotate3d(angle, axis), triangulation);
678 template <
int dim,
int spacedim>
684 transform (Scale<spacedim>(scaling_factor), triangulation);
696 const std::map<types::global_dof_index,double> &fixed_dofs,
699 const unsigned int n_dofs=S.
n();
711 SF.add_constraints(fixed_dofs);
712 SF.apply_constraints (f,
true);
713 solver.solve(SF, u, f, PF);
736 const bool solve_for_absolute_positions)
743 dof_handler.distribute_dofs(q1);
746 dof_handler.n_dofs ());
761 std::map<types::global_dof_index,double> fixed_dofs[dim];
762 typename std::map<unsigned int,Point<dim> >::const_iterator map_end=new_points.end();
766 endc=dof_handler.end();
767 for (; cell!=endc; ++cell)
771 for (
unsigned int vertex_no=0;
772 vertex_no<GeometryInfo<dim>::vertices_per_cell; ++vertex_no)
774 const unsigned int vertex_index=cell->vertex_index(vertex_no);
775 const Point<dim> &vertex_point=cell->vertex(vertex_no);
777 const typename std::map<unsigned int,Point<dim> >::const_iterator map_iter
778 = new_points.find(vertex_index);
780 if (map_iter!=map_end)
781 for (
unsigned int i=0; i<dim; ++i)
782 fixed_dofs[i].insert(std::pair<types::global_dof_index,double>
783 (cell->vertex_dof_index(vertex_no, 0),
784 (solve_for_absolute_positions ?
785 map_iter->second(i) :
786 map_iter->second(i) - vertex_point[i])
793 for (
unsigned int i=0; i<dim; ++i)
794 us[i].reinit (dof_handler.n_dofs());
798 for (
unsigned int i=0; i<dim; ++i)
800 S, fixed_dofs[i], us[i]);
805 std::vector<bool> vertex_touched (triangulation.n_vertices(),
false);
806 for (cell=dof_handler.begin_active(); cell!=endc; ++cell)
807 for (
unsigned int vertex_no=0;
808 vertex_no<GeometryInfo<dim>::vertices_per_cell; ++vertex_no)
809 if (vertex_touched[cell->vertex_index(vertex_no)] ==
false)
814 for (
unsigned int i=0; i<dim; ++i)
815 if (solve_for_absolute_positions)
816 v(i) = us[i](dof_index);
818 v(i) += us[i](dof_index);
820 vertex_touched[cell->vertex_index(vertex_no)] =
true;
824 template <
int dim,
int spacedim>
825 std::map<unsigned int, Point<spacedim> >
828 std::map<unsigned int, Point<spacedim> > vertex_map;
832 for (; cell!=endc; ++cell)
834 for (
unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
838 if (face->at_boundary())
840 for (
unsigned j = 0; j < GeometryInfo<dim>::vertices_per_face; ++j)
843 const unsigned int vertex_index = face->vertex_index(j);
844 vertex_map[vertex_index] = vertex;
856 template <
int dim,
int spacedim>
860 const bool keep_boundary)
877 double almost_infinite_length = 0;
879 cell=triangulation.begin(0); cell!=triangulation.end(0); ++cell)
880 almost_infinite_length += cell->diameter();
882 std::vector<double> minimal_length (triangulation.n_vertices(),
883 almost_infinite_length);
886 std::vector<bool> at_boundary (keep_boundary ? triangulation.n_vertices() : 0,
false);
889 const bool is_parallel_shared
892 cell=triangulation.begin_active(); cell!=triangulation.end(); ++cell)
893 if (is_parallel_shared || cell->is_locally_owned())
897 for (
unsigned int i=0; i<GeometryInfo<dim>::lines_per_cell; ++i)
899 const typename Triangulation<dim,spacedim>::line_iterator line
902 if (keep_boundary && line->at_boundary())
904 at_boundary[line->vertex_index(0)] =
true;
905 at_boundary[line->vertex_index(1)] =
true;
908 minimal_length[line->vertex_index(0)]
909 = std::min(line->diameter(),
910 minimal_length[line->vertex_index(0)]);
911 minimal_length[line->vertex_index(1)]
912 = std::min(line->diameter(),
913 minimal_length[line->vertex_index(1)]);
919 for (
unsigned int vertex=0; vertex<2; ++vertex)
920 if (cell->at_boundary(vertex) ==
true)
921 at_boundary[cell->vertex_index(vertex)] =
true;
923 minimal_length[cell->vertex_index(0)]
924 = std::min(cell->diameter(),
925 minimal_length[cell->vertex_index(0)]);
926 minimal_length[cell->vertex_index(1)]
927 = std::min(cell->diameter(),
928 minimal_length[cell->vertex_index(1)]);
937 boost::random::mt19937 rng;
938 boost::random::uniform_real_distribution<> uniform_distribution(-1,1);
946 std::vector<bool> vertex_moved (triangulation.n_vertices(),
false);
950 cell=triangulation.begin_active(); cell!=triangulation.end(); ++cell)
951 if (cell->is_locally_owned())
953 for (
unsigned int vertex_no=0; vertex_no<GeometryInfo<dim>::vertices_per_cell;
956 const unsigned global_vertex_no = cell->vertex_index(vertex_no);
961 if ((keep_boundary && at_boundary[global_vertex_no])
962 || vertex_moved[global_vertex_no]
963 || !locally_owned_vertices[global_vertex_no])
968 for (
unsigned int d=0; d<spacedim; ++d)
969 shift_vector(d) = uniform_distribution(rng);
971 shift_vector *= factor * minimal_length[global_vertex_no] /
972 std::sqrt(shift_vector.
square());
975 cell->vertex(vertex_no) += shift_vector;
976 vertex_moved[global_vertex_no] =
true;
980 #ifdef DEAL_II_WITH_P4EST 981 distributed_triangulation
982 ->communicate_locally_moved_vertices(locally_owned_vertices);
984 (void)distributed_triangulation;
995 const unsigned int n_vertices = triangulation.n_vertices();
996 std::vector<Point<spacedim> > new_vertex_locations (n_vertices);
997 const std::vector<Point<spacedim> > &old_vertex_locations
998 = triangulation.get_vertices();
1000 for (
unsigned int vertex=0; vertex<n_vertices; ++vertex)
1004 if (keep_boundary && at_boundary[vertex])
1005 new_vertex_locations[vertex] = old_vertex_locations[vertex];
1010 for (
unsigned int d=0; d<spacedim; ++d)
1011 shift_vector(d) = uniform_distribution(rng);
1013 shift_vector *= factor * minimal_length[vertex] /
1014 std::sqrt(shift_vector.
square());
1017 new_vertex_locations[vertex] = old_vertex_locations[vertex] + shift_vector;
1023 cell=triangulation.begin_active(); cell!=triangulation.end(); ++cell)
1024 for (
unsigned int vertex_no=0;
1025 vertex_no<GeometryInfo<dim>::vertices_per_cell; ++vertex_no)
1026 cell->vertex(vertex_no) = new_vertex_locations[cell->vertex_index(vertex_no)];
1040 cell = triangulation.begin_active(),
1041 endc = triangulation.end();
1042 for (; cell!=endc; ++cell)
1043 if (!cell->is_artificial())
1044 for (
unsigned int face=0;
1045 face<GeometryInfo<dim>::faces_per_cell; ++face)
1046 if (cell->face(face)->has_children() &&
1047 !cell->face(face)->at_boundary())
1051 cell->face(face)->child(0)->vertex(1)
1052 = (cell->face(face)->vertex(0) +
1053 cell->face(face)->vertex(1)) / 2;
1056 cell->face(face)->child(0)->vertex(1)
1057 = .5*(cell->face(face)->vertex(0)
1058 +cell->face(face)->vertex(1));
1059 cell->face(face)->child(0)->vertex(2)
1060 = .5*(cell->face(face)->vertex(0)
1061 +cell->face(face)->vertex(2));
1062 cell->face(face)->child(1)->vertex(3)
1063 = .5*(cell->face(face)->vertex(1)
1064 +cell->face(face)->vertex(3));
1065 cell->face(face)->child(2)->vertex(3)
1066 = .5*(cell->face(face)->vertex(2)
1067 +cell->face(face)->vertex(3));
1070 cell->face(face)->child(0)->vertex(3)
1071 = .25*(cell->face(face)->vertex(0)
1072 +cell->face(face)->vertex(1)
1073 +cell->face(face)->vertex(2)
1074 +cell->face(face)->vertex(3));
1082 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1086 const std::vector<bool> &marked_vertices)
1094 const std::vector< Point<spacedim> > &vertices = tria.
get_vertices();
1096 Assert ( tria.
get_vertices().size() == marked_vertices.size() || marked_vertices.size() ==0,
1106 Assert( marked_vertices.size()==0
1108 std::equal( marked_vertices.begin(),
1109 marked_vertices.end(),
1115 ExcMessage(
"marked_vertices should be a subset of used vertices in the triangulation " 1116 "but marked_vertices contains one or more vertices that are not used vertices!") );
1124 const std::vector<bool> &used =
1129 std::vector<bool>::const_iterator first =
1130 std::find(used.begin(), used.end(),
true);
1136 unsigned int best_vertex = std::distance(used.begin(), first);
1137 double best_dist = (p - vertices[best_vertex]).norm_square();
1141 for (
unsigned int j = best_vertex+1; j < vertices.size(); j++)
1144 double dist = (p - vertices[j]).norm_square();
1145 if (dist < best_dist)
1157 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1160 const MeshType<dim,spacedim> &mesh,
1162 const std::vector<bool> &marked_vertices)
1176 Assert ( tria.
get_vertices().size() == marked_vertices.size() || marked_vertices.size() ==0,
1186 Assert( marked_vertices.size()==0
1188 std::equal( marked_vertices.begin(),
1189 marked_vertices.end(),
1195 ExcMessage(
"marked_vertices should be a subset of used vertices in the triangulation " 1196 "but marked_vertices contains one or more vertices that are not used vertices!") );
1199 if (marked_vertices.size())
1200 for (
auto it = vertices.begin(); it != vertices.end(); )
1202 if (marked_vertices[it->first] ==
false)
1204 vertices.erase(it++);
1217 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1219 std::vector<typename MeshType<dim, spacedim>::active_cell_iterator>
1221 std::vector<typename ::internal::ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim> >::type>
1224 const unsigned int vertex)
1229 Assert(vertex < mesh.get_triangulation().n_vertices(),
1230 ExcIndexRange(0,mesh.get_triangulation().n_vertices(),vertex));
1231 Assert(mesh.get_triangulation().get_used_vertices()[vertex],
1237 std::set<typename ::internal::ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim> >::type> adjacent_cells;
1239 typename ::internal::ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim> >::type
1276 for (; cell != endc; ++cell)
1278 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell; v++)
1279 if (cell->vertex_index(v) == vertex)
1284 adjacent_cells.insert(cell);
1291 for (
unsigned int vface = 0; vface < dim; vface++)
1293 const unsigned int face =
1296 if (!cell->at_boundary(face)
1298 cell->neighbor(face)->active())
1310 adjacent_cells.insert (cell->neighbor(face));
1321 for (
unsigned int e=0; e<GeometryInfo<dim>::lines_per_cell; ++e)
1322 if (cell->line(e)->has_children())
1326 if (cell->line(e)->child(0)->vertex_index(1) == vertex)
1328 adjacent_cells.insert(cell);
1351 std::vector<typename ::internal::ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim> >::type>
1352 (adjacent_cells.begin(), adjacent_cells.end());
1359 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1360 void find_active_cell_around_point_internal
1361 (
const MeshType<dim,spacedim> &mesh,
1363 std::set<
typename MeshType<dim, spacedim>::active_cell_iterator> &searched_cells,
1364 std::set<
typename MeshType<dim, spacedim>::active_cell_iterator> &adjacent_cells)
1366 std::set<typename ::internal::ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim> >::type> &searched_cells,
1367 std::set<typename ::internal::ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim> >::type> &adjacent_cells)
1371 typedef typename MeshType<dim, spacedim>::active_cell_iterator cell_iterator;
1373 typedef typename ::internal::ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim> >::type cell_iterator;
1377 searched_cells.insert(adjacent_cells.begin(), adjacent_cells.end());
1381 std::set<cell_iterator> adjacent_cells_new;
1383 typename std::set<cell_iterator>::const_iterator
1384 cell = adjacent_cells.
begin(),
1385 endc = adjacent_cells.end();
1386 for (; cell != endc; ++cell)
1388 std::vector<cell_iterator> active_neighbors;
1389 get_active_neighbors<MeshType<dim, spacedim> >(*cell, active_neighbors);
1390 for (
unsigned int i=0; i<active_neighbors.size(); ++i)
1391 if (searched_cells.find(active_neighbors[i]) == searched_cells.end())
1392 adjacent_cells_new.insert(active_neighbors[i]);
1394 adjacent_cells.clear();
1395 adjacent_cells.insert(adjacent_cells_new.begin(), adjacent_cells_new.end());
1396 if (adjacent_cells.size() == 0)
1404 cell_iterator it = mesh.begin_active();
1405 for ( ; it!=mesh.end(); ++it)
1406 if (searched_cells.find(it) == searched_cells.end())
1408 adjacent_cells.insert(it);
1415 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1417 typename MeshType<dim, spacedim>::active_cell_iterator
1419 typename ::internal::ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim> >::type
1423 const std::vector<bool> &marked_vertices)
1426 find_active_cell_around_point<dim,MeshType,spacedim>
1428 mesh, p, marked_vertices).first;
1432 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1434 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim> >
1436 std::pair<typename ::internal::ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim> >::type,
Point<dim> >
1439 const MeshType<dim,spacedim> &mesh,
1441 const std::vector<bool> &marked_vertices)
1443 typedef typename ::internal::ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim> >::type active_cell_iterator;
1449 double best_distance = 1e-10;
1450 int best_level = -1;
1451 std::pair<active_cell_iterator, Point<dim> > best_cell;
1455 std::vector<active_cell_iterator> adjacent_cells_tmp
1464 std::set<active_cell_iterator> adjacent_cells (adjacent_cells_tmp.begin(),
1465 adjacent_cells_tmp.end());
1466 std::set<active_cell_iterator> searched_cells;
1474 const unsigned int n_active_cells = mesh.get_triangulation().n_active_cells();
1476 unsigned int cells_searched = 0;
1477 while (!found && cells_searched < n_active_cells)
1479 typename std::set<active_cell_iterator>::const_iterator
1480 cell = adjacent_cells.begin(),
1481 endc = adjacent_cells.end();
1482 for (; cell != endc; ++cell)
1496 if ((dist < best_distance)
1498 ((dist == best_distance)
1500 ((*cell)->level() > best_level)))
1503 best_distance = dist;
1504 best_level = (*cell)->level();
1505 best_cell = std::make_pair(*cell, p_cell);
1525 cells_searched += adjacent_cells.size();
1530 if (marked_vertices.size() > 0)
1531 cells_searched = n_active_cells;
1540 if (!found && cells_searched < n_active_cells)
1542 find_active_cell_around_point_internal<dim,MeshType,spacedim>
1543 (mesh, searched_cells, adjacent_cells);
1548 ExcPointNotFound<spacedim>(p));
1555 template <
int dim,
int spacedim>
1556 std::vector<std::vector<Tensor<1,spacedim> > >
1560 const std::vector<Point<spacedim> > &vertices = mesh.
get_vertices();
1561 const unsigned int n_vertices = vertex_to_cells.size();
1566 std::vector<std::vector<Tensor<1,spacedim> > > vertex_to_cell_centers(n_vertices);
1567 for (
unsigned int vertex=0; vertex<n_vertices; ++vertex)
1570 const unsigned int n_neighbor_cells = vertex_to_cells[vertex].size();
1571 vertex_to_cell_centers[vertex].resize(n_neighbor_cells);
1573 typename std::set<typename Triangulation<dim,spacedim>::active_cell_iterator>::iterator it = vertex_to_cells[vertex].begin();
1574 for (
unsigned int cell=0; cell<n_neighbor_cells; ++cell,++it)
1576 vertex_to_cell_centers[vertex][cell] = (*it)->center() - vertices[vertex];
1577 vertex_to_cell_centers[vertex][cell] /= vertex_to_cell_centers[vertex][cell].
norm();
1580 return vertex_to_cell_centers;
1586 template <
int spacedim>
1588 compare_point_association(
const unsigned int a,
1589 const unsigned int b,
1593 const double scalar_product_a = center_directions[a] * point_direction;
1594 const double scalar_product_b = center_directions[b] * point_direction;
1599 return (scalar_product_a > scalar_product_b);
1603 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
1605 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim> >
1607 std::pair<typename ::internal::ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim> >::type,
Point<dim> >
1610 const MeshType<dim,spacedim> &mesh,
1612 const std::vector<std::set<
typename MeshType<dim,spacedim>::active_cell_iterator > > &vertex_to_cells,
1614 const typename MeshType<dim, spacedim>::active_cell_iterator &cell_hint ,
1615 const std::vector<bool> &marked_vertices)
1617 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim> > cell_and_position;
1619 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim> > cell_and_position_approx;
1621 bool found_cell =
false;
1622 bool approx_cell =
false;
1624 unsigned int closest_vertex_index = 0;
1626 auto current_cell = cell_hint;
1628 while (found_cell ==
false)
1634 const unsigned int closest_vertex = find_closest_vertex_of_cell<dim,spacedim>(current_cell , p);
1635 vertex_to_point = p - current_cell ->vertex(closest_vertex);
1636 closest_vertex_index = current_cell ->vertex_index(closest_vertex);
1641 vertex_to_point = p - mesh.get_vertices()[closest_vertex_index];
1644 const double vertex_point_norm = vertex_to_point.
norm();
1645 if (vertex_point_norm > 0)
1646 vertex_to_point /= vertex_point_norm;
1648 const unsigned int n_neighbor_cells = vertex_to_cells[closest_vertex_index].size();
1651 std::vector<unsigned int> neighbor_permutation(n_neighbor_cells);
1653 for (
unsigned int i=0; i<n_neighbor_cells; ++i)
1654 neighbor_permutation[i] = i;
1656 auto comp = [&](
const unsigned int a,
const unsigned int b) ->
bool 1658 return compare_point_association<spacedim>(a,b,vertex_to_point,vertex_to_cell_centers[closest_vertex_index]);
1661 std::sort(neighbor_permutation.begin(),
1662 neighbor_permutation.end(),
1668 double best_distance = 1e-10;
1672 for (
unsigned int i=0; i<n_neighbor_cells; ++i)
1676 auto cell = vertex_to_cells[closest_vertex_index].begin();
1677 std::advance(cell,neighbor_permutation[i]);
1681 cell_and_position.first = *cell;
1682 cell_and_position.second = p_unit;
1684 approx_cell =
false;
1691 if (dist < best_distance)
1693 best_distance = dist;
1694 cell_and_position_approx.first = *cell;
1695 cell_and_position_approx.second = p_unit;
1703 if (found_cell ==
true)
1704 return cell_and_position;
1705 else if (approx_cell ==
true)
1706 return cell_and_position_approx;
1720 ExcPointNotFound<spacedim>(p));
1722 current_cell =
typename MeshType<dim,spacedim>::active_cell_iterator();
1724 return cell_and_position;
1729 template <
int dim,
int spacedim>
1735 unsigned int closest_vertex = 0;
1737 for (
unsigned int v=1; v<GeometryInfo<dim>::vertices_per_cell; ++v)
1739 const double vertex_distance = position.
distance_square(cell->vertex(v));
1740 if (vertex_distance < minimum_distance)
1743 minimum_distance = vertex_distance;
1746 return closest_vertex;
1753 namespace BoundingBoxPredicate
1755 template <
class MeshType >
1756 std::tuple< BoundingBox < MeshType::space_dimension >,
bool >
1757 compute_cell_predicate_bounding_box
1758 (
const typename MeshType::cell_iterator &parent_cell,
1759 const std::function<
bool (
const typename MeshType::active_cell_iterator &)> &predicate)
1761 bool has_predicate =
false;
1762 std::vector< typename MeshType::active_cell_iterator > active_cells;
1763 if (parent_cell->active())
1764 active_cells = {parent_cell};
1767 active_cells = get_active_child_cells < MeshType > (parent_cell);
1769 const unsigned int spacedim = MeshType::space_dimension;
1773 while ( i < active_cells.size() && !predicate(active_cells[i]) )
1777 if ( active_cells.size() == 0 || i == active_cells.size() )
1780 return std::make_tuple(bbox, has_predicate);
1787 for (; i < active_cells.size() ; ++i)
1788 if ( predicate(active_cells[i]) )
1789 for (
unsigned int v=0; v<GeometryInfo<spacedim>::vertices_per_cell; ++v)
1790 for (
unsigned int d=0;
d<spacedim; ++
d)
1792 minp[
d] = std::min( minp[d], active_cells[i]->vertex(v)[d]);
1793 maxp[
d] = std::max( maxp[d], active_cells[i]->vertex(v)[d]);
1796 has_predicate =
true;
1798 return std::make_tuple(bbox, has_predicate);
1805 template <
class MeshType >
1806 std::vector< BoundingBox<MeshType::space_dimension> >
1808 (
const MeshType &mesh,
1809 const std::function<
bool (
const typename MeshType::active_cell_iterator &)> &predicate,
1810 const unsigned int refinement_level,
1811 const bool allow_merge,
1812 const unsigned int max_boxes)
1818 Assert( refinement_level <= mesh.n_levels(),
1819 ExcMessage (
"Error: refinement level is higher then total levels in the triangulation!") );
1821 const unsigned int spacedim = MeshType::space_dimension;
1822 std::vector< BoundingBox < spacedim > > bounding_boxes;
1826 for (
unsigned int i=0; i < refinement_level; ++i)
1827 for (
typename MeshType::cell_iterator cell: mesh.active_cell_iterators_on_level(i))
1829 bool has_predicate =
false;
1831 std::tie(bbox, has_predicate) =
1832 internal::BoundingBoxPredicate::compute_cell_predicate_bounding_box <MeshType> (cell, predicate);
1834 bounding_boxes.push_back(bbox);
1838 for (
const typename MeshType::cell_iterator &cell: mesh.cell_iterators_on_level(refinement_level))
1840 bool has_predicate =
false;
1842 std::tie(bbox, has_predicate) =
1843 internal::BoundingBoxPredicate::compute_cell_predicate_bounding_box <MeshType> (cell, predicate);
1845 bounding_boxes.push_back(bbox);
1850 return bounding_boxes;
1856 std::vector<unsigned int> merged_boxes_idx;
1857 bool found_neighbors =
true;
1861 while (found_neighbors)
1863 found_neighbors =
false;
1864 for (
unsigned int i=0; i<bounding_boxes.size()-1; ++i)
1866 if ( std::find(merged_boxes_idx.begin(),merged_boxes_idx.end(),i) == merged_boxes_idx.end())
1867 for (
unsigned int j=i+1; j<bounding_boxes.size(); ++j)
1868 if ( std::find(merged_boxes_idx.begin(),merged_boxes_idx.end(),j) == merged_boxes_idx.end()
1869 && bounding_boxes[i].get_neighbor_type (bounding_boxes[j]) == NeighborType::mergeable_neighbors )
1871 bounding_boxes[i].merge_with(bounding_boxes[j]);
1872 merged_boxes_idx.push_back(j);
1873 found_neighbors =
true;
1879 std::vector< BoundingBox < spacedim > > merged_b_boxes;
1880 for (
unsigned int i=0; i<bounding_boxes.size(); ++i)
1881 if (std::find(merged_boxes_idx.begin(),merged_boxes_idx.end(),i) == merged_boxes_idx.end())
1882 merged_b_boxes.push_back(bounding_boxes[i]);
1887 if ( (merged_b_boxes.size() > max_boxes) && (spacedim > 1) )
1889 std::vector<double> volumes;
1890 for (
unsigned int i=0; i< merged_b_boxes.size(); ++i)
1891 volumes.push_back(merged_b_boxes[i].volume());
1893 while ( merged_b_boxes.size() > max_boxes)
1895 unsigned int min_idx = std::min_element(volumes.begin(),volumes.end()) -
1897 volumes.erase(volumes.begin() + min_idx);
1899 bool not_removed =
true;
1900 for (
unsigned int i=0; i<merged_b_boxes.size() && not_removed; ++i)
1903 if ( i != min_idx &&
1904 (merged_b_boxes[i].get_neighbor_type (merged_b_boxes[min_idx]) == NeighborType::attached_neighbors ||
1905 merged_b_boxes[i].get_neighbor_type (merged_b_boxes[min_idx]) == NeighborType::mergeable_neighbors) )
1907 merged_b_boxes[i].merge_with(merged_b_boxes[min_idx]);
1908 merged_b_boxes.erase(merged_b_boxes.begin() + min_idx);
1909 not_removed =
false;
1912 ExcMessage (
"Error: couldn't merge bounding boxes!") );
1915 Assert( merged_b_boxes.size() <= max_boxes,
1916 ExcMessage (
"Error: couldn't reach target number of bounding boxes!") );
1917 return merged_b_boxes;
1924 template <
int spacedim>
1925 std::tuple< std::vector< std::vector< unsigned int > >,
1926 std::map< unsigned int, unsigned int>,
1927 std::map< unsigned int, std::vector< unsigned int > > >
1932 unsigned int n_procs = global_bboxes.size();
1933 std::vector< std::vector< unsigned int > > point_owners(n_procs);
1934 std::map< unsigned int, unsigned int> map_owners_found;
1935 std::map< unsigned int, std::vector< unsigned int > > map_owners_guessed;
1937 unsigned int n_points = points.size();
1938 for (
unsigned int pt=0; pt<n_points; ++pt)
1941 std::vector< unsigned int > owners_found;
1943 for (
unsigned int rk=0; rk<n_procs; ++rk)
1946 if (bbox.point_inside(points[pt]))
1948 point_owners[rk].emplace_back(pt);
1949 owners_found.emplace_back(rk);
1953 Assert(owners_found.size() > 0,
1954 ExcMessage(
"No owners found for the point " + std::to_string(pt)));
1955 if (owners_found.size()==1)
1956 map_owners_found[pt] = owners_found[0];
1959 map_owners_guessed[pt] = owners_found;
1962 std::tuple< std::vector< std::vector< unsigned int > >,
1963 std::map< unsigned int, unsigned int>,
1964 std::map< unsigned int, std::vector< unsigned int > > >
1967 std::get<0>(output_tuple) = point_owners;
1968 std::get<1>(output_tuple) = map_owners_found;
1969 std::get<2>(output_tuple) = map_owners_guessed;
1971 return output_tuple;
1976 template <
int dim,
int spacedim>
1977 std::vector<std::set<typename Triangulation<dim,spacedim>::active_cell_iterator> >
1980 std::vector<std::set<typename Triangulation<dim,spacedim>::active_cell_iterator> >
1983 endc = triangulation.end();
1984 for (; cell!=endc; ++cell)
1985 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
1989 cell = triangulation.begin_active();
1990 for (; cell!=endc; ++cell)
1992 for (
unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
1994 if ((cell->at_boundary(i)==
false) && (cell->neighbor(i)->active()))
1998 for (
unsigned int j=0; j<GeometryInfo<dim>::vertices_per_face; ++j)
2006 for (
unsigned int i=0; i<GeometryInfo<dim>::lines_per_cell; ++i)
2007 if (cell->line(i)->has_children())
2020 template <
int dim,
int spacedim>
2021 std::map<unsigned int,types::global_vertex_index>
2025 std::map<unsigned int,types::global_vertex_index> local_to_global_vertex_index;
2027 #ifndef DEAL_II_WITH_MPI 2032 (void)triangulation;
2034 "for parallel::distributed::Triangulation " 2035 "objects if you do not have MPI enabled."));
2040 const std::vector<std::set<active_cell_iterator> > vertex_to_cell =
2045 unsigned int max_cellid_size = 0;
2046 std::set<std::pair<types::subdomain_id,types::global_vertex_index> > vertices_added;
2047 std::map<types::subdomain_id,std::set<unsigned int> > vertices_to_recv;
2050 active_cell_iterator cell = triangulation.begin_active(),
2051 endc = triangulation.end();
2052 std::set<active_cell_iterator> missing_vert_cells;
2053 std::set<unsigned int> used_vertex_index;
2054 for (; cell!=endc; ++cell)
2056 if (cell->is_locally_owned())
2058 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
2061 typename std::set<active_cell_iterator>::iterator
2062 adjacent_cell = vertex_to_cell[cell->vertex_index(i)].begin(),
2063 end_adj_cell = vertex_to_cell[cell->vertex_index(i)].end();
2064 for (; adjacent_cell!=end_adj_cell; ++adjacent_cell)
2065 lowest_subdomain_id = std::min(lowest_subdomain_id,
2066 (*adjacent_cell)->subdomain_id());
2069 if (lowest_subdomain_id==cell->subdomain_id())
2073 if (used_vertex_index.find(cell->vertex_index(i))==used_vertex_index.end())
2076 local_to_global_vertex_index[cell->vertex_index(i)] = next_index++;
2080 adjacent_cell = vertex_to_cell[cell->vertex_index(i)].begin();
2081 for (; adjacent_cell!=end_adj_cell; ++adjacent_cell)
2082 if ((*adjacent_cell)->subdomain_id()!=cell->subdomain_id())
2084 std::pair<types::subdomain_id,types::global_vertex_index>
2085 tmp((*adjacent_cell)->subdomain_id(), cell->vertex_index(i));
2086 if (vertices_added.find(tmp)==vertices_added.end())
2088 vertices_to_send[(*adjacent_cell)->subdomain_id()].emplace_back
2089 (i, cell->vertex_index(i), cell->id().to_string());
2090 if (cell->id().to_string().size() > max_cellid_size)
2091 max_cellid_size = cell->id().to_string().size();
2092 vertices_added.insert(tmp);
2095 used_vertex_index.insert(cell->vertex_index(i));
2101 vertices_to_recv[lowest_subdomain_id].insert(cell->vertex_index(i));
2102 missing_vert_cells.insert(cell);
2109 if (cell->is_ghost())
2111 for (
unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
2113 if (cell->at_boundary(i)==
false)
2115 if (cell->neighbor(i)->active())
2119 if ((adjacent_cell->is_locally_owned()))
2122 if (cell->subdomain_id()<adj_subdomain_id)
2123 for (
unsigned int j=0; j<GeometryInfo<dim>::vertices_per_face; ++j)
2125 vertices_to_recv[cell->subdomain_id()].insert(cell->face(i)->vertex_index(j));
2126 missing_vert_cells.insert(cell);
2141 std::vector<types::global_vertex_index> indices(n_cpu);
2142 int ierr = MPI_Allgather(&next_index, 1, DEAL_II_VERTEX_INDEX_MPI_TYPE, indices.data(),
2143 1, DEAL_II_VERTEX_INDEX_MPI_TYPE, triangulation.get_communicator());
2145 Assert(indices.begin() + triangulation.locally_owned_subdomain() < indices.end(),
2148 indices.begin()+triangulation.locally_owned_subdomain(),
2151 std::map<unsigned int,types::global_vertex_index>::iterator
2152 global_index_it = local_to_global_vertex_index.begin(),
2153 global_index_end = local_to_global_vertex_index.end();
2154 for (; global_index_it!=global_index_end; ++global_index_it)
2155 global_index_it->second +=
shift;
2162 std::vector<std::vector<types::global_vertex_index> > vertices_send_buffers(
2163 vertices_to_send.size());
2164 std::vector<MPI_Request> first_requests(vertices_to_send.size());
2168 vert_to_send_it = vertices_to_send.begin(),
2169 vert_to_send_end = vertices_to_send.end();
2170 for (
unsigned int i=0; vert_to_send_it!=vert_to_send_end;
2171 ++vert_to_send_it, ++i)
2173 int destination = vert_to_send_it->first;
2174 const unsigned int n_vertices = vert_to_send_it->second.size();
2175 const int buffer_size = 2*n_vertices;
2176 vertices_send_buffers[i].resize(buffer_size);
2179 for (
unsigned int j=0; j<n_vertices; ++j)
2181 vertices_send_buffers[i][2*j] = std::get<0>(vert_to_send_it->second[j]);
2182 vertices_send_buffers[i][2*j+1] =
2183 local_to_global_vertex_index[std::get<1>(vert_to_send_it->second[j])];
2187 ierr = MPI_Isend(&vertices_send_buffers[i][0],buffer_size,DEAL_II_VERTEX_INDEX_MPI_TYPE,
2188 destination, 0, triangulation.get_communicator(), &first_requests[i]);
2193 std::vector<std::vector<types::global_vertex_index> > vertices_recv_buffers(
2194 vertices_to_recv.size());
2195 typename std::map<types::subdomain_id,std::set<unsigned int> >::iterator
2196 vert_to_recv_it = vertices_to_recv.begin(),
2197 vert_to_recv_end = vertices_to_recv.end();
2198 for (
unsigned int i=0; vert_to_recv_it!=vert_to_recv_end; ++vert_to_recv_it, ++i)
2200 int source = vert_to_recv_it->first;
2201 const unsigned int n_vertices = vert_to_recv_it->second.size();
2202 const int buffer_size = 2*n_vertices;
2203 vertices_recv_buffers[i].resize(buffer_size);
2206 ierr = MPI_Recv(&vertices_recv_buffers[i][0],buffer_size,DEAL_II_VERTEX_INDEX_MPI_TYPE,
2207 source, 0, triangulation.get_communicator(), MPI_STATUS_IGNORE);
2213 std::vector<std::vector<char> > cellids_send_buffers(vertices_to_send.size());
2214 std::vector<MPI_Request> second_requests(vertices_to_send.size());
2215 vert_to_send_it = vertices_to_send.begin();
2216 for (
unsigned int i=0; vert_to_send_it!=vert_to_send_end;
2217 ++vert_to_send_it, ++i)
2219 int destination = vert_to_send_it->first;
2220 const unsigned int n_vertices = vert_to_send_it->second.size();
2221 const int buffer_size = max_cellid_size*n_vertices;
2222 cellids_send_buffers[i].resize(buffer_size);
2225 unsigned int pos = 0;
2226 for (
unsigned int j=0; j<n_vertices; ++j)
2228 std::string cell_id = std::get<2>(vert_to_send_it->second[j]);
2229 for (
unsigned int k=0; k<max_cellid_size; ++k, ++pos)
2231 if (k<cell_id.size())
2232 cellids_send_buffers[i][pos] = cell_id[k];
2236 cellids_send_buffers[i][pos] =
'-';
2241 ierr = MPI_Isend(&cellids_send_buffers[i][0], buffer_size, MPI_CHAR,
2242 destination, 0, triangulation.get_communicator(), &second_requests[i]);
2247 std::vector<std::vector<char> > cellids_recv_buffers(vertices_to_recv.size());
2248 vert_to_recv_it = vertices_to_recv.begin();
2249 for (
unsigned int i=0; vert_to_recv_it!=vert_to_recv_end; ++vert_to_recv_it, ++i)
2251 int source = vert_to_recv_it->first;
2252 const unsigned int n_vertices = vert_to_recv_it->second.size();
2253 const int buffer_size = max_cellid_size*n_vertices;
2254 cellids_recv_buffers[i].resize(buffer_size);
2257 ierr = MPI_Recv(&cellids_recv_buffers[i][0],buffer_size, MPI_CHAR,
2258 source, 0, triangulation.get_communicator(), MPI_STATUS_IGNORE);
2264 vert_to_recv_it = vertices_to_recv.begin();
2265 for (
unsigned int i=0; vert_to_recv_it!=vert_to_recv_end; ++i, ++vert_to_recv_it)
2267 for (
unsigned int j=0; j<vert_to_recv_it->second.size(); ++j)
2269 const unsigned int local_pos_recv = vertices_recv_buffers[i][2*j];
2271 const std::string cellid_recv(&cellids_recv_buffers[i][max_cellid_size*j],
2272 &cellids_recv_buffers[i][max_cellid_size*(j+1)]);
2274 typename std::set<active_cell_iterator>::iterator
2275 cell_set_it = missing_vert_cells.begin(),
2276 end_cell_set = missing_vert_cells.end();
2277 for (; (found==
false) && (cell_set_it!=end_cell_set); ++cell_set_it)
2279 typename std::set<active_cell_iterator>::iterator
2280 candidate_cell = vertex_to_cell[(*cell_set_it)->vertex_index(i)].begin(),
2281 end_cell = vertex_to_cell[(*cell_set_it)->vertex_index(i)].end();
2282 for (; candidate_cell!=end_cell; ++candidate_cell)
2284 std::string current_cellid = (*candidate_cell)->id().to_string();
2285 current_cellid.resize(max_cellid_size,
'-');
2286 if (current_cellid.compare(cellid_recv)==0)
2288 local_to_global_vertex_index[(*candidate_cell)->vertex_index(local_pos_recv)] =
2300 return local_to_global_vertex_index;
2305 template <
int dim,
int spacedim>
2310 cell_connectivity.
reinit (triangulation.n_active_cells(),
2311 triangulation.n_active_cells());
2317 std::map< std::pair<unsigned int,unsigned int>,
unsigned int >
2320 cell = triangulation.begin_active();
2321 cell != triangulation.end(); ++cell)
2322 indexmap[std::pair<unsigned int,unsigned int>(cell->level(),cell->index())] = cell->active_cell_index();
2332 cell = triangulation.begin_active();
2333 cell != triangulation.end(); ++cell)
2335 const unsigned int index = cell->active_cell_index();
2336 cell_connectivity.
add (index, index);
2337 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
2338 if ((cell->at_boundary(f) ==
false)
2340 (cell->neighbor(f)->has_children() ==
false))
2342 unsigned int other_index = indexmap.find(
2343 std::pair<unsigned int,unsigned int>(cell->neighbor(f)->level(),cell->neighbor(f)->index()))->second;
2344 cell_connectivity.
add (index, other_index);
2345 cell_connectivity.
add (other_index, index);
2352 template <
int dim,
int spacedim>
2357 std::vector<std::vector<unsigned int> > vertex_to_cell(triangulation.n_vertices());
2359 triangulation.begin_active(); cell != triangulation.end(); ++cell)
2361 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
2362 vertex_to_cell[cell->vertex_index(v)].push_back(cell->active_cell_index());
2365 cell_connectivity.
reinit (triangulation.n_active_cells(),
2366 triangulation.n_active_cells());
2368 triangulation.begin_active(); cell != triangulation.end(); ++cell)
2370 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
2371 for (
unsigned int n=0; n<vertex_to_cell[cell->vertex_index(v)].size(); ++n)
2372 cell_connectivity.
add(cell->active_cell_index(), vertex_to_cell[cell->vertex_index(v)][n]);
2377 template <
int dim,
int spacedim>
2380 const unsigned int level,
2383 std::vector<std::vector<unsigned int> > vertex_to_cell(triangulation.n_vertices());
2385 triangulation.begin(level); cell != triangulation.end(level); ++cell)
2387 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
2388 vertex_to_cell[cell->vertex_index(v)].push_back(cell->index());
2391 cell_connectivity.
reinit (triangulation.n_cells(level),
2392 triangulation.n_cells(level));
2394 triangulation.begin(level); cell != triangulation.end(level); ++cell)
2396 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
2397 for (
unsigned int n=0; n<vertex_to_cell[cell->vertex_index(v)].size(); ++n)
2398 cell_connectivity.
add(cell->index(), vertex_to_cell[cell->vertex_index(v)][n]);
2404 template <
int dim,
int spacedim>
2411 std::vector<unsigned int> cell_weights;
2414 if (!triangulation.signals.cell_weight.empty())
2416 cell_weights.resize(triangulation.n_active_cells(), std::numeric_limits<unsigned int>::max());
2420 cell = triangulation.begin_active(),
2421 endc = triangulation.end();
2422 for (; cell!=endc; ++cell, ++c)
2428 triangulation, partitioner);
2433 template <
int dim,
int spacedim>
2436 const std::vector<unsigned int> &cell_weights,
2444 ExcMessage (
"Objects of type parallel::distributed::Triangulation " 2445 "are already partitioned implicitly and can not be " 2446 "partitioned again explicitly."));
2450 if (n_partitions == 1)
2453 cell = triangulation.begin_active();
2454 cell != triangulation.end(); ++cell)
2455 cell->set_subdomain_id (0);
2469 sp_cell_connectivity.
copy_from(cell_connectivity);
2472 sp_cell_connectivity,
2480 template <
int dim,
int spacedim>
2488 std::vector<unsigned int> cell_weights;
2491 if (!triangulation.signals.cell_weight.empty())
2493 cell_weights.resize(triangulation.n_active_cells(), std::numeric_limits<unsigned int>::max());
2497 cell = triangulation.begin_active(),
2498 endc = triangulation.end();
2499 for (; cell!=endc; ++cell, ++c)
2505 cell_connection_graph,
2506 triangulation, partitioner);
2511 template <
int dim,
int spacedim>
2514 const std::vector<unsigned int> &cell_weights,
2523 ExcMessage (
"Objects of type parallel::distributed::Triangulation " 2524 "are already partitioned implicitly and can not be " 2525 "partitioned again explicitly."));
2527 Assert (cell_connection_graph.
n_rows() == triangulation.n_active_cells(),
2528 ExcMessage (
"Connectivity graph has wrong size"));
2529 Assert (cell_connection_graph.
n_cols() == triangulation.n_active_cells(),
2530 ExcMessage (
"Connectivity graph has wrong size"));
2533 if (n_partitions == 1)
2536 cell = triangulation.begin_active();
2537 cell != triangulation.end(); ++cell)
2538 cell->set_subdomain_id (0);
2546 std::vector<unsigned int> partition_indices (triangulation.n_active_cells());
2548 partition_indices, partitioner);
2553 cell = triangulation.begin_active();
2554 cell != triangulation.end(); ++cell)
2555 cell->set_subdomain_id (partition_indices[cell->active_cell_index()]);
2565 void set_subdomain_id_in_zorder_recursively(IT cell,
2566 unsigned int ¤t_proc_idx,
2567 unsigned int ¤t_cell_idx,
2568 const unsigned int n_active_cells,
2569 const unsigned int n_partitions)
2573 while (current_cell_idx >= floor((
long)n_active_cells*(current_proc_idx+1)/n_partitions))
2575 cell->set_subdomain_id(current_proc_idx);
2580 for (
unsigned int n=0; n<cell->n_children(); ++n)
2581 set_subdomain_id_in_zorder_recursively(cell->child(n),
2590 template <
int dim,
int spacedim>
2598 ExcMessage (
"Objects of type parallel::distributed::Triangulation " 2599 "are already partitioned implicitly and can not be " 2600 "partitioned again explicitly."));
2604 if (n_partitions == 1)
2607 cell = triangulation.begin_active();
2608 cell != triangulation.end(); ++cell)
2609 cell->set_subdomain_id (0);
2615 std::vector<types::global_dof_index> coarse_cell_to_p4est_tree_permutation;
2616 std::vector<types::global_dof_index> p4est_tree_to_coarse_cell_permutation;
2620 coarse_cell_to_p4est_tree_permutation.resize (triangulation.n_cells(0));
2622 coarse_cell_to_p4est_tree_permutation);
2624 p4est_tree_to_coarse_cell_permutation
2627 unsigned int current_proc_idx=0;
2628 unsigned int current_cell_idx=0;
2629 const unsigned int n_active_cells = triangulation.n_active_cells();
2633 for (
unsigned int idx=0; idx<triangulation.n_cells(0); ++idx)
2635 const unsigned int coarse_cell_idx = p4est_tree_to_coarse_cell_permutation[idx];
2637 coarse_cell (&triangulation, 0, coarse_cell_idx);
2639 set_subdomain_id_in_zorder_recursively(coarse_cell,
2655 cell = triangulation.begin(),
2656 endc = triangulation.end();
2657 for (; cell!=endc; ++cell)
2661 bool all_children_active =
true;
2662 std::map<unsigned int, unsigned int> map_cpu_n_cells;
2663 for (
unsigned int n=0; n<cell->n_children(); ++n)
2664 if (!cell->child(n)->active())
2666 all_children_active =
false;
2670 ++map_cpu_n_cells[cell->child(n)->subdomain_id()];
2672 if (!all_children_active)
2675 unsigned int new_owner = cell->child(0)->subdomain_id();
2676 for (std::map<unsigned int, unsigned int>::iterator it = map_cpu_n_cells.begin();
2677 it != map_cpu_n_cells.end();
2679 if (it->second > map_cpu_n_cells[new_owner])
2680 new_owner = it->first;
2682 for (
unsigned int n=0; n<cell->n_children(); ++n)
2683 cell->child(n)->set_subdomain_id(new_owner);
2689 template <
int dim,
int spacedim>
2693 unsigned int n_levels = triangulation.n_levels();
2694 for (
int lvl = n_levels-1; lvl>=0; --lvl)
2697 cell = triangulation.begin(lvl),
2698 endc = triangulation.end(lvl);
2699 for (; cell!=endc; ++cell)
2701 if (!cell->has_children())
2702 cell->set_level_subdomain_id(cell->subdomain_id());
2705 Assert(cell->child(0)->level_subdomain_id()
2707 cell->set_level_subdomain_id(cell->child(0)->level_subdomain_id());
2714 template <
int dim,
int spacedim>
2717 std::vector<types::subdomain_id> &subdomain)
2719 Assert (subdomain.size() == triangulation.n_active_cells(),
2721 triangulation.n_active_cells()));
2723 cell = triangulation.begin_active(); cell!=triangulation.end(); ++cell)
2724 subdomain[cell->active_cell_index()] = cell->subdomain_id();
2729 template <
int dim,
int spacedim>
2734 unsigned int count = 0;
2736 cell = triangulation.begin_active();
2737 cell!=triangulation.end(); ++cell)
2738 if (cell->subdomain_id() == subdomain)
2746 template <
int dim,
int spacedim>
2751 std::vector<bool> locally_owned_vertices = triangulation.get_used_vertices();
2761 cell = triangulation.begin_active();
2762 cell != triangulation.end(); ++cell)
2763 if (cell->is_artificial()
2765 (cell->is_ghost() &&
2766 (cell->subdomain_id() < tr->locally_owned_subdomain())))
2767 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
2768 locally_owned_vertices[cell->vertex_index(v)] =
false;
2770 return locally_owned_vertices;
2777 template <
int dim,
int spacedim>
2785 return (vertices[1]-vertices[0]).norm();
2787 return std::max((vertices[3]-vertices[0]).norm(),
2788 (vertices[2]-vertices[1]).norm());
2790 return std::max( std::max((vertices[7]-vertices[0]).norm(),
2791 (vertices[6]-vertices[1]).norm()),
2792 std::max((vertices[2]-vertices[5]).norm(),
2793 (vertices[3]-vertices[4]).norm()) );
2802 template <
int dim,
int spacedim>
2807 double min_diameter = std::numeric_limits<double>::max();
2808 for (
const auto &cell: triangulation.active_cell_iterators())
2809 if (!cell->is_artificial())
2810 min_diameter = std::min (min_diameter,
2811 diameter<dim,spacedim>(cell, mapping));
2813 double global_min_diameter = 0;
2815 #ifdef DEAL_II_WITH_MPI 2821 global_min_diameter = min_diameter;
2823 return global_min_diameter;
2828 template <
int dim,
int spacedim>
2833 double max_diameter = 0.;
2834 for (
const auto &cell: triangulation.active_cell_iterators())
2835 if (!cell->is_artificial())
2836 max_diameter = std::max (max_diameter,
2839 double global_max_diameter = 0;
2841 #ifdef DEAL_II_WITH_MPI 2847 global_max_diameter = max_diameter;
2849 return global_max_diameter;
2856 namespace FixUpDistortedChildCells
2879 template <
typename Iterator,
int spacedim>
2881 objective_function (
const Iterator &
object,
2884 const unsigned int structdim = Iterator::AccessorType::structure_dimension;
2885 Assert (spacedim == Iterator::AccessorType::dimension,
2898 Tensor<spacedim-structdim,spacedim> parent_alternating_forms
2901 for (
unsigned int i=0; i<GeometryInfo<structdim>::vertices_per_cell; ++i)
2902 parent_vertices[i] = object->vertex(i);
2905 parent_alternating_forms);
2907 const Tensor<spacedim-structdim,spacedim>
2908 average_parent_alternating_form
2909 = std::accumulate (parent_alternating_forms,
2923 Tensor<spacedim-structdim,spacedim> child_alternating_forms
2927 for (
unsigned int c=0; c<
object->n_children(); ++c)
2928 for (
unsigned int i=0; i<GeometryInfo<structdim>::vertices_per_cell; ++i)
2929 child_vertices[c][i] = object->child(c)->vertex(i);
2937 for (
unsigned int c=0; c<
object->n_children(); ++c)
2941 for (
unsigned int c=0; c<
object->n_children(); ++c)
2943 child_alternating_forms[c]);
2955 double objective = 0;
2956 for (
unsigned int c=0; c<
object->n_children(); ++c)
2957 for (
unsigned int i=0; i<GeometryInfo<structdim>::vertices_per_cell; ++i)
2958 objective += (child_alternating_forms[c][i] -
2959 average_parent_alternating_form/std::pow(2.,1.*structdim))
2971 template <
typename Iterator>
2973 get_face_midpoint (
const Iterator &
object,
2974 const unsigned int f,
2975 std::integral_constant<int, 1>)
2977 return object->vertex(f);
2987 template <
typename Iterator>
2989 get_face_midpoint (
const Iterator &
object,
2990 const unsigned int f,
2991 std::integral_constant<int, 2>)
2993 return object->line(f)->center();
3003 template <
typename Iterator>
3005 get_face_midpoint (
const Iterator &
object,
3006 const unsigned int f,
3007 std::integral_constant<int, 3>)
3009 return object->face(f)->center();
3037 template <
typename Iterator>
3039 minimal_diameter (
const Iterator &
object)
3042 structdim = Iterator::AccessorType::structure_dimension;
3044 double diameter =
object->diameter();
3045 for (
unsigned int f=0;
3046 f<GeometryInfo<structdim>::faces_per_cell;
3048 for (
unsigned int e=f+1;
3049 e<GeometryInfo<structdim>::faces_per_cell;
3054 std::integral_constant<int, structdim>())
3055 .distance (get_face_midpoint
3058 std::integral_constant<int, structdim>())));
3068 template <
typename Iterator>
3070 fix_up_object (
const Iterator &
object)
3072 const unsigned int structdim = Iterator::AccessorType::structure_dimension;
3073 const unsigned int spacedim = Iterator::AccessorType::space_dimension;
3089 unsigned int iteration = 0;
3090 const double diameter = minimal_diameter (
object);
3093 double current_value = objective_function (
object, object_mid_point);
3094 double initial_delta = 0;
3101 const double step_length =
diameter / 4 / (iteration + 1);
3106 for (
unsigned int d=0;
d<spacedim; ++
d)
3108 const double eps = step_length/10;
3113 gradient[
d] = (objective_function (
object,
3116 objective_function (
object,
3122 if (gradient.
norm() == 0)
3130 object_mid_point -= std::min(2 * current_value / (gradient*gradient),
3131 step_length / gradient.norm()) * gradient;
3135 const double previous_value = current_value;
3136 current_value = objective_function (
object, object_mid_point);
3139 initial_delta = (previous_value - current_value);
3142 if ((iteration >= 1) &&
3143 ((previous_value - current_value < 0)
3145 (std::fabs (previous_value - current_value)
3147 0.001 * initial_delta)))
3152 while (iteration < 20);
3168 double old_min_product, new_min_product;
3172 for (
unsigned int i=0; i<GeometryInfo<structdim>::vertices_per_cell; ++i)
3173 parent_vertices[i] = object->vertex(i);
3175 Tensor<spacedim-structdim,spacedim> parent_alternating_forms
3178 parent_alternating_forms);
3184 for (
unsigned int c=0; c<
object->n_children(); ++c)
3185 for (
unsigned int i=0; i<GeometryInfo<structdim>::vertices_per_cell; ++i)
3186 child_vertices[c][i] = object->child(c)->vertex(i);
3188 Tensor<spacedim-structdim,spacedim> child_alternating_forms
3192 for (
unsigned int c=0; c<
object->n_children(); ++c)
3194 child_alternating_forms[c]);
3196 old_min_product = child_alternating_forms[0][0] * parent_alternating_forms[0];
3197 for (
unsigned int c=0; c<
object->n_children(); ++c)
3198 for (
unsigned int i=0; i<GeometryInfo<structdim>::vertices_per_cell; ++i)
3199 for (
unsigned int j=0; j<GeometryInfo<structdim>::vertices_per_cell; ++j)
3201 std::min<double> (old_min_product,
3202 child_alternating_forms[c][i] *
3203 parent_alternating_forms[j]);
3211 for (
unsigned int c=0; c<
object->n_children(); ++c)
3215 for (
unsigned int c=0; c<
object->n_children(); ++c)
3217 child_alternating_forms[c]);
3219 new_min_product = child_alternating_forms[0][0] * parent_alternating_forms[0];
3220 for (
unsigned int c=0; c<
object->n_children(); ++c)
3221 for (
unsigned int i=0; i<GeometryInfo<structdim>::vertices_per_cell; ++i)
3222 for (
unsigned int j=0; j<GeometryInfo<structdim>::vertices_per_cell; ++j)
3224 std::min<double> (new_min_product,
3225 child_alternating_forms[c][i] *
3226 parent_alternating_forms[j]);
3234 if (new_min_product >= old_min_product)
3241 return (std::max (new_min_product, old_min_product) > 0);
3246 void fix_up_faces (const ::Triangulation<1,1>::cell_iterator &,
3247 std::integral_constant<int, 1>,
3248 std::integral_constant<int, 1>)
3256 template <
int dim,
int spacedim>
3257 void fix_up_faces (
const typename ::Triangulation<dim,spacedim>::cell_iterator &cell,
3258 std::integral_constant<int, dim>,
3259 std::integral_constant<int, spacedim>)
3269 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
3272 Assert (cell->face(f)->refinement_case() ==
3276 bool subface_is_more_refined =
false;
3277 for (
unsigned int g=0; g<GeometryInfo<dim>::max_children_per_face; ++g)
3278 if (cell->face(f)->child(g)->has_children())
3280 subface_is_more_refined =
true;
3284 if (subface_is_more_refined ==
true)
3288 fix_up_object (cell->face(f));
3295 template <
int dim,
int spacedim>
3312 ExcMessage(
"This function is only valid for a list of cells that " 3313 "have children (i.e., no cell in the list may be active)."));
3315 internal::FixUpDistortedChildCells
3316 ::fix_up_faces (cell,
3317 std::integral_constant<int, dim>(),
3318 std::integral_constant<int, spacedim>());
3321 if (!internal::FixUpDistortedChildCells::fix_up_object (cell))
3325 return unfixable_subset;
3330 template <
int dim,
int spacedim>
3332 const bool reset_boundary_ids)
3335 std::vector<types::manifold_id> dst_manifold_ids(src_boundary_ids.size());
3336 auto m_it = dst_manifold_ids.begin();
3337 for (
auto b : src_boundary_ids)
3342 const std::vector<types::boundary_id> reset_boundary_id =
3343 reset_boundary_ids ?
3344 std::vector<types::boundary_id>(src_boundary_ids.size(), 0) : src_boundary_ids;
3350 template <
int dim,
int spacedim>
3352 const std::vector<types::manifold_id> &dst_manifold_ids,
3354 const std::vector<types::boundary_id> &reset_boundary_ids_)
3357 const auto reset_boundary_ids = reset_boundary_ids_.size() ?
3358 reset_boundary_ids_ : src_boundary_ids;
3368 cell != tria.
end(); ++cell)
3369 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
3370 if (cell->face(f)->at_boundary())
3371 for (
unsigned int e=0; e<GeometryInfo<dim>::lines_per_face; ++e)
3373 const auto bid = cell->face(f)->line(e)->boundary_id();
3374 const auto ind = std::find(src_boundary_ids.begin(), src_boundary_ids.end(), bid)-
3375 src_boundary_ids.begin();
3376 if ((
unsigned int)ind < src_boundary_ids.size())
3377 cell->face(f)->line(e)->set_manifold_id(dst_manifold_ids[ind]);
3383 cell != tria.
end(); ++cell)
3384 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
3385 if (cell->face(f)->at_boundary())
3387 const auto bid = cell->face(f)->boundary_id();
3388 const auto ind = std::find(src_boundary_ids.begin(), src_boundary_ids.end(), bid)-
3389 src_boundary_ids.begin();
3391 if ((
unsigned int)ind < src_boundary_ids.size())
3394 cell->face(f)->set_manifold_id(dst_manifold_ids[ind]);
3396 cell->face(f)->set_boundary_id(reset_boundary_ids[ind]);
3400 for (
unsigned int e=0; e<GeometryInfo<dim>::lines_per_face; ++e)
3402 const auto bid = cell->face(f)->line(e)->boundary_id();
3403 const auto ind = std::find(src_boundary_ids.begin(), src_boundary_ids.end(), bid)-
3404 src_boundary_ids.begin();
3405 if ((
unsigned int)ind < src_boundary_ids.size())
3406 cell->face(f)->line(e)->set_boundary_id(reset_boundary_ids[ind]);
3412 template <
int dim,
int spacedim>
3414 const bool compute_face_ids)
3419 for (; cell != endc; ++cell)
3421 cell->set_manifold_id(cell->material_id());
3422 if (compute_face_ids ==
true)
3424 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
3426 if (cell->at_boundary(f) ==
false)
3427 cell->face(f)->set_manifold_id
3428 (std::min(cell->material_id(),
3429 cell->neighbor(f)->material_id()));
3431 cell->face(f)->set_manifold_id(cell->material_id());
3437 template <
int dim,
int spacedim>
3438 std::pair<unsigned int, double>
3441 double max_ratio = 1;
3442 unsigned int index = 0;
3444 for (
unsigned int i = 0; i < dim; ++i)
3445 for (
unsigned int j = i+1; j < dim; ++j)
3447 unsigned int ax = i % dim;
3448 unsigned int next_ax = j % dim;
3450 double ratio = cell->extent_in_direction(ax) / cell->extent_in_direction(next_ax);
3452 if ( ratio > max_ratio )
3457 else if ( 1.0 /ratio > max_ratio )
3459 max_ratio = 1.0 /ratio;
3463 return std::make_pair(index, max_ratio);
3467 template <
int dim,
int spacedim>
3470 const bool isotropic,
3471 const unsigned int max_iterations)
3473 unsigned int iter = 0;
3474 bool continue_refinement =
true;
3480 while ( continue_refinement && (iter < max_iterations) )
3483 continue_refinement =
false;
3486 for (
unsigned int j = 0; j < GeometryInfo<dim>::faces_per_cell; j++)
3487 if (cell->at_boundary(j)==
false && cell->neighbor(j)->has_children())
3491 cell->set_refine_flag();
3492 continue_refinement =
true;
3495 continue_refinement |= cell->flag_for_face_refinement(j);
3502 template <
int dim,
int spacedim>
3505 const double max_ratio,
3506 const unsigned int max_iterations)
3508 unsigned int iter = 0;
3509 bool continue_refinement =
true;
3515 while ( continue_refinement && (iter<max_iterations) )
3518 continue_refinement =
false;
3521 std::pair<unsigned int, double> info = GridTools::get_longest_direction<dim,spacedim>(cell);
3522 if (info.second > max_ratio)
3525 continue_refinement =
true;
3533 template <
int dim,
int spacedim>
3535 const double limit_angle_fraction)
3542 "have hanging nodes."));
3545 bool has_cells_with_more_than_dim_faces_on_boundary =
true;
3546 bool has_cells_with_dim_faces_on_boundary =
false;
3548 unsigned int refinement_cycles = 0;
3550 while (has_cells_with_more_than_dim_faces_on_boundary)
3552 has_cells_with_more_than_dim_faces_on_boundary =
false;
3556 unsigned int boundary_face_counter = 0;
3557 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
3558 if (cell->face(f)->at_boundary())
3559 boundary_face_counter++;
3560 if (boundary_face_counter > dim)
3562 has_cells_with_more_than_dim_faces_on_boundary =
true;
3565 else if (boundary_face_counter == dim)
3566 has_cells_with_dim_faces_on_boundary =
true;
3568 if (has_cells_with_more_than_dim_faces_on_boundary)
3571 refinement_cycles++;
3575 if (has_cells_with_dim_faces_on_boundary)
3578 refinement_cycles++;
3582 while (refinement_cycles>0)
3585 cell->set_coarsen_flag();
3587 refinement_cycles--;
3593 std::vector<Point<spacedim> > vertices = tria.
get_vertices();
3595 std::vector<bool> faces_to_remove(tria.
n_raw_faces(),
false);
3597 std::vector<CellData<dim> > cells_to_add;
3603 v2 = (dim > 1 ? 2:0), v3 = (dim > 1 ? 3:0);
3607 double angle_fraction = 0;
3613 p0[spacedim > 1 ? 1 : 0] = 1;
3617 if (cell->face(v0)->at_boundary() && cell->face(v3)->at_boundary())
3619 p0 = cell->vertex(v0) - cell->vertex(v2);
3620 p1 = cell->vertex(v3) - cell->vertex(v2);
3621 vertex_at_corner = v2;
3623 else if (cell->face(v3)->at_boundary() && cell->face(v1)->at_boundary())
3625 p0 = cell->vertex(v2) - cell->vertex(v3);
3626 p1 = cell->vertex(v1) - cell->vertex(v3);
3627 vertex_at_corner = v3;
3629 else if (cell->face(1)->at_boundary() && cell->face(2)->at_boundary())
3631 p0 = cell->vertex(v0) - cell->vertex(v1);
3632 p1 = cell->vertex(v3) - cell->vertex(v1);
3633 vertex_at_corner = v1;
3635 else if (cell->face(2)->at_boundary() && cell->face(0)->at_boundary())
3637 p0 = cell->vertex(v2) - cell->vertex(v0);
3638 p1 = cell->vertex(v1) - cell->vertex(v0);
3639 vertex_at_corner = v0;
3651 if (angle_fraction > limit_angle_fraction)
3654 auto flags_removal = [&](
unsigned int f1,
unsigned int f2,
3655 unsigned int n1,
unsigned int n2) ->
void 3657 cells_to_remove[cell->active_cell_index()] =
true;
3658 cells_to_remove[cell->neighbor(n1)->active_cell_index()] =
true;
3659 cells_to_remove[cell->neighbor(n2)->active_cell_index()] =
true;
3661 faces_to_remove[cell->face(f1)->index()] =
true;
3662 faces_to_remove[cell->face(f2)->index()] =
true;
3664 faces_to_remove[cell->neighbor(n1)->face(f1)->index()] =
true;
3665 faces_to_remove[cell->neighbor(n2)->face(f2)->index()] =
true;
3668 auto cell_creation = [&](
3669 const unsigned int vv0,
3670 const unsigned int vv1,
3671 const unsigned int f0,
3672 const unsigned int f1,
3674 const unsigned int n0,
3675 const unsigned int v0n0,
3676 const unsigned int v1n0,
3678 const unsigned int n1,
3679 const unsigned int v0n1,
3680 const unsigned int v1n1)
3685 c1.
vertices[v0] = cell->vertex_index(vv0);
3686 c1.
vertices[v1] = cell->vertex_index(vv1);
3687 c1.
vertices[v2] = cell->neighbor(n0)->vertex_index(v0n0);
3688 c1.
vertices[v3] = cell->neighbor(n0)->vertex_index(v1n0);
3693 c2.
vertices[v0] = cell->vertex_index(vv0);
3694 c2.
vertices[v1] = cell->neighbor(n1)->vertex_index(v0n1);
3695 c2.
vertices[v2] = cell->vertex_index(vv1);
3696 c2.
vertices[v3] = cell->neighbor(n1)->vertex_index(v1n1);
3701 l1.
vertices[0] = cell->vertex_index(vv0);
3702 l1.
vertices[1] = cell->neighbor(n0)->vertex_index(v0n0);
3708 l2.
vertices[0] = cell->vertex_index(vv0);
3709 l2.
vertices[1] = cell->neighbor(n1)->vertex_index(v0n1);
3715 cells_to_add.push_back(c1);
3716 cells_to_add.push_back(c2);
3721 switch (vertex_at_corner)
3724 flags_removal(0,2,3,1);
3725 cell_creation(0,3, 0,2, 3,2,3, 1,1,3);
3728 flags_removal(1,2,3,0);
3729 cell_creation(1,2, 2,1, 0,0,2, 3,3,2);
3732 flags_removal(3,0,1,2);
3733 cell_creation(2,1, 3,0, 1,3,1, 2,0,1);
3736 flags_removal(3,1,0,2);
3737 cell_creation(3,0, 1,3, 2,1,0, 0,2,0);
3750 if (cells_to_add.size() == 0)
3752 while (refinement_cycles>0)
3755 cell->set_coarsen_flag();
3757 refinement_cycles--;
3765 if (cells_to_remove[cell->active_cell_index()] ==
false)
3768 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
3769 c.
vertices[v] = cell->vertex_index(v);
3772 cells_to_add.push_back(c);
3780 for (; face != endf; ++face)
3782 && faces_to_remove[face->index()] ==
false)
3784 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_face; ++l)
3789 for (
unsigned int v=0; v<GeometryInfo<1>::vertices_per_cell; ++v)
3790 line.
vertices[v] = face->vertex_index(v);
3796 for (
unsigned int v=0; v<GeometryInfo<1>::vertices_per_cell; ++v)
3797 line.
vertices[v] = face->line(l)->vertex_index(v);
3806 for (
unsigned int v=0; v<GeometryInfo<2>::vertices_per_cell; ++v)
3807 quad.
vertices[v] = face->vertex_index(v);
3818 std::map<types::manifold_id, std::unique_ptr<Manifold<dim,spacedim> > > manifolds;
3820 for (
auto manifold_id: manifold_ids)
3822 manifolds[manifold_id] = tria.
get_manifold(manifold_id).clone();
3829 for (
auto manifold_id: manifold_ids)
3831 tria.
set_manifold(manifold_id, *manifolds[manifold_id]);
3836 template <
int dim,
int spacedim>
3838 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator >,
3839 std::vector< std::vector< Point<dim> > >,
3840 std::vector< std::vector<unsigned int> > >
3846 const unsigned int np = points.size();
3849 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator >,
3850 std::vector< std::vector< Point<dim> > >,
3851 std::vector< std::vector<unsigned int> > >
3855 if (np==0)
return cell_qpoint_map;
3858 std::pair<typename Triangulation<dim, spacedim>::active_cell_iterator,
Point<dim> >
3862 (cache, points[0],cell_hint);
3867 std::get<0>(cell_qpoint_map).emplace_back(my_pair.first);
3868 std::get<1>(cell_qpoint_map).emplace_back(1, my_pair.second);
3869 std::get<2>(cell_qpoint_map).emplace_back(1, 0);
3872 if (np==1)
return cell_qpoint_map;
3874 Point<spacedim> cell_center = std::get<0>(cell_qpoint_map)[0]->center();
3875 double cell_diameter = std::get<0>(cell_qpoint_map)[0]->
diameter()*
3876 (0.5 + std::numeric_limits<double>::epsilon() );
3879 for (
unsigned int p=1; p< np; ++p)
3883 if ( cell_center.
distance(points[p]) < cell_diameter )
3885 (cache, points[p],std::get<0>(cell_qpoint_map).back());
3891 if ( my_pair.first == std::get<0>(cell_qpoint_map).back() )
3894 std::get<1>(cell_qpoint_map).back().emplace_back(my_pair.second);
3895 std::get<2>(cell_qpoint_map).back().emplace_back(p);
3900 typename std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>::iterator
3901 cells_it = std::find(std::get<0>(cell_qpoint_map).begin(),std::get<0>(cell_qpoint_map).end()-1,my_pair.first);
3903 if ( cells_it == std::get<0>(cell_qpoint_map).end()-1 )
3906 std::get<0>(cell_qpoint_map).emplace_back(my_pair.first);
3907 std::get<1>(cell_qpoint_map).emplace_back(1, my_pair.second);
3908 std::get<2>(cell_qpoint_map).emplace_back(1, p);
3910 cell_center = std::get<0>(cell_qpoint_map).back()->center();
3911 cell_diameter = std::get<0>(cell_qpoint_map).back()->diameter()*
3912 (0.5 + std::numeric_limits<double>::epsilon() );
3916 unsigned int current_cell = cells_it - std::get<0>(cell_qpoint_map).begin();
3918 std::get<1>(cell_qpoint_map)[current_cell].emplace_back(my_pair.second);
3919 std::get<2>(cell_qpoint_map)[current_cell].emplace_back(p);
3925 Assert(std::get<0>(cell_qpoint_map).size() == std::get<2>(cell_qpoint_map).size(),
3926 ExcDimensionMismatch(std::get<0>(cell_qpoint_map).size(), std::get<2>(cell_qpoint_map).size()));
3928 Assert(std::get<0>(cell_qpoint_map).size() == std::get<1>(cell_qpoint_map).size(),
3929 ExcDimensionMismatch(std::get<0>(cell_qpoint_map).size(), std::get<1>(cell_qpoint_map).size()));
3932 unsigned int c = std::get<0>(cell_qpoint_map).size();
3933 unsigned int qps = 0;
3938 for (
unsigned int n=0; n<c; ++n)
3940 Assert(std::get<1>(cell_qpoint_map)[n].size() ==
3941 std::get<2>(cell_qpoint_map)[n].size(),
3943 std::get<2>(cell_qpoint_map)[n].size()));
3944 qps += std::get<1>(cell_qpoint_map)[n].size();
3950 return cell_qpoint_map;
3958 namespace distributed_cptloc
3961 template <
int dim,
int spacedim>
3967 return k->active_cell_index();
3975 template <
int dim,
int spacedim>
3976 std::unordered_map< typename Triangulation<dim, spacedim>::active_cell_iterator,
3977 std::pair<std::vector<Point<dim> >,std::vector<unsigned int> >, cell_hash<dim,spacedim> >
3982 const unsigned int np = points.size();
3984 std::unordered_map< typename Triangulation<dim, spacedim>::active_cell_iterator,
3985 std::pair<std::vector<Point<dim> >,std::vector<unsigned int> >, cell_hash<dim,spacedim> >
3989 if (np==0)
return cell_qpoint_map;
3994 auto last_cell = cell_qpoint_map.emplace(
3995 std::make_pair(my_pair.first, std::make_pair(
3997 std::vector<unsigned int> {0})));
3999 if (np==1)
return cell_qpoint_map;
4002 double cell_diameter = my_pair.first->diameter()*
4003 (0.5 + std::numeric_limits<double>::epsilon() );
4006 for (
unsigned int p=1; p< np; ++p)
4010 if ( cell_center.
distance(points[p]) < cell_diameter )
4012 (cache, points[p],last_cell.first->first);
4017 if ( last_cell.first->first == my_pair.first)
4019 last_cell.first->second.first.emplace_back(my_pair.second);
4020 last_cell.first->second.second.emplace_back(p);
4025 last_cell = cell_qpoint_map.emplace(std::make_pair(my_pair.first, std::make_pair(
4027 std::vector<unsigned int> {p})));
4029 if ( last_cell.second ==
false )
4032 last_cell.first->second.first.emplace_back(my_pair.second);
4033 last_cell.first->second.second.emplace_back(p);
4038 cell_center = my_pair.first->center();
4039 cell_diameter = my_pair.first->diameter()*
4040 (0.5 + std::numeric_limits<double>::epsilon() );
4046 unsigned int qps = 0;
4051 for (
const auto &m: cell_qpoint_map)
4053 Assert(m.second.second.size() ==
4054 m.second.first.size(),
4056 m.second.first.size()));
4057 qps += m.second.second.size();
4062 return cell_qpoint_map;
4075 template <
int dim,
int spacedim>
4077 merge_cptloc_outputs(
4081 std::vector< unsigned int >,
4083 std::vector< unsigned int >
4085 cell_hash<dim,spacedim>> &output_unmap,
4087 const std::vector< std::vector<
Point<dim> > > &in_qpoints,
4088 const std::vector< std::vector<unsigned int> > &in_maps,
4090 const unsigned int in_rank
4094 for (
unsigned int c=0; c< in_cells.size(); ++c)
4097 auto current_c = output_unmap.emplace(
4098 std::make_pair(in_cells[c],
4099 std::make_tuple(in_qpoints[c],
4102 std::vector<unsigned int>
4103 (in_points[c].size(),in_rank))));
4105 if ( current_c.second ==
false )
4109 auto &cell_qpts = std::get<0>(current_c.first->second);
4110 auto &cell_maps = std::get<1>(current_c.first->second);
4111 auto &cell_pts = std::get<2>(current_c.first->second);
4112 auto &cell_ranks = std::get<3>(current_c.first->second);
4113 cell_qpts.insert(cell_qpts.end(),
4114 in_qpoints[c].begin(),
4115 in_qpoints[c].end());
4116 cell_maps.insert(cell_maps.end(),
4119 cell_pts.insert(cell_pts.end(),
4120 in_points[c].begin(),
4121 in_points[c].end());
4122 std::vector< unsigned int > ranks_tmp(in_points[c].size(),in_rank);
4123 cell_ranks.insert(cell_ranks.end(),
4139 template <
int dim,
int spacedim>
4141 compute_and_classify_points(
4144 const std::vector< unsigned int > &local_points_idx,
4149 std::vector< unsigned int >,
4151 std::vector< unsigned int >
4153 cell_hash<dim,spacedim>> &output_unmap,
4154 std::map<
unsigned int,
4156 std::vector< CellId >,
4158 std::vector< std::vector< unsigned int > >,
4161 std::vector< unsigned int > &classified_pts
4164 auto cpt_loc_pts = compute_point_locations_unmap(cache,local_points);
4168 for (
auto const &cell_tuples : cpt_loc_pts)
4170 auto &cell_loc = cell_tuples.first;
4171 auto &q_loc = std::get<0>(cell_tuples.second);
4172 auto &indices_loc = std::get<1>(cell_tuples.second);
4173 if (cell_loc->is_locally_owned() )
4176 std::vector < Point<spacedim> > cell_points(indices_loc.size());
4177 std::vector < unsigned int > cell_points_idx(indices_loc.size());
4178 for (
unsigned int i=0; i< indices_loc.size(); ++i)
4181 cell_points[i] = local_points[indices_loc[i]];
4186 cell_points_idx[i] = local_points_idx[indices_loc[i]];
4187 classified_pts.emplace_back(local_points_idx[indices_loc[i]]);
4189 output_unmap.emplace(std::make_pair(cell_loc,
4190 std::make_tuple(q_loc,
4193 std::vector<unsigned int>
4194 (indices_loc.size(),cell_loc->subdomain_id()))));
4196 else if ( cell_loc->is_ghost() )
4200 std::vector < Point<spacedim> > cell_points(indices_loc.size());
4201 std::vector < unsigned int > cell_points_idx(indices_loc.size());
4202 for (
unsigned int i=0; i< indices_loc.size(); ++i)
4204 cell_points[i] = local_points[indices_loc[i]];
4205 cell_points_idx[i] = local_points_idx[indices_loc[i]];
4206 classified_pts.emplace_back(local_points_idx[indices_loc[i]]);
4211 auto &map_tuple_owner = ghost_loc_pts[cell_loc->subdomain_id()];
4213 std::get<0>(map_tuple_owner).emplace_back(cell_loc->id());
4214 std::get<1>(map_tuple_owner).emplace_back(q_loc);
4215 std::get<2>(map_tuple_owner).emplace_back(cell_points_idx);
4216 std::get<3>(map_tuple_owner).emplace_back(cell_points);
4229 template <
int dim,
int spacedim>
4231 compute_and_merge_from_map(
4233 const std::map<
unsigned int,
4236 std::vector < unsigned int > >
4241 std::vector< unsigned int >,
4243 std::vector< unsigned int >
4245 cell_hash<dim,spacedim>> &output_unmap,
4246 const bool &check_owned
4249 bool no_check = !check_owned;
4252 for (
auto const &rank_and_points : map_pts)
4255 const auto &received_process = rank_and_points.first;
4256 const auto &received_points = rank_and_points.second.first;
4257 const auto &received_map = rank_and_points.second.second;
4260 std::vector< typename Triangulation<dim, spacedim>::active_cell_iterator > in_cell;
4261 std::vector< std::vector< Point<dim> > > in_qpoints;
4262 std::vector< std::vector< unsigned int > > in_maps;
4263 std::vector< std::vector< Point<spacedim> > > in_points;
4265 auto cpt_loc_pts = compute_point_locations_unmap(cache,rank_and_points.second.first);
4266 for (
const auto &map_c_pt_idx: cpt_loc_pts)
4269 const auto &proc_cell = map_c_pt_idx.first;
4270 const auto &proc_qpoints = map_c_pt_idx.second.first;
4271 const auto &proc_maps = map_c_pt_idx.second.second;
4275 if ( no_check || proc_cell->is_locally_owned() )
4277 in_cell.emplace_back(proc_cell);
4278 in_qpoints.emplace_back(proc_qpoints);
4280 unsigned int loc_size = proc_qpoints.size();
4281 std::vector< unsigned int > cell_maps(loc_size);
4282 std::vector< Point<spacedim> > cell_points(loc_size);
4283 for (
unsigned int pt=0; pt<loc_size; ++pt)
4285 cell_maps[pt] = received_map[proc_maps[pt]];
4286 cell_points[pt] = received_points[proc_maps[pt]];
4288 in_maps.emplace_back(cell_maps);
4289 in_points.emplace_back(cell_points);
4294 internal::distributed_cptloc::merge_cptloc_outputs(output_unmap,
4307 template <
int dim,
int spacedim>
4309 std::vector< typename Triangulation<dim, spacedim>::active_cell_iterator >,
4310 std::vector< std::vector< Point<dim> > >,
4311 std::vector< std::vector< unsigned int > >,
4312 std::vector< std::vector< Point<spacedim> > >,
4313 std::vector< std::vector< unsigned int > >
4320 #ifndef DEAL_II_WITH_MPI 4323 (void)global_bboxes;
4324 Assert(
false,
ExcMessage(
"GridTools::distributed_compute_point_locations() requires MPI."));
4326 std::vector< typename Triangulation<dim, spacedim>::active_cell_iterator >,
4327 std::vector< std::vector< Point<dim> > >,
4328 std::vector< std::vector< unsigned int > >,
4329 std::vector< std::vector< Point<spacedim> > >,
4330 std::vector< std::vector< unsigned int > >
4335 const auto &tria_mpi =
4338 Assert(tria_mpi,
ExcMessage(
"GridTools::distributed_compute_point_locations() requires a parallel triangulation."));
4339 auto mpi_communicator = tria_mpi->get_communicator();
4342 std::vector< typename Triangulation<dim, spacedim>::active_cell_iterator >,
4343 std::vector< std::vector< Point<dim> > >,
4344 std::vector< std::vector< unsigned int > >,
4345 std::vector< std::vector< Point<spacedim> > >,
4346 std::vector< std::vector< unsigned int > >
4350 std::unordered_map< typename Triangulation<dim, spacedim>::active_cell_iterator,
4352 std::vector< Point<dim> >,
4353 std::vector< unsigned int >,
4354 std::vector< Point<spacedim> >,
4355 std::vector< unsigned int >
4357 internal::distributed_cptloc::cell_hash<dim,spacedim> >
4365 std::tuple< std::vector< std::vector< unsigned int > >, std::map< unsigned int, unsigned int >,
4366 std::map< unsigned int, std::vector< unsigned int > > > guessed_points;
4372 const auto &guess_loc_idx = std::get<0>(guessed_points)[my_rank];
4373 const unsigned int n_local_guess = guess_loc_idx.size();
4375 std::vector< Point<spacedim> > guess_local_pts(n_local_guess);
4376 for (
unsigned int i=0; i<n_local_guess; ++i)
4377 guess_local_pts[i] = local_points[ guess_loc_idx[i] ];
4380 std::map<
unsigned int,
4382 std::vector< CellId >,
4383 std::vector< std::vector< Point<dim> > >,
4384 std::vector< std::vector< unsigned int > >,
4385 std::vector< std::vector< Point<spacedim> > > > > ghost_loc_pts;
4388 std::vector< unsigned int > classified_pts;
4394 &internal::distributed_cptloc::compute_and_classify_points<dim,spacedim>,
4404 const auto &other_owned_idx = std::get<1>(guessed_points);
4407 std::pair< std::vector<Point<spacedim>> , std::vector<unsigned int > > >
4410 for (
const auto &indices: other_owned_idx)
4411 if (indices.second != my_rank)
4414 auto ¤t_pts = other_owned_pts[indices.second];
4416 current_pts.first.emplace_back(local_points[indices.first]);
4417 current_pts.second.emplace_back(indices.first);
4428 =
Threads::new_task (&internal::distributed_cptloc::compute_and_merge_from_map<dim,spacedim>,
4443 std::pair< std::vector < Point<spacedim> >,
4444 std::vector< unsigned int > > >
4449 const std::map< unsigned int, std::vector< unsigned int > >
4450 &other_check_idx = std::get<2>(guessed_points);
4457 std::sort (classified_pts.begin(), classified_pts.end());
4459 for (
const auto &pt_to_guesses: other_check_idx)
4461 const auto &point_idx = pt_to_guesses.first;
4462 const auto &probable_owners_rks = pt_to_guesses.second;
4463 if ( !std::binary_search(
4464 classified_pts.begin(), classified_pts.end(), point_idx) )
4466 for (
unsigned int i=0; i<probable_owners_rks.size(); ++i)
4467 if (probable_owners_rks[i] != my_rank)
4470 auto ¤t_pts = other_check_pts[probable_owners_rks[i]];
4472 current_pts.first.emplace_back(local_points[point_idx]);
4474 current_pts.second.emplace_back(point_idx);
4481 owned_pts_tsk.
join();
4484 for (
const auto &rank_vals: cpt_ghost)
4487 const auto &cell_ids = std::get<0>(rank_vals.second);
4488 unsigned int n_cells = cell_ids.size();
4489 std::vector< typename Triangulation<dim, spacedim>::active_cell_iterator >
4491 for (
unsigned int c=0; c<n_cells; ++c)
4494 internal::distributed_cptloc::merge_cptloc_outputs(temporary_unmap,
4496 std::get<1>(rank_vals.second),
4497 std::get<2>(rank_vals.second),
4498 std::get<3>(rank_vals.second),
4504 internal::distributed_cptloc::compute_and_merge_from_map(
4512 unsigned int size_output = temporary_unmap.size();
4513 auto &out_cells = std::get<0>(output_tuple);
4514 auto &out_qpoints = std::get<1>(output_tuple);
4515 auto &out_maps = std::get<2>(output_tuple);
4516 auto &out_points = std::get<3>(output_tuple);
4517 auto &out_ranks = std::get<4>(output_tuple);
4519 out_cells.resize(size_output);
4520 out_qpoints.resize(size_output);
4521 out_maps.resize(size_output);
4522 out_points.resize(size_output);
4523 out_ranks.resize(size_output);
4526 for (
const auto &rank_and_tuple: temporary_unmap)
4528 out_cells[c] = rank_and_tuple.first;
4529 out_qpoints[c] = std::get<0>(rank_and_tuple.second);
4530 out_maps[c] = std::get<1>(rank_and_tuple.second);
4531 out_points[c] = std::get<2>(rank_and_tuple.second);
4532 out_ranks[c] = std::get<3>(rank_and_tuple.second);
4536 return output_tuple;
4541 template<
int dim,
int spacedim>
4542 std::map<unsigned int, Point<spacedim> >
4546 std::map<unsigned int, Point<spacedim> > result;
4550 for (
unsigned int i=0; i<vs.size(); ++i)
4551 result[cell->vertex_index(i)]=vs[i];
4559 template<
int spacedim>
4564 auto id_and_v = std::min_element(vertices.begin(), vertices.end(),
4565 [&](
const std::pair<const unsigned int, Point<spacedim>> &p1,
4566 const std::pair<const unsigned int, Point<spacedim>> &p2) ->
bool 4568 return p1.second.
distance(p) < p2.second.distance(p);
4571 return id_and_v->first;
4575 template<
int dim,
int spacedim>
4576 std::pair<typename Triangulation<dim,spacedim>::active_cell_iterator,
Point<dim> >
4580 const std::vector<bool> &marked_vertices)
4589 vertex_to_cell_centers,
4594 template<
int spacedim>
4595 std::vector< std::vector< BoundingBox<spacedim> > >
4597 MPI_Comm mpi_communicator)
4599 #ifndef DEAL_II_WITH_MPI 4601 (void)mpi_communicator;
4602 Assert(
false,
ExcMessage(
"GridTools::exchange_local_bounding_boxes() requires MPI."));
4606 unsigned int n_bboxes = local_bboxes.size();
4608 int n_local_data = 2*spacedim*n_bboxes;
4610 std::vector<double> loc_data_array(n_local_data);
4611 for (
unsigned int i=0; i<n_bboxes; ++i)
4612 for (
unsigned int d=0;
d < spacedim; ++
d)
4615 loc_data_array[2*i*spacedim +
d] = local_bboxes[i].get_boundary_points().first[
d];
4616 loc_data_array[2*i*spacedim + spacedim +
d] = local_bboxes[i].get_boundary_points().second[
d];
4623 std::vector<int> size_all_data(n_procs);
4626 int ierr = MPI_Allgather(&n_local_data, 1, MPI_INT,
4627 &(size_all_data[0]), 1, MPI_INT,
4633 std::vector<int> rdispls(n_procs);
4635 for (
unsigned int i=1; i < n_procs; ++i)
4636 rdispls[i] = rdispls[i-1] + size_all_data[i-1];
4640 std::vector<double> data_array(rdispls.back() + size_all_data.back());
4642 ierr = MPI_Allgatherv(&(loc_data_array[0]), n_local_data, MPI_DOUBLE,
4643 &(data_array[0]), &(size_all_data[0]),
4644 &(rdispls[0]), MPI_DOUBLE, mpi_communicator);
4648 std::vector< std::vector< BoundingBox<spacedim> > > global_bboxes(n_procs);
4649 unsigned int begin_idx = 0;
4650 for (
unsigned int i=0; i < n_procs; ++i)
4653 unsigned int n_bbox_i = size_all_data[i]/(spacedim*2);
4654 global_bboxes[i].resize(n_bbox_i);
4655 for (
unsigned int bbox=0; bbox<n_bbox_i; ++bbox)
4658 for (
unsigned int d=0;
d<spacedim; ++
d)
4660 p1[
d] = data_array[ begin_idx + 2*bbox*spacedim +
d];
4661 p2[
d] = data_array[ begin_idx + 2*bbox*spacedim + spacedim +
d];
4664 global_bboxes[i][bbox] = loc_bbox;
4667 begin_idx += size_all_data[i];
4669 return global_bboxes;
4670 #endif // DEAL_II_WITH_MPI 4678 #include "grid_tools.inst" 4680 DEAL_II_NAMESPACE_CLOSE
void map_boundary_to_manifold_ids(const std::vector< types::boundary_id > &src_boundary_ids, const std::vector< types::manifold_id > &dst_manifold_ids, Triangulation< dim, spacedim > &tria, const std::vector< types::boundary_id > &reset_boundary_ids={})
std::vector< CellData< 1 > > boundary_lines
Transformed quadrature weights.
static ::ExceptionBase & ExcScalingFactorNotPositive(double arg1)
unsigned int n_active_cells() const
unsigned int n_used_vertices() const
static const unsigned int invalid_unsigned_int
#define AssertDimension(dim1, dim2)
void copy_boundary_to_manifold_id(Triangulation< dim, spacedim > &tria, const bool reset_boundary_ids=false)
active_face_iterator begin_active_face() const
cell_iterator begin(const unsigned int level=0) const
virtual bool has_hanging_nodes() const
unsigned int n_cells() const
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
Task< RT > new_task(const std::function< RT()> &function)
std::list< typename Triangulation< dim, spacedim >::cell_iterator > distorted_cells
void add(const size_type i, const size_type j)
IteratorRange< active_cell_iterator > active_cell_iterators() const
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
active_cell_iterator begin_active(const unsigned int level=0) const
void create_laplace_matrix(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, SparseMatrix< double > &matrix, const Function< spacedim > *const a=nullptr, const ConstraintMatrix &constraints=ConstraintMatrix())
#define AssertThrow(cond, exc)
numbers::NumberTraits< Number >::real_type norm() const
types::boundary_id boundary_id
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const =0
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
void set_manifold(const types::manifold_id number, const Manifold< dim, spacedim > &manifold_object)
static double distance_to_unit_cell(const Point< dim > &p)
cell_iterator end() const
active_cell_iterator begin_active(const unsigned int level=0) const
virtual void execute_coarsening_and_refinement()
static ::ExceptionBase & ExcInvalidNumberOfPartitions(int arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
bool check_consistency(const unsigned int dim) const
T sum(const T &t, const MPI_Comm &mpi_communicator)
unsigned int global_dof_index
#define Assert(cond, exc)
void reinit(const size_type m, const size_type n, const IndexSet &rowset=IndexSet())
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
types::material_id material_id
const std::vector< Point< spacedim > > & get_vertices() const
void initialize(const MatrixType &A, const AdditionalData ¶meters=AdditionalData())
static void reorder_cells(std::vector< CellData< dim > > &original_cells, const bool use_new_style_ordering=false)
unsigned long long int global_vertex_index
void copy_material_to_manifold_id(Triangulation< dim, spacedim > &tria, const bool compute_face_ids=false)
void copy_from(const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end)
std::vector< unsigned int > invert_permutation(const std::vector< unsigned int > &permutation)
unsigned int subdomain_id
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
unsigned int size() const
types::manifold_id manifold_id
unsigned int n_raw_faces() const
void swap(Vector< Number > &u, Vector< Number > &v)
static void alternating_form_at_vertices(const Point< spacedim >(&vertices)[vertices_per_cell], Tensor< spacedim-dim, spacedim >(&forms)[vertices_per_cell])
const types::subdomain_id artificial_subdomain_id
std::vector< types::boundary_id > get_boundary_ids() const
const types::manifold_id invalid_manifold_id
#define AssertThrowMPI(error_code)
void make_sparsity_pattern(const DoFHandlerType &dof_handler, SparsityPatternType &sparsity_pattern, const ConstraintMatrix &constraints=ConstraintMatrix(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
double JxW(const unsigned int quadrature_point) const
const std::vector< bool > & get_used_vertices() const
T min(const T &t, const MPI_Comm &mpi_communicator)
std::vector< CellData< 2 > > boundary_quads
numbers::NumberTraits< Number >::real_type square() const
void refine_global(const unsigned int times=1)
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
virtual bool preserves_vertex_locations() const =0
static ::ExceptionBase & ExcNotImplemented()
Iterator points to a valid object.
static ::ExceptionBase & ExcVertexNotUsed(unsigned int arg1)
face_iterator end_face() const
IteratorState::IteratorStates state() const
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access > > &cell)
std::map< unsigned int, T > some_to_some(const MPI_Comm &comm, const std::map< unsigned int, T > &objects_to_send)
std::vector< types::manifold_id > get_manifold_ids() const
bool vertex_used(const unsigned int index) const
const Manifold< dim, spacedim > & get_manifold(const types::manifold_id number) const
T max(const T &t, const MPI_Comm &mpi_communicator)
numbers::NumberTraits< Number >::real_type distance_square(const Point< dim, Number > &p) const
unsigned int vertices[GeometryInfo< structdim >::vertices_per_cell]
static ::ExceptionBase & ExcInternalError()
Triangulation< dim, spacedim > & get_triangulation()