From: Wolfgang Bangerth Date: Tue, 10 Sep 2013 18:13:31 +0000 (+0000) Subject: Implement the hp::DoFHandler cell dof indices cache in the same way as we already... X-Git-Tag: v8.1.0~836 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=21c1d982fb58ce6c78dd1c598cb949ab5489b5e0;p=dealii.git Implement the hp::DoFHandler cell dof indices cache in the same way as we already have for the ::DoFHandler. git-svn-id: https://svn.dealii.org/trunk@30666 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/include/deal.II/dofs/dof_accessor.templates.h b/deal.II/include/deal.II/dofs/dof_accessor.templates.h index a2a4c03f10..5fee314d9e 100644 --- a/deal.II/include/deal.II/dofs/dof_accessor.templates.h +++ b/deal.II/include/deal.II/dofs/dof_accessor.templates.h @@ -2556,15 +2556,44 @@ namespace internal template static void - update_cell_dof_indices_cache (const DoFCellAccessor, lda> &) + update_cell_dof_indices_cache (const DoFCellAccessor, lda> &accessor) { -//TODO[WB]: should implement a dof indices cache for hp as well + // caches are only for cells with DoFs, i.e., for active ones + if (accessor.has_children()) + return; - // not implemented, but should also - // not be called - Assert (false, ExcNotImplemented()); + const unsigned int dofs_per_cell = accessor.get_fe().dofs_per_cell; + + // make sure the cache is at least + // as big as we need it when + // writing to the last element of + // this cell + Assert (static_cast(accessor.present_index) + < + accessor.dof_handler->levels[accessor.present_level] + ->cell_cache_offsets.size(), + ExcInternalError()); + Assert (accessor.dof_handler->levels[accessor.present_level] + ->cell_cache_offsets[accessor.present_index] + <= + accessor.dof_handler->levels[accessor.present_level] + ->cell_dof_indices_cache.size(), + ExcInternalError()); + + std::vector dof_indices (dofs_per_cell); + accessor.dealii::DoFAccessor,lda>::get_dof_indices + (dof_indices, accessor.active_fe_index()); + + types::global_dof_index *next_dof_index + = &accessor.dof_handler->levels[accessor.present_level] + ->cell_dof_indices_cache[accessor.dof_handler->levels[accessor.present_level] + ->cell_cache_offsets[accessor.present_index]]; + for (unsigned int i=0; i, lda> &accessor, std::vector &dof_indices) { - // no caching for - // hp::DoFHandler - // implemented - typedef - dealii::DoFAccessor, lda> - DoFAccessor; - accessor.DoFAccessor::get_dof_indices (dof_indices, - accessor.active_fe_index()); + AssertDimension (dof_indices.size(), accessor.get_fe().dofs_per_cell); + + Assert (!accessor.has_children(), + ExcMessage ("Cell must be active.")); + + types::global_dof_index *cache = &accessor.dof_handler->levels[accessor.level()] + ->cell_dof_indices_cache[accessor.dof_handler->levels[accessor.level()]->cell_cache_offsets[accessor.present_index]]; + for (unsigned int i=0; i dof_indices; + /** + * The offsets for each cell of the cache that holds all DoF indices. + */ + std::vector cell_cache_offsets; + + /** + * Cache for the DoF indices + * on cells. The size of this + * array equals the sum over all cells of selected_fe[active_fe_index[cell]].dofs_per_cell. + */ + std::vector cell_dof_indices_cache; + public: /** @@ -275,6 +291,7 @@ namespace internal */ template friend class dealii::hp::DoFHandler; friend struct dealii::internal::hp::DoFHandler::Implementation; + friend struct dealii::internal::DoFCellAccessor::Implementation; }; diff --git a/deal.II/source/hp/dof_handler.cc b/deal.II/source/hp/dof_handler.cc index eeec5b43e1..d6afadf43c 100644 --- a/deal.II/source/hp/dof_handler.cc +++ b/deal.II/source/hp/dof_handler.cc @@ -580,8 +580,13 @@ namespace internal = std::vector ( dof_handler.tria->n_raw_lines(level), (DoFLevel::offset_type)(-1)); + dof_handler.levels[level]->cell_cache_offsets + = std::vector ( + dof_handler.tria->n_raw_lines(level), + (DoFLevel::offset_type)(-1)); types::global_dof_index next_free_dof = 0; + types::global_dof_index cache_size = 0; for (typename DoFHandler::active_cell_iterator cell=dof_handler.begin_active(level); cell!=dof_handler.end_active(level); ++cell) @@ -589,11 +594,17 @@ namespace internal { dof_handler.levels[level]->dof_offsets[cell->index()] = next_free_dof; next_free_dof += cell->get_fe().dofs_per_line; + + dof_handler.levels[level]->cell_cache_offsets[cell->index()] = cache_size; + cache_size += cell->get_fe().dofs_per_cell; } dof_handler.levels[level]->dof_indices = std::vector (next_free_dof, DoFHandler::invalid_dof_index); + dof_handler.levels[level]->cell_dof_indices_cache + = std::vector (cache_size, + DoFHandler::invalid_dof_index); } // safety check: make sure that @@ -697,8 +708,13 @@ namespace internal = std::vector ( dof_handler.tria->n_raw_quads(level), (DoFLevel::offset_type)(-1)); + dof_handler.levels[level]->cell_cache_offsets + = std::vector ( + dof_handler.tria->n_raw_quads(level), + (DoFLevel::offset_type)(-1)); types::global_dof_index next_free_dof = 0; + types::global_dof_index cache_size = 0; for (typename DoFHandler::active_cell_iterator cell=dof_handler.begin_active(level); cell!=dof_handler.end_active(level); ++cell) @@ -706,11 +722,17 @@ namespace internal { dof_handler.levels[level]->dof_offsets[cell->index()] = next_free_dof; next_free_dof += cell->get_fe().dofs_per_quad; + + dof_handler.levels[level]->cell_cache_offsets[cell->index()] = cache_size; + cache_size += cell->get_fe().dofs_per_cell; } dof_handler.levels[level]->dof_indices = std::vector (next_free_dof, DoFHandler::invalid_dof_index); + dof_handler.levels[level]->cell_dof_indices_cache + = std::vector (cache_size, + DoFHandler::invalid_dof_index); } // safety check: make sure that @@ -1059,11 +1081,16 @@ namespace internal for (unsigned int level=0; leveln_levels(); ++level) { dof_handler.levels[level]->dof_offsets - = std::vector - (dof_handler.tria->n_raw_hexs(level), - (DoFLevel::offset_type)(-1)); + = std::vector ( + dof_handler.tria->n_raw_hexs(level), + (DoFLevel::offset_type)(-1)); + dof_handler.levels[level]->cell_cache_offsets + = std::vector ( + dof_handler.tria->n_raw_hexs(level), + (DoFLevel::offset_type)(-1)); types::global_dof_index next_free_dof = 0; + types::global_dof_index cache_size = 0; for (typename DoFHandler::active_cell_iterator cell=dof_handler.begin_active(level); cell!=dof_handler.end_active(level); ++cell) @@ -1071,11 +1098,17 @@ namespace internal { dof_handler.levels[level]->dof_offsets[cell->index()] = next_free_dof; next_free_dof += cell->get_fe().dofs_per_hex; + + dof_handler.levels[level]->cell_cache_offsets[cell->index()] = cache_size; + cache_size += cell->get_fe().dofs_per_cell; } dof_handler.levels[level]->dof_indices = std::vector (next_free_dof, DoFHandler::invalid_dof_index); + dof_handler.levels[level]->cell_dof_indices_cache + = std::vector (cache_size, + DoFHandler::invalid_dof_index); } // safety check: make sure that @@ -2689,6 +2722,11 @@ namespace hp = std::vector (1, number_cache.locally_owned_dofs); + // update the cache used for cell dof indices and compress the data on the levels + for (active_cell_iterator cell = begin_active(); + cell != end(); ++cell) + cell->update_cell_dof_indices_cache (); + for (unsigned int level=0; levelcompress_data (*finite_elements); @@ -2736,6 +2774,11 @@ namespace hp // do the renumbering renumber_dofs_internal (new_numbers, dealii::internal::int2type()); + // update the cache used for cell dof indices + for (active_cell_iterator cell = begin_active(); + cell != end(); ++cell) + cell->update_cell_dof_indices_cache (); + // now re-compress the dof indices for (unsigned int level=0; levelcompress_data (*finite_elements);