16 #include <deal.II/base/memory_consumption.h> 17 #include <deal.II/dofs/dof_handler.h> 18 #include <deal.II/dofs/dof_handler_policy.h> 19 #include <deal.II/dofs/dof_levels.h> 20 #include <deal.II/dofs/dof_faces.h> 21 #include <deal.II/dofs/dof_accessor.h> 22 #include <deal.II/grid/tria_accessor.h> 23 #include <deal.II/grid/tria_iterator.h> 24 #include <deal.II/grid/tria_levels.h> 25 #include <deal.II/grid/tria.h> 26 #include <deal.II/base/geometry_info.h> 27 #include <deal.II/fe/fe.h> 28 #include <deal.II/distributed/shared_tria.h> 29 #include <deal.II/distributed/tria.h> 34 DEAL_II_NAMESPACE_OPEN
41 template<
int dim,
int spacedim>
44 template<
int dim,
int spacedim>
47 template <
int dim,
int spacedim>
50 template <
int dim,
int spacedim>
58 template <
int dim,
int spacedim>
69 template<
int dim,
int spacedim>
70 std::string policy_to_string(const ::internal::DoFHandler::Policy::PolicyBase<dim,spacedim> &policy)
72 std::string policy_name;
73 if (
dynamic_cast<const typename ::internal::DoFHandler::Policy::Sequential<dim,spacedim>*
>(&policy))
74 policy_name =
"Policy::Sequential<";
75 else if (
dynamic_cast<const typename ::internal::DoFHandler::Policy::ParallelDistributed<dim,spacedim>*
>(&policy))
76 policy_name =
"Policy::ParallelDistributed<";
77 else if (
dynamic_cast<const typename ::internal::DoFHandler::Policy::ParallelShared<dim,spacedim>*
>(&policy))
78 policy_name =
"Policy::ParallelShared<";
105 template <
int spacedim>
110 return std::min(static_cast<types::global_dof_index>(3*dof_handler.
selected_fe->dofs_per_vertex +
117 template <
int spacedim>
145 switch (dof_handler.
tria->max_adjacent_cells())
148 max_couplings=19*dof_handler.
selected_fe->dofs_per_vertex +
153 max_couplings=21*dof_handler.
selected_fe->dofs_per_vertex +
158 max_couplings=28*dof_handler.
selected_fe->dofs_per_vertex +
163 max_couplings=30*dof_handler.
selected_fe->dofs_per_vertex +
168 max_couplings=37*dof_handler.
selected_fe->dofs_per_vertex +
178 max_couplings=39*dof_handler.
selected_fe->dofs_per_vertex +
183 max_couplings=46*dof_handler.
selected_fe->dofs_per_vertex +
188 max_couplings=48*dof_handler.
selected_fe->dofs_per_vertex +
193 max_couplings=55*dof_handler.
selected_fe->dofs_per_vertex +
198 max_couplings=57*dof_handler.
selected_fe->dofs_per_vertex +
203 max_couplings=63*dof_handler.
selected_fe->dofs_per_vertex +
208 max_couplings=65*dof_handler.
selected_fe->dofs_per_vertex +
213 max_couplings=72*dof_handler.
selected_fe->dofs_per_vertex +
222 return std::min(max_couplings,dof_handler.
n_dofs());
226 template <
int spacedim>
244 const unsigned int max_adjacent_cells
245 = dof_handler.
tria->max_adjacent_cells();
248 if (max_adjacent_cells <= 8)
249 max_couplings=7*7*7*dof_handler.
selected_fe->dofs_per_vertex +
259 return std::min(max_couplings,dof_handler.
n_dofs());
272 template <
int spacedim>
277 .resize(dof_handler.
tria->n_vertices() *
281 for (
unsigned int i=0; i<dof_handler.
tria->n_levels(); ++i)
286 dof_handler.
levels.back()->dof_object.dofs
287 .resize (dof_handler.
tria->n_raw_cells(i) *
291 dof_handler.
levels.back()->cell_dof_indices_cache
292 .resize (dof_handler.
tria->n_raw_cells(i) *
299 template <
int spacedim>
304 .resize(dof_handler.
tria->n_vertices() *
308 for (
unsigned int i=0; i<dof_handler.
tria->n_levels(); ++i)
312 dof_handler.
levels.back()->dof_object.dofs
313 .resize (dof_handler.
tria->n_raw_cells(i) *
317 dof_handler.
levels.back()->cell_dof_indices_cache
318 .resize (dof_handler.
tria->n_raw_cells(i) *
325 if (dof_handler.
tria->n_cells() > 0)
327 dof_handler.
faces->lines.dofs
328 .resize (dof_handler.
tria->n_raw_lines() *
335 template <
int spacedim>
340 .resize(dof_handler.
tria->n_vertices() *
344 for (
unsigned int i=0; i<dof_handler.
tria->n_levels(); ++i)
348 dof_handler.
levels.back()->dof_object.dofs
349 .resize (dof_handler.
tria->n_raw_cells(i) *
353 dof_handler.
levels.back()->cell_dof_indices_cache
354 .resize (dof_handler.
tria->n_raw_cells(i) *
361 if (dof_handler.
tria->n_cells() > 0)
363 dof_handler.
faces->lines.dofs
364 .resize (dof_handler.
tria->n_raw_lines() *
367 dof_handler.
faces->quads.dofs
368 .resize (dof_handler.
tria->n_raw_quads() *
374 template<
int spacedim>
379 dof_handler.clear_mg_space ();
382 const unsigned int &dofs_per_line = dof_handler.
get_fe ().dofs_per_line;
383 const unsigned int &n_levels = tria.n_levels ();
385 for (
unsigned int i = 0; i < n_levels; ++i)
388 dof_handler.mg_levels.back ()->dof_object.dofs = std::vector<types::global_dof_index> (tria.n_raw_lines (i) * dofs_per_line,
DoFHandler<1>::invalid_dof_index);
391 const unsigned int &n_vertices = tria.n_vertices ();
395 std::vector<unsigned int> max_level (n_vertices, 0);
396 std::vector<unsigned int> min_level (n_vertices, n_levels);
398 for (typename ::Triangulation<1, spacedim>::cell_iterator cell = tria.begin (); cell != tria.end (); ++cell)
400 const unsigned int level = cell->level ();
402 for (
unsigned int vertex = 0; vertex < GeometryInfo<1>::vertices_per_cell; ++vertex)
404 const unsigned int vertex_index = cell->vertex_index (vertex);
406 if (min_level[vertex_index] > level)
407 min_level[vertex_index] = level;
409 if (max_level[vertex_index] < level)
410 max_level[vertex_index] = level;
414 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
415 if (tria.vertex_used (vertex))
419 dof_handler.
mg_vertex_dofs[vertex].init (min_level[vertex], max_level[vertex], dof_handler.
get_fe ().dofs_per_vertex);
430 template<
int spacedim>
435 dof_handler.clear_mg_space ();
437 const ::FiniteElement<2, spacedim> &fe = dof_handler.
get_fe ();
439 const unsigned int &n_levels = tria.n_levels ();
441 for (
unsigned int i = 0; i < n_levels; ++i)
444 dof_handler.mg_levels.back ()->dof_object.dofs = std::vector<types::global_dof_index> (tria.n_raw_quads (i) * fe.dofs_per_quad,
DoFHandler<2>::invalid_dof_index);
450 const unsigned int &n_vertices = tria.n_vertices ();
454 std::vector<unsigned int> max_level (n_vertices, 0);
455 std::vector<unsigned int> min_level (n_vertices, n_levels);
457 for (typename ::Triangulation<2, spacedim>::cell_iterator cell = tria.begin (); cell != tria.end (); ++cell)
459 const unsigned int level = cell->level ();
461 for (
unsigned int vertex = 0; vertex < GeometryInfo<2>::vertices_per_cell; ++vertex)
463 const unsigned int vertex_index = cell->vertex_index (vertex);
465 if (min_level[vertex_index] > level)
466 min_level[vertex_index] = level;
468 if (max_level[vertex_index] < level)
469 max_level[vertex_index] = level;
473 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
474 if (tria.vertex_used (vertex))
478 dof_handler.
mg_vertex_dofs[vertex].init (min_level[vertex], max_level[vertex], fe.dofs_per_vertex);
489 template<
int spacedim>
494 dof_handler.clear_mg_space ();
496 const ::FiniteElement<3, spacedim> &fe = dof_handler.
get_fe ();
498 const unsigned int &n_levels = tria.n_levels ();
500 for (
unsigned int i = 0; i < n_levels; ++i)
503 dof_handler.mg_levels.back ()->dof_object.dofs = std::vector<types::global_dof_index> (tria.n_raw_hexs (i) * fe.dofs_per_hex,
DoFHandler<3>::invalid_dof_index);
510 const unsigned int &n_vertices = tria.n_vertices ();
514 std::vector<unsigned int> max_level (n_vertices, 0);
515 std::vector<unsigned int> min_level (n_vertices, n_levels);
517 for (typename ::Triangulation<3, spacedim>::cell_iterator cell = tria.begin (); cell != tria.end (); ++cell)
519 const unsigned int level = cell->level ();
521 for (
unsigned int vertex = 0; vertex < GeometryInfo<3>::vertices_per_cell; ++vertex)
523 const unsigned int vertex_index = cell->vertex_index (vertex);
525 if (min_level[vertex_index] > level)
526 min_level[vertex_index] = level;
528 if (max_level[vertex_index] < level)
529 max_level[vertex_index] = level;
533 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
534 if (tria.vertex_used (vertex))
538 dof_handler.
mg_vertex_dofs[vertex].init (min_level[vertex], max_level[vertex], fe.dofs_per_vertex);
551 template<
int spacedim>
558 const unsigned int obj_index,
559 const unsigned int fe_index,
560 const unsigned int local_index,
563 return mg_level.
dof_object.get_dof_index (dof_handler, obj_index, fe_index, local_index);
566 template<
int spacedim>
574 template<
int spacedim>
579 return mg_level.
dof_object.get_dof_index (dof_handler, obj_index, fe_index, local_index);
582 template<
int spacedim>
590 template<
int spacedim>
598 template<
int spacedim>
603 return mg_level.
dof_object.get_dof_index (dof_handler, obj_index, fe_index, local_index);
606 template<
int spacedim>
610 mg_level.
dof_object.set_dof_index (dof_handler, obj_index, fe_index, local_index, global_index);
613 template<
int spacedim>
617 mg_faces.
lines.
set_dof_index (dof_handler, obj_index, fe_index, local_index, global_index);
620 template<
int spacedim>
624 mg_level.
dof_object.set_dof_index (dof_handler, obj_index, fe_index, local_index, global_index);
627 template<
int spacedim>
631 mg_faces.
lines.
set_dof_index (dof_handler, obj_index, fe_index, local_index, global_index);
634 template<
int spacedim>
638 mg_faces.
quads.
set_dof_index (dof_handler, obj_index, fe_index, local_index, global_index);
641 template<
int spacedim>
645 mg_level.
dof_object.set_dof_index (dof_handler, obj_index, fe_index, local_index, global_index);
653 template<
int dim,
int spacedim>
656 tria(&tria, typeid(*this).name()),
657 selected_fe(0, typeid(*this).name()),
677 template<
int dim,
int spacedim>
680 tria(0, typeid(*this).name()),
681 selected_fe(0, typeid(*this).name()),
687 template <
int dim,
int spacedim>
695 template<
int dim,
int spacedim>
703 number_cache.n_global_dofs = 0;
726 template <
int dim,
int spacedim>
731 if (cell == this->get_triangulation().end(level))
738 template <
int dim,
int spacedim>
746 while (i->has_children())
754 template <
int dim,
int spacedim>
765 template <
int dim,
int spacedim>
776 template <
int dim,
int spacedim>
788 template <
int dim,
int spacedim>
789 typename DoFHandler<dim, spacedim>::level_cell_iterator
795 if (cell == this->get_triangulation().end(level))
796 return end_mg(level);
797 return level_cell_iterator (*cell,
this);
801 template <
int dim,
int spacedim>
802 typename DoFHandler<dim, spacedim>::level_cell_iterator
810 return level_cell_iterator (*cell,
this);
814 template <
int dim,
int spacedim>
815 typename DoFHandler<dim, spacedim>::level_cell_iterator
818 return level_cell_iterator (&this->get_triangulation(), -1, -1,
this);
823 template <
int dim,
int spacedim>
833 template <
int dim,
int spacedim>
839 (begin_active(), end());
844 template <
int dim,
int spacedim>
850 (begin_mg(), end_mg());
856 template <
int dim,
int spacedim>
862 (begin(level), end(level));
867 template <
int dim,
int spacedim>
873 (begin_active(level), end_active(level));
878 template <
int dim,
int spacedim>
884 (begin_mg(level), end_mg(level));
896 return 2*get_fe().dofs_per_vertex;
902 template <
typename number>
909 i!=boundary_ids.end(); ++i)
910 Assert ((i->first == 0) || (i->first == 1),
913 return boundary_ids.size()*get_fe().dofs_per_vertex;
924 for (std::set<types::boundary_id>::const_iterator i=boundary_ids.begin();
925 i!=boundary_ids.end(); ++i)
926 Assert ((*i == 0) || (*i == 1),
929 return boundary_ids.size()*get_fe().dofs_per_vertex;
936 return 2*get_fe().dofs_per_vertex;
942 template <
typename number>
949 i!=boundary_ids.end(); ++i)
950 Assert ((i->first == 0) || (i->first == 1),
953 return boundary_ids.size()*get_fe().dofs_per_vertex;
964 for (std::set<types::boundary_id>::const_iterator i=boundary_ids.begin();
965 i!=boundary_ids.end(); ++i)
966 Assert ((*i == 0) || (*i == 1),
969 return boundary_ids.size()*get_fe().dofs_per_vertex;
974 template<
int dim,
int spacedim>
977 std::set<int> boundary_dofs;
979 const unsigned int dofs_per_face = get_fe().dofs_per_face;
980 std::vector<types::global_dof_index> dofs_on_face(dofs_per_face);
998 for (; cell!=endc; ++cell)
999 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
1000 if (cell->at_boundary(f))
1002 cell->face(f)->get_dof_indices (dofs_on_face);
1003 for (
unsigned int i=0; i<dofs_per_face; ++i)
1004 boundary_dofs.insert(dofs_on_face[i]);
1007 return boundary_dofs.size();
1012 template <
int dim,
int spacedim>
1019 std::set<types::global_dof_index> boundary_dofs;
1021 const unsigned int dofs_per_face = get_fe().dofs_per_face;
1022 std::vector<types::global_dof_index> dofs_on_face(dofs_per_face);
1029 for (; cell!=endc; ++cell)
1030 for (
unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
1031 if (cell->at_boundary(f)
1033 (std::find (boundary_ids.begin(),
1035 cell->face(f)->boundary_id()) !=
1036 boundary_ids.end()))
1038 cell->face(f)->get_dof_indices (dofs_on_face);
1039 for (
unsigned int i=0; i<dofs_per_face; ++i)
1040 boundary_dofs.insert(dofs_on_face[i]);
1043 return boundary_dofs.size();
1048 template<
int dim,
int spacedim>
1058 sizeof (number_cache) +
1061 for (
unsigned int i=0; i<levels.size(); ++i)
1064 for (
unsigned int level = 0; level < mg_levels.size (); ++level)
1065 mem += mg_levels[level]->memory_consumption ();
1070 for (
unsigned int i = 0; i < mg_vertex_dofs.size (); ++i)
1071 mem +=
sizeof (MGVertexDoFs) + (1 + mg_vertex_dofs[i].get_finest_level () - mg_vertex_dofs[i].get_coarsest_level ()) *
sizeof (
types::global_dof_index);
1078 template<
int dim,
int spacedim>
1103 policy->distribute_dofs (*
this,number_cache);
1110 block_info_object.initialize(*
this,
false,
true);
1114 template<
int dim,
int spacedim>
1118 Assert(levels.size()>0,
ExcMessage(
"Distribute active DoFs using distribute_dofs() before calling distribute_mg_dofs()."));
1122 Assert(old_fe == &fe,
ExcMessage(
"You are required to use the same FE for level and active DoFs!") );
1126 ExcMessage(
"The mesh smoothing requirement 'limit_level_difference_at_vertices' has to be set for using multigrid!"));
1130 internal::DoFHandler::Implementation::reserve_space_mg (*
this);
1133 mg_number_cache.resize((*tria).n_levels());
1137 policy->distribute_mg_dofs (*
this, mg_number_cache);
1144 block_info_object.initialize (*
this,
true,
false);
1147 template<
int dim,
int spacedim>
1153 template<
int dim,
int spacedim>
1156 for (
unsigned int i = 0; i < mg_levels.size (); ++i)
1157 delete mg_levels[i];
1163 std::vector<MGVertexDoFs> tmp;
1165 std::swap (mg_vertex_dofs, tmp);
1167 mg_number_cache.clear();
1171 template<
int dim,
int spacedim>
1174 block_info_object.initialize_local(*
this);
1179 template<
int dim,
int spacedim>
1192 template <
int dim,
int spacedim>
1196 Assert(levels.size()>0,
ExcMessage(
"You need to distribute DoFs before you can renumber them."));
1198 Assert (new_numbers.size() == n_locally_owned_dofs(),
1211 if (n_locally_owned_dofs() == n_dofs())
1213 std::vector<types::global_dof_index> tmp(new_numbers);
1214 std::sort (tmp.begin(), tmp.end());
1215 std::vector<types::global_dof_index>::const_iterator p = tmp.begin();
1217 for (; p!=tmp.end(); ++p, ++i)
1218 Assert (*p == i, ExcNewNumbersNotConsecutive(i));
1222 Assert (new_numbers[i] < n_dofs(),
1223 ExcMessage (
"New DoF index is not less than the total number of dofs."));
1226 policy->renumber_dofs (new_numbers, *
this,number_cache);
1230 template <
int dim,
int spacedim>
1233 const std::vector<types::global_dof_index> &)
1241 const std::vector<types::global_dof_index> &new_numbers)
1243 Assert(mg_levels.size()>0 && levels.size()>0,
1244 ExcMessage(
"You need to distribute active and level DoFs before you can renumber level DoFs."));
1253 for (std::vector<MGVertexDoFs>::iterator i=mg_vertex_dofs.begin();
1254 i!=mg_vertex_dofs.end(); ++i)
1257 if ((i->get_coarsest_level() <= level) &&
1258 (i->get_finest_level() >= level))
1259 for (
unsigned int d=0; d<this->get_fe().dofs_per_vertex; ++d)
1260 i->set_index (level, d, new_numbers[i->get_index (level, d)]);
1262 for (std::vector<types::global_dof_index>::iterator i=mg_levels[level]->dof_object.dofs.begin();
1263 i!=mg_levels[level]->dof_object.dofs.end(); ++i)
1268 *i = new_numbers[*i];
1277 const std::vector<types::global_dof_index> &new_numbers)
1279 Assert(mg_levels.size()>0 && levels.size()>0,
1280 ExcMessage(
"You need to distribute active and level DoFs before you can renumber level DoFs."));
1281 Assert (new_numbers.size() == n_dofs(level),
1284 for (std::vector<MGVertexDoFs>::iterator i=mg_vertex_dofs.begin();
1285 i!=mg_vertex_dofs.end(); ++i)
1288 if ((i->get_coarsest_level() <= level) &&
1289 (i->get_finest_level() >= level))
1290 for (
unsigned int d=0;
d<this->get_fe().dofs_per_vertex; ++
d)
1291 i->set_index (level, d, new_numbers[i->get_index (level, d)]);
1293 if (this->get_fe().dofs_per_line > 0)
1296 std::vector<bool> user_flags;
1297 this->get_triangulation().save_user_flags(user_flags);
1298 const_cast<Triangulation<2> &
>(this->get_triangulation()).clear_user_flags ();
1304 level_cell_iterator cell, endc = end(level);
1305 for (cell = begin(level); cell != endc; ++cell)
1306 for (
unsigned int line=0; line < GeometryInfo<2>::faces_per_cell; ++line)
1307 cell->face(line)->set_user_flag();
1309 for (cell_iterator cell = begin(); cell != end(); ++cell)
1310 for (
unsigned int l=0; l<GeometryInfo<2>::lines_per_cell; ++
l)
1311 if (cell->line(l)->user_flag_set())
1313 for (
unsigned int d=0;
d<this->get_fe().dofs_per_line; ++
d)
1314 cell->line(l)->set_mg_dof_index (level, d,
1315 new_numbers[cell->line(l)->mg_dof_index(level, d)]);
1316 cell->line(l)->clear_user_flag();
1319 const_cast<Triangulation<2> &
>(this->get_triangulation()).load_user_flags (user_flags);
1322 for (std::vector<types::global_dof_index>::iterator i=mg_levels[level]->dof_object.dofs.begin();
1323 i!=mg_levels[level]->dof_object.dofs.end(); ++i)
1328 *i = new_numbers[*i];
1337 const std::vector<types::global_dof_index> &new_numbers)
1339 Assert(mg_levels.size()>0 && levels.size()>0,
1340 ExcMessage(
"You need to distribute active and level DoFs before you can renumber level DoFs."));
1341 Assert (new_numbers.size() == n_dofs(level),
1344 for (std::vector<MGVertexDoFs>::iterator i=mg_vertex_dofs.begin();
1345 i!=mg_vertex_dofs.end(); ++i)
1348 if ((i->get_coarsest_level() <= level) &&
1349 (i->get_finest_level() >= level))
1350 for (
unsigned int d=0;
d<this->get_fe().dofs_per_vertex; ++
d)
1351 i->set_index (level, d, new_numbers[i->get_index (level, d)]);
1354 if (this->get_fe().dofs_per_line > 0)
1357 std::vector<bool> user_flags;
1358 this->get_triangulation().save_user_flags(user_flags);
1359 const_cast<Triangulation<3> &
>(this->get_triangulation()).clear_user_flags ();
1365 level_cell_iterator cell, endc = end(level);
1366 for (cell = begin(level) ; cell != endc ; ++cell)
1367 for (
unsigned int line=0; line < GeometryInfo<3>::lines_per_cell; ++line)
1368 cell->line(line)->set_user_flag();
1371 for (cell = begin(level); cell != endc; ++cell)
1372 for (
unsigned int l=0; l<GeometryInfo<3>::lines_per_cell; ++
l)
1373 if (cell->line(l)->user_flag_set())
1375 for (
unsigned int d=0;
d<this->get_fe().dofs_per_line; ++
d)
1376 cell->line(l)->set_mg_dof_index (level, d,
1377 new_numbers[cell->line(l)->mg_dof_index(level, d)]);
1378 cell->line(l)->clear_user_flag();
1381 const_cast<Triangulation<3> &
>(this->get_triangulation()).load_user_flags (user_flags);
1385 if (this->get_fe().dofs_per_quad > 0)
1388 std::vector<bool> user_flags;
1389 this->get_triangulation().save_user_flags(user_flags);
1390 const_cast<Triangulation<3> &
>(this->get_triangulation()).clear_user_flags ();
1396 level_cell_iterator cell, endc = end(level);
1397 for (cell = begin(level) ; cell != endc; ++cell)
1398 for (
unsigned int quad=0; quad < GeometryInfo<3>::faces_per_cell; ++quad)
1399 cell->face(quad)->set_user_flag();
1401 for (cell = begin(level); cell != endc; ++cell)
1402 for (
unsigned int q=0; q<GeometryInfo<3>::quads_per_cell; ++q)
1403 if (cell->quad(q)->user_flag_set())
1405 for (
unsigned int d=0;
d<this->get_fe().dofs_per_quad; ++
d)
1406 cell->quad(q)->set_mg_dof_index (level, d,
1407 new_numbers[cell->quad(q)->mg_dof_index(level, d)]);
1408 cell->quad(q)->clear_user_flag();
1411 const_cast<Triangulation<3> &
>(this->get_triangulation()).load_user_flags (user_flags);
1415 for (std::vector<types::global_dof_index>::iterator i=mg_levels[level]->dof_object.dofs.begin();
1416 i!=mg_levels[level]->dof_object.dofs.end(); ++i)
1421 *i = new_numbers[*i];
1429 template <
int dim,
int spacedim>
1438 template <
int dim,
int spacedim>
1445 return get_fe().dofs_per_vertex;
1447 return (3*get_fe().dofs_per_vertex +
1448 2*get_fe().dofs_per_line);
1466 return (19*get_fe().dofs_per_vertex +
1467 28*get_fe().dofs_per_line +
1468 8*get_fe().dofs_per_quad);
1477 template<
int dim,
int spacedim>
1480 for (
unsigned int i=0; i<levels.size(); ++i)
1487 std::vector<types::global_dof_index> tmp;
1488 std::swap (vertex_dofs, tmp);
1490 number_cache.clear ();
1493 template<
int dim,
int spacedim>
1494 template<
int structdim>
1497 const unsigned int obj_level,
1498 const unsigned int obj_index,
1499 const unsigned int fe_index,
1500 const unsigned int local_index)
const 1502 return internal::DoFHandler::Implementation::get_dof_index (*
this, *this->mg_levels[obj_level],
1503 *this->mg_faces, obj_index,
1504 fe_index, local_index,
1508 template<
int dim,
int spacedim>
1509 template<
int structdim>
1512 internal::DoFHandler::Implementation::set_dof_index (*
this, *this->mg_levels[obj_level], *this->mg_faces, obj_index, fe_index, local_index, global_index,
internal::int2type<structdim> ());
1516 template<
int dim,
int spacedim>
1522 template<
int dim,
int spacedim>
1526 delete[] indices_offset;
1529 template<
int dim,
int spacedim>
1538 if (indices_offset != 0)
1540 delete[] indices_offset;
1544 coarsest_level = cl;
1550 const unsigned int n_levels = finest_level - coarsest_level + 1;
1551 const unsigned int n_indices = n_levels * dofs_per_vertex;
1554 Assert (indices != 0, ExcNoMemory ());
1556 for (
unsigned int i = 0; i < n_indices; ++i)
1560 Assert (indices != 0, ExcNoMemory ());
1562 for (
unsigned int i = 0; i < n_levels; ++i)
1563 indices_offset[i] = i * dofs_per_vertex;
1566 template<
int dim,
int spacedim>
1569 return coarsest_level;
1572 template<
int dim,
int spacedim>
1575 return finest_level;
1580 #include "dof_handler.inst" 1583 DEAL_II_NAMESPACE_CLOSE
IteratorRange< active_cell_iterator > active_cell_iterators_on_level(const unsigned int level) const
unsigned int max_couplings_between_dofs() const
unsigned int get_coarsest_level() const
std::vector< MGVertexDoFs > mg_vertex_dofs
level_cell_iterator end_mg() const
static const unsigned int invalid_unsigned_int
cell_iterator begin(const unsigned int level=0) const
std_cxx11::shared_ptr<::internal::DoFHandler::Policy::PolicyBase< dim, spacedim > > policy
static ::ExceptionBase & ExcRenumberingIncomplete()
cell_iterator end() const
SmartPointer< const FiniteElement< dim, spacedim >, DoFHandler< dim, spacedim > > selected_fe
IteratorRange< level_cell_iterator > mg_cell_iterators_on_level(const unsigned int level) const
std::vector<::internal::DoFHandler::DoFLevel< dim > * > levels
ActiveSelector::active_cell_iterator active_cell_iterator
SmartPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > tria
#define AssertThrow(cond, exc)
virtual std::size_t memory_consumption() const
void initialize(const Triangulation< dim, spacedim > &tria, const FiniteElement< dim, spacedim > &fe)
const FiniteElement< dim, spacedim > & get_fe() const
virtual void distribute_mg_dofs(const FiniteElement< dim, spacedim > &fe)
active_cell_iterator begin_active(const unsigned int level=0) const
static ::ExceptionBase & ExcMessage(std::string arg1)
ActiveSelector::cell_iterator cell_iterator
unsigned int global_dof_index
#define Assert(cond, exc)
IteratorRange< active_cell_iterator > active_cell_iterators() const
void set_dof_index(const ::DoFHandler< dh_dim, spacedim > &dof_handler, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index)
active_cell_iterator end_active(const unsigned int level) const
types::global_dof_index n_boundary_dofs() const
types::global_dof_index n_dofs() const
level_cell_iterator begin_mg(const unsigned int level=0) const
static ::ExceptionBase & ExcRenumberingIncomplete()
::internal::DoFHandler::DoFFaces< dim > * faces
internal::DoFHandler::DoFObjects< 1 > lines
DoFObjects< dim > dof_object
unsigned int get_finest_level() const
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)
virtual unsigned int n_global_levels() const
void renumber_dofs(const std::vector< types::global_dof_index > &new_numbers)
void initialize_local_block_info()
static unsigned int max_couplings_between_dofs(const DoFHandler< 1, spacedim > &dof_handler)
internal::DoFHandler::DoFObjects< 1 > lines
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
IteratorRange< level_cell_iterator > mg_cell_iterators() const
internal::DoFHandler::DoFObjects< 2 > quads
void init(const unsigned int coarsest_level, const unsigned int finest_level, const unsigned int dofs_per_vertex)
virtual void distribute_dofs(const FiniteElement< dim, spacedim > &fe)
types::global_dof_index get_dof_index(const ::DoFHandler< dh_dim, spacedim > &dof_handler, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index) const
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
static ::ExceptionBase & ExcNotImplemented()
const Triangulation< dim, spacedim > & get_triangulation() const
Iterator points to a valid object.
unsigned char boundary_id
const types::boundary_id internal_face_boundary_id
unsigned int max_couplings_between_boundary_dofs() const
IteratorState::IteratorStates state() const
static void reserve_space(DoFHandler< 1, spacedim > &dof_handler)
IteratorRange< cell_iterator > cell_iterators() const
std::vector< types::global_dof_index > vertex_dofs
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static ::ExceptionBase & ExcInternalError()