From: Luca Heltai Date: Tue, 24 Apr 2018 11:27:56 +0000 (+0200) Subject: Added implementation of MappingFEField::get_vertices X-Git-Tag: v9.0.0-rc1~68^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=c35f04452455a8002b55f3b9787b676d96e8f167;p=dealii.git Added implementation of MappingFEField::get_vertices --- diff --git a/include/deal.II/fe/mapping_fe_field.h b/include/deal.II/fe/mapping_fe_field.h index c7235ef21e..f27680e7d6 100644 --- a/include/deal.II/fe/mapping_fe_field.h +++ b/include/deal.II/fe/mapping_fe_field.h @@ -24,6 +24,7 @@ #include #include #include +#include #include @@ -125,8 +126,6 @@ public: virtual std::unique_ptr> clone () const override; - - /** * See the documentation of Mapping::preserves_vertex_locations() * for the purpose of this function. The implementation in this @@ -135,6 +134,16 @@ public: virtual bool preserves_vertex_locations () const override; + /** + * Return the mapped vertices of a cell. + * + * This mapping ignores the vertices of the Triangulation it is associated to, + * and constructs the position of the vertices according to the @p euler_vector + * that was passed at construction time. + */ + virtual + std::array, GeometryInfo::vertices_per_cell> + get_vertices (const typename Triangulation::cell_iterator &cell) const override; /** * @name Mapping points between reference and real cells @@ -206,7 +215,6 @@ public: * @} */ - /** * Return the degree of the mapping, i.e. the value which was passed to the * constructor. @@ -574,6 +582,17 @@ private: */ std::vector fe_to_real; + /** + * FEValues object used to query the given finite element field at the + * support points in the reference configuration. + */ + mutable FEValues fe_values; + + /** + * A variable to guard access to the fe_values variable. + */ + mutable Threads::Mutex fe_values_mutex; + void compute_data (const UpdateFlags update_flags, const Quadrature &q, @@ -599,9 +618,6 @@ private: #ifndef DOXYGEN - - - #endif // DOXYGEN DEAL_II_NAMESPACE_CLOSE diff --git a/source/fe/mapping_fe_field.cc b/source/fe/mapping_fe_field.cc index 0d339c5fad..3009f8cc8f 100644 --- a/source/fe/mapping_fe_field.cc +++ b/source/fe/mapping_fe_field.cc @@ -202,6 +202,27 @@ MappingFEField::InternalData::fourth_der +namespace +{ + // Construct a quadrature formula containing the vertices of the reference + // cell in dimension dim (with invalid weights) + template + Quadrature &get_vertex_quadrature() + { + static Quadrature quad; + if (quad.size() == 0) + { + std::vector > points(GeometryInfo::vertices_per_cell); + for (unsigned int i=0; i::vertices_per_cell; ++i) + points[i] = GeometryInfo::unit_cell_vertex(i); + quad = Quadrature(points); + } + return quad; + } +} + + + template MappingFEField::MappingFEField (const DoFHandlerType &euler_dof_handler, @@ -212,7 +233,8 @@ MappingFEField::MappingFEField euler_dof_handler(&euler_dof_handler), fe_mask(mask.size() ? mask : ComponentMask(euler_dof_handler.get_fe().get_nonzero_components(0).size(), true)), - fe_to_real(fe_mask.size(), numbers::invalid_unsigned_int) + fe_to_real(fe_mask.size(), numbers::invalid_unsigned_int), + fe_values(this->euler_dof_handler->get_fe(), get_vertex_quadrature(), update_values) { unsigned int size = 0; for (unsigned int i=0; i::MappingFEField euler_vector(mapping.euler_vector), euler_dof_handler(mapping.euler_dof_handler), fe_mask(mapping.fe_mask), - fe_to_real(mapping.fe_to_real) + fe_to_real(mapping.fe_to_real), + fe_values(mapping.euler_dof_handler->get_fe(), get_vertex_quadrature(), update_values) {} -template +template inline const double & -MappingFEField::InternalData::shape +MappingFEField::InternalData::shape (const unsigned int qpoint, const unsigned int shape_nr) const { @@ -251,15 +274,51 @@ MappingFEField::InternalData::shape -template +template bool -MappingFEField::preserves_vertex_locations () const +MappingFEField::preserves_vertex_locations () const { return false; } +template +std::array, GeometryInfo::vertices_per_cell> +MappingFEField:: +get_vertices +(const typename Triangulation::cell_iterator &cell) const +{ + // we transform our tria iterator into a dof iterator so we can access + // data not associated with triangulations + const typename DoFHandler::cell_iterator dof_cell(*cell, + euler_dof_handler); + + Assert (dof_cell->active() == true, ExcInactiveCell()); + AssertDimension(GeometryInfo::vertices_per_cell, fe_values.n_quadrature_points); + AssertDimension(fe_to_real.size(), euler_dof_handler->get_fe().n_components()); + + std::vector > + values(fe_values.n_quadrature_points, + Vector (euler_dof_handler->get_fe().n_components())); + + { + Threads::Mutex::ScopedLock lock(fe_values_mutex); + fe_values.reinit(dof_cell); + fe_values.get_function_values(*euler_vector, values); + } + std::array, GeometryInfo::vertices_per_cell> vertices; + + for (unsigned int i=0; i::vertices_per_cell; ++i) + for (unsigned int j=0; j void MappingFEField:: diff --git a/tests/grid/grid_tools_min_max_diameter_02.cc b/tests/grid/grid_tools_min_max_diameter_02.cc index 8dd263bad8..f653eae67b 100644 --- a/tests/grid/grid_tools_min_max_diameter_02.cc +++ b/tests/grid/grid_tools_min_max_diameter_02.cc @@ -13,6 +13,9 @@ // // --------------------------------------------------------------------- +// test GridTools::minimal_cell_diameter and GridTools::maximal_cell_diameter with a +// MappingFEField + #include "../tests.h" #include #include