From: kronbichler Date: Mon, 9 Mar 2009 17:44:59 +0000 (+0000) Subject: A few updates of the CellSimilarity concept. X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=392d06ac896847bbea75724afde6f0cc278354db;p=dealii-svn.git A few updates of the CellSimilarity concept. git-svn-id: https://svn.dealii.org/trunk@18469 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 a7e8bc797a..33fe2b8384 100644 --- a/deal.II/deal.II/include/fe/fe_update_flags.h +++ b/deal.II/deal.II/include/fe/fe_update_flags.h @@ -354,7 +354,8 @@ namespace CellSimilarity enum Similarity { none, - translation + translation, + invalid_next_cell }; } diff --git a/deal.II/deal.II/include/fe/mapping.h b/deal.II/deal.II/include/fe/mapping.h index 290a9b0757..334d2bff0e 100644 --- a/deal.II/deal.II/include/fe/mapping.h +++ b/deal.II/deal.II/include/fe/mapping.h @@ -546,7 +546,7 @@ class Mapping : public Subscriptor virtual void fill_fe_values (const typename Triangulation::cell_iterator &cell, const Quadrature &quadrature, - const enum CellSimilarity::Similarity cell_similarity, + enum CellSimilarity::Similarity &cell_similarity, InternalDataBase &internal, std::vector > &quadrature_points, std::vector &JxW_values, diff --git a/deal.II/deal.II/include/fe/mapping_cartesian.h b/deal.II/deal.II/include/fe/mapping_cartesian.h index 9b094b6ef8..2645ad8d3f 100644 --- a/deal.II/deal.II/include/fe/mapping_cartesian.h +++ b/deal.II/deal.II/include/fe/mapping_cartesian.h @@ -61,7 +61,7 @@ class MappingCartesian : public Mapping virtual void fill_fe_values (const typename Triangulation::cell_iterator &cell, const Quadrature &quadrature, - const enum CellSimilarity::Similarity cell_similarity, + enum CellSimilarity::Similarity &cell_similarity, typename Mapping::InternalDataBase &mapping_data, std::vector > &quadrature_points, std::vector &JxW_values, diff --git a/deal.II/deal.II/include/fe/mapping_q.h b/deal.II/deal.II/include/fe/mapping_q.h index 038144041b..9ca3716412 100644 --- a/deal.II/deal.II/include/fe/mapping_q.h +++ b/deal.II/deal.II/include/fe/mapping_q.h @@ -205,7 +205,7 @@ class MappingQ : public MappingQ1 virtual void fill_fe_values (const typename Triangulation::cell_iterator &cell, const Quadrature &quadrature, - enum CellSimilarity::Similarity cell_similarity, + enum CellSimilarity::Similarity &cell_similarity, typename Mapping::InternalDataBase &mapping_data, typename std::vector > &quadrature_points, std::vector &JxW_values, diff --git a/deal.II/deal.II/include/fe/mapping_q1.h b/deal.II/deal.II/include/fe/mapping_q1.h index f21b8eee9a..8a8ee2f0a6 100644 --- a/deal.II/deal.II/include/fe/mapping_q1.h +++ b/deal.II/deal.II/include/fe/mapping_q1.h @@ -295,7 +295,7 @@ class MappingQ1 : public Mapping virtual void fill_fe_values (const typename Triangulation::cell_iterator &cell, const Quadrature &quadrature, - const enum CellSimilarity::Similarity cell_similarity, + enum CellSimilarity::Similarity &cell_similarity, typename Mapping::InternalDataBase &mapping_data, typename std::vector > &quadrature_points, std::vector &JxW_values, diff --git a/deal.II/deal.II/include/fe/mapping_q1_eulerian.h b/deal.II/deal.II/include/fe/mapping_q1_eulerian.h index 540c377bfc..ff2635e7cf 100644 --- a/deal.II/deal.II/include/fe/mapping_q1_eulerian.h +++ b/deal.II/deal.II/include/fe/mapping_q1_eulerian.h @@ -152,7 +152,7 @@ class MappingQ1Eulerian : public MappingQ1 virtual void fill_fe_values (const typename Triangulation::cell_iterator &cell, const Quadrature &quadrature, - enum CellSimilarity::Similarity cell_similarity, + enum CellSimilarity::Similarity &cell_similarity, typename Mapping::InternalDataBase &mapping_data, typename std::vector > &quadrature_points, std::vector &JxW_values, diff --git a/deal.II/deal.II/include/fe/mapping_q_eulerian.h b/deal.II/deal.II/include/fe/mapping_q_eulerian.h index 809bab67f8..c1c5454e63 100644 --- a/deal.II/deal.II/include/fe/mapping_q_eulerian.h +++ b/deal.II/deal.II/include/fe/mapping_q_eulerian.h @@ -155,7 +155,7 @@ class MappingQEulerian : public MappingQ virtual void fill_fe_values (const typename Triangulation::cell_iterator &cell, const Quadrature &quadrature, - enum CellSimilarity::Similarity cell_similarity, + enum CellSimilarity::Similarity &cell_similarity, typename Mapping::InternalDataBase &mapping_data, typename std::vector > &quadrature_points, std::vector &JxW_values, diff --git a/deal.II/deal.II/source/fe/fe_values.cc b/deal.II/deal.II/source/fe/fe_values.cc index c66288f97f..b3ff227694 100644 --- a/deal.II/deal.II/source/fe/fe_values.cc +++ b/deal.II/deal.II/source/fe/fe_values.cc @@ -3249,6 +3249,16 @@ FEValuesBase::check_cell_similarity // case that there has not been any cell // before if (&*this->present_cell == 0) + { + cell_similarity = CellSimilarity::none; + return; + } + + // in MappingQ, data can have been + // modified during the last call. Then, + // we can't use that data on the new + // cell. + if (cell_similarity == CellSimilarity::invalid_next_cell) { cell_similarity = CellSimilarity::none; return; @@ -3273,7 +3283,7 @@ FEValuesBase::check_cell_similarity 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) + if (dist_new.norm_square() > 1e-30) { is_translation = false; break; diff --git a/deal.II/deal.II/source/fe/mapping_cartesian.cc b/deal.II/deal.II/source/fe/mapping_cartesian.cc index db2c77f370..c5e55981e2 100644 --- a/deal.II/deal.II/source/fe/mapping_cartesian.cc +++ b/deal.II/deal.II/source/fe/mapping_cartesian.cc @@ -185,22 +185,25 @@ MappingCartesian::compute_fill (const typename Triangulationvertex(1)(0) - start(0); - break; - case 2: - data.length[0] = cell->vertex(1)(0) - start(0); - data.length[1] = cell->vertex(2)(1) - start(1); - break; - case 3: - data.length[0] = cell->vertex(1)(0) - start(0); - data.length[1] = cell->vertex(2)(1) - start(1); - data.length[2] = cell->vertex(4)(2) - start(2); - break; - default: - Assert(false, ExcNotImplemented()); + switch (dim) + { + case 1: + data.length[0] = cell->vertex(1)(0) - start(0); + break; + case 2: + data.length[0] = cell->vertex(1)(0) - start(0); + data.length[1] = cell->vertex(2)(1) - start(1); + break; + case 3: + data.length[0] = cell->vertex(1)(0) - start(0); + data.length[1] = cell->vertex(2)(1) - start(1); + data.length[2] = cell->vertex(4)(2) - start(2); + break; + default: + Assert(false, ExcNotImplemented()); + } } @@ -300,7 +303,7 @@ void MappingCartesian:: fill_fe_values (const typename Triangulation::cell_iterator& cell, const Quadrature& q, - const enum CellSimilarity::Similarity cell_similarity, + enum CellSimilarity::Similarity& cell_similarity, typename Mapping::InternalDataBase& mapping_data, std::vector >& quadrature_points, std::vector& JxW_values, @@ -347,7 +350,7 @@ fill_fe_values (const typename Triangulation::cell_iterator& cell, { jacobians[i]=Tensor<2,dim>(); for (unsigned int j=0; j::fill_fe_values ( const typename Triangulation::cell_iterator &cell, const Quadrature &q, - enum CellSimilarity::Similarity cell_similarity, + enum CellSimilarity::Similarity &cell_similarity, typename Mapping::InternalDataBase &mapping_data, std::vector > &quadrature_points, std::vector &JxW_values, @@ -341,17 +341,17 @@ MappingQ::fill_fe_values ( // recomputed anyway since then the // mapping changes the data. this needs // to be known also for later operations, - // so modify the variable here. TODO: - // still need to check how this will work - // when the previous cell disabled this - // flag. + // so modify the variable here. this also + // affects the calculation of the next + // cell -- if we use Q1 data on the next + // cell, the data will still be invalid. typename MappingQ1::InternalData *p_data=0; if (data.use_mapping_q1_on_current_cell) p_data=&data.mapping_q1_data; else { p_data=&data; - cell_similarity = CellSimilarity::none; + cell_similarity = CellSimilarity::invalid_next_cell; } 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 c5985e9b5d..90741b0ad7 100644 --- a/deal.II/deal.II/source/fe/mapping_q1.cc +++ b/deal.II/deal.II/source/fe/mapping_q1.cc @@ -739,7 +739,7 @@ void MappingQ1::fill_fe_values ( const typename Triangulation::cell_iterator &cell, const Quadrature &q, - const enum CellSimilarity::Similarity cell_similarity, + enum CellSimilarity::Similarity &cell_similarity, typename Mapping::InternalDataBase &mapping_data, std::vector > &quadrature_points, std::vector &JxW_values, 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 3e1bd7f832..d6dfaf9e66 100644 --- a/deal.II/deal.II/source/fe/mapping_q1_eulerian.cc +++ b/deal.II/deal.II/source/fe/mapping_q1_eulerian.cc @@ -120,7 +120,7 @@ void MappingQ1Eulerian::fill_fe_values ( const typename Triangulation::cell_iterator &cell, const Quadrature &q, - enum CellSimilarity::Similarity cell_similarity, + enum CellSimilarity::Similarity &cell_similarity, typename Mapping::InternalDataBase &mapping_data, std::vector > &quadrature_points, std::vector &JxW_values, @@ -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::none; + cell_similarity = CellSimilarity::invalid_next_cell; 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 cff88e4ba0..1431787d12 100644 --- a/deal.II/deal.II/source/fe/mapping_q_eulerian.cc +++ b/deal.II/deal.II/source/fe/mapping_q_eulerian.cc @@ -187,7 +187,7 @@ void MappingQEulerian::fill_fe_values ( const typename Triangulation::cell_iterator &cell, const Quadrature &q, - enum CellSimilarity::Similarity cell_similarity, + enum CellSimilarity::Similarity &cell_similarity, typename Mapping::InternalDataBase &mapping_data, std::vector > &quadrature_points, std::vector &JxW_values, @@ -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::none; + cell_similarity = CellSimilarity::invalid_next_cell; MappingQ::fill_fe_values (cell, q, cell_similarity, mapping_data, quadrature_points, JxW_values, jacobians, jacobian_grads, inverse_jacobians,