From: kronbichler Date: Mon, 9 Mar 2009 11:15:09 +0000 (+0000) Subject: Improved some details in CellSimilarity concept. X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=94b32cc8b4bf1037f2c4efe6d36acdbf8e3badec;p=dealii-svn.git Improved some details in CellSimilarity concept. git-svn-id: https://svn.dealii.org/trunk@18467 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/deal.II/include/fe/fe_update_flags.h b/deal.II/deal.II/include/fe/fe_update_flags.h index 9eb404b7a7..a7e8bc797a 100644 --- a/deal.II/deal.II/include/fe/fe_update_flags.h +++ b/deal.II/deal.II/include/fe/fe_update_flags.h @@ -340,8 +340,29 @@ operator &= (UpdateFlags &f1, UpdateFlags f2) +/** + * This enum definition is used for storing similarities of the current cell + * to the previously visited cell. This information is used for reusing data + * when calling the method FEValues::reinit() (like derivatives, which do + * not change if one cell is just a translation of the previous). Currently, + * this variable does only recognize a translation. However, this concept + * makes it easy to add additional staties to be detected in + * FEValues/FEFaceValues for making use of these similarities as well. + */ +namespace CellSimilarity +{ + enum Similarity + { + none, + translation + }; +} + + /*@}*/ + + DEAL_II_NAMESPACE_CLOSE #endif diff --git a/deal.II/deal.II/include/fe/fe_values.h b/deal.II/deal.II/include/fe/fe_values.h index 4d77e95177..73ee333fcf 100644 --- a/deal.II/deal.II/include/fe/fe_values.h +++ b/deal.II/deal.II/include/fe/fe_values.h @@ -2244,7 +2244,6 @@ class FEValuesBase : protected FEValuesData, const std::vector > & get_cell_normal_vectors () const; - /** * Return the outward normal vector to * the cell at the ith quadrature @@ -2253,6 +2252,17 @@ class FEValuesBase : protected FEValuesData, */ const Point & cell_normal_vector (const unsigned int i) const; + /** + * Return the relation of the current + * cell to the previous cell. This + * allows re-use of some cell data + * (like local matrices for equations + * with constant coefficients) if the + * result is + * CellSimilarity::translation. + */ + const enum CellSimilarity::Similarity get_cell_similarity () const; + /** * Determine an estimate for the * memory consumption (in bytes) @@ -2603,6 +2613,28 @@ class FEValuesBase : protected FEValuesData, */ UpdateFlags compute_update_flags (const UpdateFlags update_flags) const; + /** + * An enum variable that can store + * different states of the current cell + * in comparison to the previously + * visited cell. If wanted, additional + * states can be checked here and used + * in one of the methods used during + * reinit. + */ + enum CellSimilarity::Similarity cell_similarity; + + /** + * A function that checks whether the + * new cell is similar to the one + * previously used. Then, a significant + * amount of the data can be reused, + * e.g. the derivatives of the basis + * functions in real space, shape_grad. + */ + void + check_cell_similarity (const typename Triangulation::cell_iterator &cell); + private: /** * Copy constructor. Since @@ -2679,16 +2711,16 @@ class FEValues : public FEValuesBase */ FEValues (const Mapping &mapping, const FiniteElement &fe, - const Quadrature &quadrature, - const UpdateFlags update_flags); + const Quadrature &quadrature, + const UpdateFlags update_flags); /** * Constructor. Uses MappingQ1 * implicitly. */ FEValues (const FiniteElement &fe, - const Quadrature &quadrature, - const UpdateFlags update_flags); + const Quadrature &quadrature, + const UpdateFlags update_flags); /** * Reinitialize the gradients, @@ -2815,28 +2847,6 @@ class FEValues : public FEValuesBase */ void initialize (const UpdateFlags update_flags); - /** - * An enum variable that can store - * different states of the current cell - * in comparison to the previously - * visited cell. If wanted, additional - * states can be checked here and used - * in one of the methods used during - * reinit. - */ - enum CellSimilarity::Similarity cell_similarity; - - /** - * A function that checks whether the - * new cell is similar to the one - * previously used. Then, a significant - * amount of the data can be reused, - * e.g. the derivatives of the basis - * functions in real space, shape_grad. - */ - void - check_cell_similarity (const typename Triangulation::cell_iterator &cell); - /** * The reinit() functions do * only that part of the work @@ -2895,12 +2905,12 @@ class FEFaceValuesBase : public FEValuesBase * of faces times the number of subfaces * per face. */ - FEFaceValuesBase (const unsigned int n_q_points, - const unsigned int dofs_per_cell, - const UpdateFlags update_flags, + FEFaceValuesBase (const unsigned int n_q_points, + const unsigned int dofs_per_cell, + const UpdateFlags update_flags, const Mapping &mapping, const FiniteElement &fe, - const Quadrature& quadrature); + const Quadrature& quadrature); /** * Return the outward normal vector to @@ -3026,16 +3036,16 @@ class FEFaceValues : public FEFaceValuesBase */ FEFaceValues (const Mapping &mapping, const FiniteElement &fe, - const Quadrature &quadrature, - const UpdateFlags update_flags); + const Quadrature &quadrature, + const UpdateFlags update_flags); /** * Constructor. Uses MappingQ1 * implicitly. */ FEFaceValues (const FiniteElement &fe, - const Quadrature &quadrature, - const UpdateFlags update_flags); + const Quadrature &quadrature, + const UpdateFlags update_flags); /** * Reinitialize the gradients, Jacobi @@ -3044,7 +3054,7 @@ class FEFaceValues : public FEFaceValuesBase * and the given finite element. */ void reinit (const typename DoFHandler::cell_iterator &cell, - const unsigned int face_no); + const unsigned int face_no); /** * Reinitialize the gradients, @@ -3060,7 +3070,7 @@ class FEFaceValues : public FEFaceValuesBase * object. */ void reinit (const typename hp::DoFHandler::cell_iterator &cell, - const unsigned int face_no); + const unsigned int face_no); /** * Reinitialize the gradients, @@ -3076,7 +3086,7 @@ class FEFaceValues : public FEFaceValuesBase * object. */ void reinit (const typename MGDoFHandler::cell_iterator &cell, - const unsigned int face_no); + const unsigned int face_no); /** * Reinitialize the gradients, @@ -3105,7 +3115,7 @@ class FEFaceValues : public FEFaceValuesBase * handler type objects. */ void reinit (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no); + const unsigned int face_no); /** * Return a reference to this diff --git a/deal.II/deal.II/include/fe/mapping.h b/deal.II/deal.II/include/fe/mapping.h index cfbaebbc29..290a9b0757 100644 --- a/deal.II/deal.II/include/fe/mapping.h +++ b/deal.II/deal.II/include/fe/mapping.h @@ -32,24 +32,6 @@ template class FEValues; template class FEFaceValues; template class FESubfaceValues; - /** - * An enum variable that can store - * different states of the current cell - * in comparison to the previously - * visited cell. If wanted, additional - * states can be checked here and used - * in one of the methods used during - * reinit. - */ -namespace CellSimilarity -{ - enum Similarity - { - no_similarity, - translation - }; -} - /** * The transformation type used * for the Mapping::transform() functions. diff --git a/deal.II/deal.II/source/fe/fe_system.cc b/deal.II/deal.II/source/fe/fe_system.cc index 9989fb6ea8..53917b353d 100644 --- a/deal.II/deal.II/source/fe/fe_system.cc +++ b/deal.II/deal.II/source/fe/fe_system.cc @@ -690,7 +690,7 @@ fill_fe_face_values (const Mapping &mapping, FEValuesData &data) const { compute_fill (mapping, cell, face_no, invalid_face_number, quadrature, - CellSimilarity::no_similarity, mapping_data, fe_data, data); + CellSimilarity::none, mapping_data, fe_data, data); } @@ -709,7 +709,7 @@ fill_fe_subface_values (const Mapping &mappin FEValuesData &data) const { compute_fill (mapping, cell, face_no, sub_no, quadrature, - CellSimilarity::no_similarity, mapping_data, fe_data, data); + CellSimilarity::none, mapping_data, fe_data, data); } diff --git a/deal.II/deal.II/source/fe/fe_values.cc b/deal.II/deal.II/source/fe/fe_values.cc index 0cb2bd6f11..c66288f97f 100644 --- a/deal.II/deal.II/source/fe/fe_values.cc +++ b/deal.II/deal.II/source/fe/fe_values.cc @@ -3239,6 +3239,70 @@ FEValuesBase::compute_update_flags (const UpdateFlags update_flags } + +template +inline +void +FEValuesBase::check_cell_similarity + (const typename Triangulation::cell_iterator &cell) +{ + // case that there has not been any cell + // before + if (&*this->present_cell == 0) + { + cell_similarity = CellSimilarity::none; + return; + } + + const typename Triangulation::cell_iterator & present_cell = + *this->present_cell; + + // test for translation + { + // otherwise, go through the vertices and + // check... The cell is a translation of + // the previous one in case the distance + // between the individual vertices in the + // two cell is the same for all the + // vertices. So do the check by first + // getting the distance on the first + // vertex, and then checking whether all + // others have the same. + bool is_translation = true; + const Point dist = cell->vertex(0) - present_cell->vertex(0); + for (unsigned int i=1; i::vertices_per_cell; ++i) + { + Point dist_new = cell->vertex(i) - present_cell->vertex(i) - dist; + if (dist_new.norm_square() > 1e-28) + { + is_translation = false; + break; + } + } + + cell_similarity = (is_translation == true + ? + CellSimilarity::translation + : + CellSimilarity::none); + return; + } + + // TODO: here, one could implement other + // checks for similarity, e.g. for + // children of a parallelogram. +} + + + +template +const enum CellSimilarity::Similarity +FEValuesBase::get_cell_similarity () const +{ + return cell_similarity; +} + + /*------------------------------- FEValues -------------------------------*/ @@ -3469,53 +3533,6 @@ void FEValues::reinit (const typename Triangulation: -template -inline -void -FEValues::check_cell_similarity (const typename Triangulation::cell_iterator &cell) -{ - // case that there has not been any cell - // before - if (&*this->present_cell == 0) - { - cell_similarity = CellSimilarity::no_similarity; - return; - } - - const typename Triangulation::cell_iterator & present_cell = - *this->present_cell; - - // test for translation - { - // otherwise, go through the vertices and - // check... - bool is_translation = true; - const Point dist = cell->vertex(0) - present_cell->vertex(0); - for (unsigned int i=1; i::vertices_per_cell; ++i) - { - Point dist_new = cell->vertex(i) - present_cell->vertex(i) - dist; - if (dist_new.norm_square() > 1e-28) - { - is_translation = false; - break; - } - } - - cell_similarity = (is_translation == true - ? - CellSimilarity::translation - : - CellSimilarity::no_similarity); - return; - } - - // TODO: here, one could implement other - // checks for similarity, e.g. for - // children of a parallelogram. -} - - - template void FEValues::do_reinit () { diff --git a/deal.II/deal.II/source/fe/mapping_cartesian.cc b/deal.II/deal.II/source/fe/mapping_cartesian.cc index ec021ded6c..db2c77f370 100644 --- a/deal.II/deal.II/source/fe/mapping_cartesian.cc +++ b/deal.II/deal.II/source/fe/mapping_cartesian.cc @@ -431,7 +431,7 @@ MappingCartesian::fill_fe_face_values ( InternalData &data = static_cast (mapping_data); compute_fill (cell, face_no, invalid_face_number, - CellSimilarity::no_similarity, + CellSimilarity::none, data, quadrature_points, normal_vectors); @@ -484,7 +484,7 @@ MappingCartesian::fill_fe_subface_values ( Assert (dynamic_cast (&mapping_data) != 0, ExcInternalError()); InternalData &data = static_cast (mapping_data); - compute_fill (cell, face_no, sub_no, CellSimilarity::no_similarity, + compute_fill (cell, face_no, sub_no, CellSimilarity::none, data, quadrature_points, normal_vectors); diff --git a/deal.II/deal.II/source/fe/mapping_q.cc b/deal.II/deal.II/source/fe/mapping_q.cc index fc6c6336b5..106bd50532 100644 --- a/deal.II/deal.II/source/fe/mapping_q.cc +++ b/deal.II/deal.II/source/fe/mapping_q.cc @@ -351,7 +351,7 @@ MappingQ::fill_fe_values ( else { p_data=&data; - cell_similarity = CellSimilarity::no_similarity; + cell_similarity = CellSimilarity::none; } MappingQ1::fill_fe_values(cell, q, cell_similarity, *p_data, diff --git a/deal.II/deal.II/source/fe/mapping_q1.cc b/deal.II/deal.II/source/fe/mapping_q1.cc index d7b3aa9f00..c5985e9b5d 100644 --- a/deal.II/deal.II/source/fe/mapping_q1.cc +++ b/deal.II/deal.II/source/fe/mapping_q1.cc @@ -962,7 +962,7 @@ MappingQ1::compute_fill_face ( std::vector > &boundary_forms, std::vector > &normal_vectors) const { - compute_fill (cell, n_q_points, data_set, CellSimilarity::no_similarity, + compute_fill (cell, n_q_points, data_set, CellSimilarity::none, data, quadrature_points); const UpdateFlags update_flags(data.current_update_flags()); diff --git a/deal.II/deal.II/source/fe/mapping_q1_eulerian.cc b/deal.II/deal.II/source/fe/mapping_q1_eulerian.cc index 287a809e10..3e1bd7f832 100644 --- a/deal.II/deal.II/source/fe/mapping_q1_eulerian.cc +++ b/deal.II/deal.II/source/fe/mapping_q1_eulerian.cc @@ -132,7 +132,7 @@ MappingQ1Eulerian::fill_fe_values ( // disable any previously detected // similarity and hand on to the // respective function of the base class. - cell_similarity = CellSimilarity::no_similarity; + cell_similarity = CellSimilarity::none; MappingQ1::fill_fe_values (cell, q, cell_similarity, mapping_data, quadrature_points, JxW_values, jacobians, jacobian_grads, inverse_jacobians, diff --git a/deal.II/deal.II/source/fe/mapping_q_eulerian.cc b/deal.II/deal.II/source/fe/mapping_q_eulerian.cc index cfa9ce1930..cff88e4ba0 100644 --- a/deal.II/deal.II/source/fe/mapping_q_eulerian.cc +++ b/deal.II/deal.II/source/fe/mapping_q_eulerian.cc @@ -199,7 +199,7 @@ MappingQEulerian::fill_fe_values ( // disable any previously detected // similarity and hand on to the respective // function of the base class. - cell_similarity = CellSimilarity::no_similarity; + cell_similarity = CellSimilarity::none; MappingQ::fill_fe_values (cell, q, cell_similarity, mapping_data, quadrature_points, JxW_values, jacobians, jacobian_grads, inverse_jacobians, diff --git a/deal.II/lac/source/trilinos_sparse_matrix.cc b/deal.II/lac/source/trilinos_sparse_matrix.cc index 9ef30e4069..496c11b175 100755 --- a/deal.II/lac/source/trilinos_sparse_matrix.cc +++ b/deal.II/lac/source/trilinos_sparse_matrix.cc @@ -406,7 +406,9 @@ namespace TrilinosWrappers for (int col=0; col < row_length; ++col) row_indices[col] = sparsity_pattern.column_number (row, col); - graph->InsertGlobalIndices (row, row_length, &row_indices[0]); + graph->InsertGlobalIndices (static_cast(row), + row_length, &row_indices[0]); + } // Now, fill the graph (sort indices, make