From: Marc Fehling Date: Fri, 6 Nov 2020 19:23:12 +0000 (-0700) Subject: Deprecated DoFHandlerType in MappingFEField. X-Git-Tag: v9.3.0-rc1~864^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=7076d15470abe416bc1c00097a0c8e770a08fff9;p=dealii.git Deprecated DoFHandlerType in MappingFEField. --- diff --git a/doc/news/changes/incompatibilities/20201115Fehling b/doc/news/changes/incompatibilities/20201115Fehling new file mode 100644 index 0000000000..da8f19c876 --- /dev/null +++ b/doc/news/changes/incompatibilities/20201115Fehling @@ -0,0 +1,5 @@ +Deprecated: The DoFHandlerType template argument for the MappingFEField +class has been deprecated. Use MappingFEField +instead. +
+(Marc Fehling, 2020/11/15) diff --git a/examples/step-60/step-60.cc b/examples/step-60/step-60.cc index d9885590b0..acfe56d454 100644 --- a/examples/step-60/step-60.cc +++ b/examples/step-60/step-60.cc @@ -677,10 +677,7 @@ namespace Step60 embedded_configuration); else embedded_mapping = - std::make_unique, - DoFHandler>>( + std::make_unique>>( *embedded_configuration_dh, embedded_configuration); setup_embedded_dofs(); diff --git a/include/deal.II/fe/mapping_fe_field.h b/include/deal.II/fe/mapping_fe_field.h index 5c2eae1d94..8814ae17d4 100644 --- a/include/deal.II/fe/mapping_fe_field.h +++ b/include/deal.II/fe/mapping_fe_field.h @@ -39,6 +39,56 @@ DEAL_II_NAMESPACE_OPEN /*!@addtogroup mapping */ /*@{*/ +/** + * @deprecated Use MappingFEField instead. + */ +template , + typename DoFHandlerType = void> +class MappingFEField; + +#ifndef DOXYGEN +// prevent doxygen from complaining about potential recursive class relations +template +class MappingFEField : public MappingFEField +{ +public: + DEAL_II_DEPRECATED + MappingFEField(const DoFHandlerType &euler_dof_handler, + const VectorType & euler_vector, + const ComponentMask & mask = ComponentMask()) + : MappingFEField(euler_dof_handler, + euler_vector, + mask) + {} + + DEAL_II_DEPRECATED + MappingFEField(const DoFHandlerType & euler_dof_handler, + const std::vector &euler_vector, + const ComponentMask & mask = ComponentMask()) + : MappingFEField(euler_dof_handler, + euler_vector, + mask) + {} + + DEAL_II_DEPRECATED + MappingFEField(const DoFHandlerType & euler_dof_handler, + const MGLevelObject &euler_vector, + const ComponentMask & mask = ComponentMask()) + : MappingFEField(euler_dof_handler, + euler_vector, + mask) + {} + + DEAL_II_DEPRECATED + MappingFEField( + const MappingFEField &mapping) + : MappingFEField(mapping) + {} +}; +#endif // DOXYGEN + /** * The MappingFEField is a generalization of the MappingQEulerian class, for * arbitrary vector finite elements. The two main differences are that this @@ -52,7 +102,7 @@ DEAL_II_NAMESPACE_OPEN * can be arbitrarily selected at construction time. * * The idea is to consider the Triangulation as a parameter configuration - * space, on which we construct an arbitrary geometrical mapping, using the + * space, on which we construct an arbitrary geometrical mapping, using the * instruments of the deal.II library: a vector of degrees of freedom, a * DoFHandler associated to the geometry of the problem and a ComponentMask * that tells us which components of the FiniteElement to use for the mapping. @@ -74,11 +124,9 @@ DEAL_II_NAMESPACE_OPEN * MappingFEField map(dhq, eulerq, mask); * @endcode */ -template , - typename DoFHandlerType = DoFHandler> -class MappingFEField : public Mapping +template +class MappingFEField + : public Mapping { public: /** @@ -113,9 +161,9 @@ public: * * If an incompatible mask is passed, an exception is thrown. */ - MappingFEField(const DoFHandlerType &euler_dof_handler, - const VectorType & euler_vector, - const ComponentMask & mask = ComponentMask()); + MappingFEField(const DoFHandler &euler_dof_handler, + const VectorType & euler_vector, + const ComponentMask & mask = ComponentMask()); /** * Constructor taking vectors on the multigrid levels rather than the active @@ -126,9 +174,9 @@ public: * has been called. Apart from the level vectors, the same arguments as in * the other constructor need to be provided. */ - MappingFEField(const DoFHandlerType & euler_dof_handler, - const std::vector &euler_vector, - const ComponentMask & mask = ComponentMask()); + MappingFEField(const DoFHandler &euler_dof_handler, + const std::vector & euler_vector, + const ComponentMask & mask = ComponentMask()); /** * Constructor with MGLevelObject instead of std::vector, otherwise the same @@ -137,7 +185,7 @@ public: * zero or more — it only needs to be consistent between what is set * here and later used for evaluation of the mapping. */ - MappingFEField(const DoFHandlerType & euler_dof_handler, + MappingFEField(const DoFHandler &euler_dof_handler, const MGLevelObject &euler_vector, const ComponentMask & mask = ComponentMask()); @@ -145,7 +193,7 @@ public: * Copy constructor. */ MappingFEField( - const MappingFEField &mapping); + const MappingFEField &mapping); /** * Return a pointer to a copy of the present object. The caller of this copy @@ -537,16 +585,15 @@ private: /** * Reference to the vector of shifts. */ - std::vector< - SmartPointer>> + std::vector>> euler_vector; /** * Pointer to the DoFHandler to which the mapping vector is associated. */ - SmartPointer> + SmartPointer, + MappingFEField> euler_dof_handler; private: @@ -596,8 +643,8 @@ private: void update_internal_dofs( const typename Triangulation::cell_iterator &cell, - const typename MappingFEField:: - InternalData &data) const; + const typename MappingFEField::InternalData + &data) const; /** * See the documentation of the base class for detailed information. @@ -605,8 +652,8 @@ private: virtual void compute_shapes_virtual( const std::vector> &unit_points, - typename MappingFEField:: - InternalData &data) const; + typename MappingFEField::InternalData + &data) const; /* * Which components to use for the mapping. diff --git a/source/fe/mapping_fe_field.cc b/source/fe/mapping_fe_field.cc index 50667b2f51..2da4c2ff17 100644 --- a/source/fe/mapping_fe_field.cc +++ b/source/fe/mapping_fe_field.cc @@ -58,10 +58,10 @@ DEAL_II_NAMESPACE_OPEN -template -MappingFEField::InternalData:: - InternalData(const FiniteElement &fe, - const ComponentMask & mask) +template +MappingFEField::InternalData::InternalData( + const FiniteElement &fe, + const ComponentMask & mask) : unit_tangentials() , n_shape_functions(fe.n_dofs_per_cell()) , mask(mask) @@ -71,9 +71,9 @@ MappingFEField::InternalData:: -template +template std::size_t -MappingFEField::InternalData:: +MappingFEField::InternalData:: memory_consumption() const { Assert(false, ExcNotImplemented()); @@ -82,9 +82,9 @@ MappingFEField::InternalData:: -template +template double & -MappingFEField::InternalData::shape( +MappingFEField::InternalData::shape( const unsigned int qpoint, const unsigned int shape_nr) { @@ -93,10 +93,11 @@ MappingFEField::InternalData::shape( } -template +template const Tensor<1, dim> & -MappingFEField::InternalData:: - derivative(const unsigned int qpoint, const unsigned int shape_nr) const +MappingFEField::InternalData::derivative( + const unsigned int qpoint, + const unsigned int shape_nr) const { AssertIndexRange(qpoint * n_shape_functions + shape_nr, shape_derivatives.size()); @@ -105,10 +106,11 @@ MappingFEField::InternalData:: -template +template Tensor<1, dim> & -MappingFEField::InternalData:: - derivative(const unsigned int qpoint, const unsigned int shape_nr) +MappingFEField::InternalData::derivative( + const unsigned int qpoint, + const unsigned int shape_nr) { AssertIndexRange(qpoint * n_shape_functions + shape_nr, shape_derivatives.size()); @@ -116,9 +118,9 @@ MappingFEField::InternalData:: } -template +template const Tensor<2, dim> & -MappingFEField::InternalData:: +MappingFEField::InternalData:: second_derivative(const unsigned int qpoint, const unsigned int shape_nr) const { @@ -129,9 +131,9 @@ MappingFEField::InternalData:: -template +template Tensor<2, dim> & -MappingFEField::InternalData:: +MappingFEField::InternalData:: second_derivative(const unsigned int qpoint, const unsigned int shape_nr) { AssertIndexRange(qpoint * n_shape_functions + shape_nr, @@ -140,10 +142,11 @@ MappingFEField::InternalData:: } -template +template const Tensor<3, dim> & -MappingFEField::InternalData:: - third_derivative(const unsigned int qpoint, const unsigned int shape_nr) const +MappingFEField::InternalData::third_derivative( + const unsigned int qpoint, + const unsigned int shape_nr) const { AssertIndexRange(qpoint * n_shape_functions + shape_nr, shape_third_derivatives.size()); @@ -152,10 +155,11 @@ MappingFEField::InternalData:: -template +template Tensor<3, dim> & -MappingFEField::InternalData:: - third_derivative(const unsigned int qpoint, const unsigned int shape_nr) +MappingFEField::InternalData::third_derivative( + const unsigned int qpoint, + const unsigned int shape_nr) { AssertIndexRange(qpoint * n_shape_functions + shape_nr, shape_third_derivatives.size()); @@ -163,9 +167,9 @@ MappingFEField::InternalData:: } -template +template const Tensor<4, dim> & -MappingFEField::InternalData:: +MappingFEField::InternalData:: fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) const { @@ -176,9 +180,9 @@ MappingFEField::InternalData:: -template +template Tensor<4, dim> & -MappingFEField::InternalData:: +MappingFEField::InternalData:: fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) { AssertIndexRange(qpoint * n_shape_functions + shape_nr, @@ -210,11 +214,11 @@ namespace -template -MappingFEField::MappingFEField( - const DoFHandlerType &euler_dof_handler, - const VectorType & euler_vector, - const ComponentMask & mask) +template +MappingFEField::MappingFEField( + const DoFHandler &euler_dof_handler, + const VectorType & euler_vector, + const ComponentMask & mask) : uses_level_dofs(false) , euler_vector({&euler_vector}) , euler_dof_handler(&euler_dof_handler) @@ -239,11 +243,11 @@ MappingFEField::MappingFEField( -template -MappingFEField::MappingFEField( - const DoFHandlerType & euler_dof_handler, - const std::vector &euler_vector, - const ComponentMask & mask) +template +MappingFEField::MappingFEField( + const DoFHandler &euler_dof_handler, + const std::vector & euler_vector, + const ComponentMask & mask) : uses_level_dofs(true) , euler_dof_handler(&euler_dof_handler) , fe_mask(mask.size() ? @@ -278,9 +282,9 @@ MappingFEField::MappingFEField( -template -MappingFEField::MappingFEField( - const DoFHandlerType & euler_dof_handler, +template +MappingFEField::MappingFEField( + const DoFHandler &euler_dof_handler, const MGLevelObject &euler_vector, const ComponentMask & mask) : uses_level_dofs(true) @@ -319,9 +323,9 @@ MappingFEField::MappingFEField( -template -MappingFEField::MappingFEField( - const MappingFEField &mapping) +template +MappingFEField::MappingFEField( + const MappingFEField &mapping) : uses_level_dofs(mapping.uses_level_dofs) , euler_vector(mapping.euler_vector) , euler_dof_handler(mapping.euler_dof_handler) @@ -334,9 +338,9 @@ MappingFEField::MappingFEField( -template +template inline const double & -MappingFEField::InternalData::shape( +MappingFEField::InternalData::shape( const unsigned int qpoint, const unsigned int shape_nr) const { @@ -346,20 +350,20 @@ MappingFEField::InternalData::shape( -template +template bool -MappingFEField:: - preserves_vertex_locations() const +MappingFEField::preserves_vertex_locations() + const { return false; } -template +template boost::container::small_vector, GeometryInfo::vertices_per_cell> -MappingFEField::get_vertices( +MappingFEField::get_vertices( const typename Triangulation::cell_iterator &cell) const { // we transform our tria iterator into a dof iterator so we can access @@ -420,13 +424,12 @@ MappingFEField::get_vertices( -template +template void -MappingFEField:: - compute_shapes_virtual( - const std::vector> &unit_points, - typename MappingFEField:: - InternalData &data) const +MappingFEField::compute_shapes_virtual( + const std::vector> &unit_points, + typename MappingFEField::InternalData &data) + const { const auto fe = &euler_dof_handler->get_fe(); const unsigned int n_points = unit_points.size(); @@ -459,10 +462,10 @@ MappingFEField:: } -template +template UpdateFlags -MappingFEField:: - requires_update_flags(const UpdateFlags in) const +MappingFEField::requires_update_flags( + const UpdateFlags in) const { // add flags if the respective quantities are necessary to compute // what we need. note that some flags appear in both conditions and @@ -514,9 +517,9 @@ MappingFEField:: -template +template void -MappingFEField::compute_data( +MappingFEField::compute_data( const UpdateFlags update_flags, const Quadrature &q, const unsigned int n_original_q_points, @@ -565,9 +568,9 @@ MappingFEField::compute_data( } -template +template void -MappingFEField::compute_face_data( +MappingFEField::compute_face_data( const UpdateFlags update_flags, const Quadrature &q, const unsigned int n_original_q_points, @@ -606,9 +609,9 @@ MappingFEField::compute_face_data( } -template +template typename std::unique_ptr::InternalDataBase> -MappingFEField::get_data( +MappingFEField::get_data( const UpdateFlags update_flags, const Quadrature &quadrature) const { @@ -622,9 +625,9 @@ MappingFEField::get_data( -template +template std::unique_ptr::InternalDataBase> -MappingFEField::get_face_data( +MappingFEField::get_face_data( const UpdateFlags update_flags, const Quadrature &quadrature) const { @@ -640,9 +643,9 @@ MappingFEField::get_face_data( } -template +template std::unique_ptr::InternalDataBase> -MappingFEField::get_subface_data( +MappingFEField::get_subface_data( const UpdateFlags update_flags, const Quadrature &quadrature) const { @@ -671,16 +674,12 @@ namespace internal * have already been set), but only if the update_flags of the @p data * argument indicate so. */ - template + template void maybe_compute_q_points( const typename dealii::QProjector::DataSetDescriptor data_set, - const typename dealii:: - MappingFEField:: - InternalData & data, + const typename dealii::MappingFEField:: + InternalData & data, const FiniteElement &fe, const ComponentMask & fe_mask, const std::vector & fe_to_real, @@ -717,16 +716,12 @@ namespace internal * * Skip the computation if possible as indicated by the first argument. */ - template + template void maybe_update_Jacobians( const typename dealii::QProjector::DataSetDescriptor data_set, - const typename dealii:: - MappingFEField:: - InternalData & data, + const typename dealii::MappingFEField:: + InternalData & data, const FiniteElement &fe, const ComponentMask & fe_mask, const std::vector & fe_to_real) @@ -790,16 +785,12 @@ namespace internal * * Skip the computation if possible as indicated by the first argument. */ - template + template void maybe_update_jacobian_grads( const typename dealii::QProjector::DataSetDescriptor data_set, - const typename dealii:: - MappingFEField:: - InternalData & data, + const typename dealii::MappingFEField:: + InternalData & data, const FiniteElement & fe, const ComponentMask & fe_mask, const std::vector & fe_to_real, @@ -844,16 +835,12 @@ namespace internal * * Skip the computation if possible as indicated by the first argument. */ - template + template void maybe_update_jacobian_pushed_forward_grads( const typename dealii::QProjector::DataSetDescriptor data_set, - const typename dealii:: - MappingFEField:: - InternalData & data, + const typename dealii::MappingFEField:: + InternalData & data, const FiniteElement &fe, const ComponentMask & fe_mask, const std::vector & fe_to_real, @@ -921,16 +908,12 @@ namespace internal * * Skip the computation if possible as indicated by the first argument. */ - template + template void maybe_update_jacobian_2nd_derivatives( const typename dealii::QProjector::DataSetDescriptor data_set, - const typename dealii:: - MappingFEField:: - InternalData & data, + const typename dealii::MappingFEField:: + InternalData & data, const FiniteElement & fe, const ComponentMask & fe_mask, const std::vector & fe_to_real, @@ -979,16 +962,12 @@ namespace internal * * Skip the computation if possible as indicated by the first argument. */ - template + template void maybe_update_jacobian_pushed_forward_2nd_derivatives( const typename dealii::QProjector::DataSetDescriptor data_set, - const typename dealii:: - MappingFEField:: - InternalData & data, + const typename dealii::MappingFEField:: + InternalData & data, const FiniteElement &fe, const ComponentMask & fe_mask, const std::vector & fe_to_real, @@ -1079,16 +1058,12 @@ namespace internal * * Skip the computation if possible as indicated by the first argument. */ - template + template void maybe_update_jacobian_3rd_derivatives( const typename dealii::QProjector::DataSetDescriptor data_set, - const typename dealii:: - MappingFEField:: - InternalData & data, + const typename dealii::MappingFEField:: + InternalData & data, const FiniteElement & fe, const ComponentMask & fe_mask, const std::vector & fe_to_real, @@ -1141,16 +1116,12 @@ namespace internal * * Skip the computation if possible as indicated by the first argument. */ - template + template void maybe_update_jacobian_pushed_forward_3rd_derivatives( const typename dealii::QProjector::DataSetDescriptor data_set, - const typename dealii:: - MappingFEField:: - InternalData & data, + const typename dealii::MappingFEField:: + InternalData & data, const FiniteElement &fe, const ComponentMask & fe_mask, const std::vector & fe_to_real, @@ -1271,10 +1242,7 @@ namespace internal * * The resulting data is put into the @p output_data argument. */ - template + template void maybe_compute_face_data( const dealii::Mapping &mapping, @@ -1283,9 +1251,8 @@ namespace internal const unsigned int face_no, const unsigned int subface_no, const std::vector &weights, - const typename dealii:: - MappingFEField:: - InternalData &data, + const typename dealii::MappingFEField:: + InternalData &data, internal::FEValuesImplementation::MappingRelatedData &output_data) { @@ -1432,10 +1399,7 @@ namespace internal * 'data_set' to differentiate whether we will work on a face (and if so, * which one) or subface. */ - template + template void do_fill_fe_face_values( const dealii::Mapping &mapping, @@ -1445,16 +1409,15 @@ namespace internal const unsigned int subface_no, const typename dealii::QProjector::DataSetDescriptor data_set, const Quadrature & quadrature, - const typename dealii:: - MappingFEField:: - InternalData & data, + const typename dealii::MappingFEField:: + InternalData & data, const FiniteElement &fe, const ComponentMask & fe_mask, const std::vector & fe_to_real, internal::FEValuesImplementation::MappingRelatedData &output_data) { - maybe_compute_q_points( + maybe_compute_q_points( data_set, data, fe, @@ -1462,16 +1425,13 @@ namespace internal fe_to_real, output_data.quadrature_points); - maybe_update_Jacobians( + maybe_update_Jacobians( data_set, data, fe, fe_mask, fe_to_real); - maybe_update_jacobian_grads( + maybe_update_jacobian_grads( data_set, data, fe, fe_mask, fe_to_real, output_data.jacobian_grads); - maybe_update_jacobian_pushed_forward_grads( + maybe_update_jacobian_pushed_forward_grads( data_set, data, fe, @@ -1479,10 +1439,7 @@ namespace internal fe_to_real, output_data.jacobian_pushed_forward_grads); - maybe_update_jacobian_2nd_derivatives( + maybe_update_jacobian_2nd_derivatives( data_set, data, fe, @@ -1492,8 +1449,7 @@ namespace internal maybe_update_jacobian_pushed_forward_2nd_derivatives( + VectorType>( data_set, data, fe, @@ -1501,10 +1457,7 @@ namespace internal fe_to_real, output_data.jacobian_pushed_forward_2nd_derivatives); - maybe_update_jacobian_3rd_derivatives( + maybe_update_jacobian_3rd_derivatives( data_set, data, fe, @@ -1514,8 +1467,7 @@ namespace internal maybe_update_jacobian_pushed_forward_3rd_derivatives( + VectorType>( data_set, data, fe, @@ -1523,7 +1475,7 @@ namespace internal fe_to_real, output_data.jacobian_pushed_forward_3rd_derivatives); - maybe_compute_face_data( + maybe_compute_face_data( mapping, cell, face_no, @@ -1539,9 +1491,9 @@ namespace internal // Note that the CellSimilarity flag is modifiable, since MappingFEField can // need to recalculate data even when cells are similar. -template +template CellSimilarity::Similarity -MappingFEField::fill_fe_values( +MappingFEField::fill_fe_values( const typename Triangulation::cell_iterator &cell, const CellSimilarity::Similarity, const Quadrature & quadrature, @@ -1560,7 +1512,7 @@ MappingFEField::fill_fe_values( update_internal_dofs(cell, data); internal::MappingFEFieldImplementation:: - maybe_compute_q_points( + maybe_compute_q_points( QProjector::DataSetDescriptor::cell(), data, euler_dof_handler->get_fe(), @@ -1569,7 +1521,7 @@ MappingFEField::fill_fe_values( output_data.quadrature_points); internal::MappingFEFieldImplementation:: - maybe_update_Jacobians( + maybe_update_Jacobians( QProjector::DataSetDescriptor::cell(), data, euler_dof_handler->get_fe(), @@ -1668,7 +1620,7 @@ MappingFEField::fill_fe_values( // calculate derivatives of the Jacobians internal::MappingFEFieldImplementation:: - maybe_update_jacobian_grads( + maybe_update_jacobian_grads( QProjector::DataSetDescriptor::cell(), data, euler_dof_handler->get_fe(), @@ -1679,10 +1631,7 @@ MappingFEField::fill_fe_values( // calculate derivatives of the Jacobians pushed forward to real cell // coordinates internal::MappingFEFieldImplementation:: - maybe_update_jacobian_pushed_forward_grads( + maybe_update_jacobian_pushed_forward_grads( QProjector::DataSetDescriptor::cell(), data, euler_dof_handler->get_fe(), @@ -1691,23 +1640,20 @@ MappingFEField::fill_fe_values( output_data.jacobian_pushed_forward_grads); // calculate hessians of the Jacobians - internal::MappingFEFieldImplementation::maybe_update_jacobian_2nd_derivatives< - dim, - spacedim, - VectorType, - DoFHandlerType>(QProjector::DataSetDescriptor::cell(), - data, - euler_dof_handler->get_fe(), - fe_mask, - fe_to_real, - output_data.jacobian_2nd_derivatives); + internal::MappingFEFieldImplementation:: + maybe_update_jacobian_2nd_derivatives( + QProjector::DataSetDescriptor::cell(), + data, + euler_dof_handler->get_fe(), + fe_mask, + fe_to_real, + output_data.jacobian_2nd_derivatives); // calculate hessians of the Jacobians pushed forward to real cell coordinates internal::MappingFEFieldImplementation:: maybe_update_jacobian_pushed_forward_2nd_derivatives( + VectorType>( QProjector::DataSetDescriptor::cell(), data, euler_dof_handler->get_fe(), @@ -1716,24 +1662,21 @@ MappingFEField::fill_fe_values( output_data.jacobian_pushed_forward_2nd_derivatives); // calculate gradients of the hessians of the Jacobians - internal::MappingFEFieldImplementation::maybe_update_jacobian_3rd_derivatives< - dim, - spacedim, - VectorType, - DoFHandlerType>(QProjector::DataSetDescriptor::cell(), - data, - euler_dof_handler->get_fe(), - fe_mask, - fe_to_real, - output_data.jacobian_3rd_derivatives); + internal::MappingFEFieldImplementation:: + maybe_update_jacobian_3rd_derivatives( + QProjector::DataSetDescriptor::cell(), + data, + euler_dof_handler->get_fe(), + fe_mask, + fe_to_real, + output_data.jacobian_3rd_derivatives); // calculate gradients of the hessians of the Jacobians pushed forward to real // cell coordinates internal::MappingFEFieldImplementation:: maybe_update_jacobian_pushed_forward_3rd_derivatives( + VectorType>( QProjector::DataSetDescriptor::cell(), data, euler_dof_handler->get_fe(), @@ -1746,9 +1689,9 @@ MappingFEField::fill_fe_values( -template +template void -MappingFEField::fill_fe_face_values( +MappingFEField::fill_fe_face_values( const typename Triangulation::cell_iterator &cell, const unsigned int face_no, const Quadrature & quadrature, @@ -1764,39 +1707,38 @@ MappingFEField::fill_fe_face_values( update_internal_dofs(cell, data); - internal::MappingFEFieldImplementation:: - do_fill_fe_face_values( - *this, - cell, - face_no, - numbers::invalid_unsigned_int, - QProjector::DataSetDescriptor::face(ReferenceCell::get_hypercube( - dim), - face_no, - cell->face_orientation(face_no), - cell->face_flip(face_no), - cell->face_rotation(face_no), - quadrature.size()), - quadrature, - data, - euler_dof_handler->get_fe(), - fe_mask, - fe_to_real, - output_data); + internal::MappingFEFieldImplementation::do_fill_fe_face_values( + *this, + cell, + face_no, + numbers::invalid_unsigned_int, + QProjector::DataSetDescriptor::face(ReferenceCell::get_hypercube(dim), + face_no, + cell->face_orientation(face_no), + cell->face_flip(face_no), + cell->face_rotation(face_no), + quadrature.size()), + quadrature, + data, + euler_dof_handler->get_fe(), + fe_mask, + fe_to_real, + output_data); } -template +template void -MappingFEField:: - fill_fe_subface_values( - const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const unsigned int subface_no, - const Quadrature & quadrature, - const typename Mapping::InternalDataBase & internal_data, - internal::FEValuesImplementation::MappingRelatedData - &output_data) const +MappingFEField::fill_fe_subface_values( + const typename Triangulation::cell_iterator &cell, + const unsigned int face_no, + const unsigned int subface_no, + const Quadrature & quadrature, + const typename Mapping::InternalDataBase & internal_data, + internal::FEValuesImplementation::MappingRelatedData + &output_data) const { // convert data object to internal data for this class. fails with an // exception if that is not possible @@ -1806,27 +1748,28 @@ MappingFEField:: update_internal_dofs(cell, data); - internal::MappingFEFieldImplementation:: - do_fill_fe_face_values( - *this, - cell, - face_no, - numbers::invalid_unsigned_int, - QProjector::DataSetDescriptor::subface( - ReferenceCell::get_hypercube(dim), - face_no, - subface_no, - cell->face_orientation(face_no), - cell->face_flip(face_no), - cell->face_rotation(face_no), - quadrature.size(), - cell->subface_case(face_no)), - quadrature, - data, - euler_dof_handler->get_fe(), - fe_mask, - fe_to_real, - output_data); + internal::MappingFEFieldImplementation::do_fill_fe_face_values( + *this, + cell, + face_no, + numbers::invalid_unsigned_int, + QProjector::DataSetDescriptor::subface(ReferenceCell::get_hypercube( + dim), + face_no, + subface_no, + cell->face_orientation(face_no), + cell->face_flip(face_no), + cell->face_rotation(face_no), + quadrature.size(), + cell->subface_case(face_no)), + quadrature, + data, + euler_dof_handler->get_fe(), + fe_mask, + fe_to_real, + output_data); } @@ -1836,11 +1779,7 @@ namespace internal { namespace { - template + template void transform_fields( const ArrayView> & input, @@ -1849,19 +1788,17 @@ namespace internal const ArrayView> & output) { AssertDimension(input.size(), output.size()); - Assert((dynamic_cast< - const typename dealii:: - MappingFEField:: - InternalData *>(&mapping_data) != nullptr), - ExcInternalError()); - const typename dealii::MappingFEField::InternalData - &data = static_cast< + Assert( + (dynamic_cast< + const typename dealii:: + MappingFEField::InternalData *>( + &mapping_data) != nullptr), + ExcInternalError()); + const typename dealii::MappingFEField:: + InternalData &data = static_cast< const typename dealii:: - MappingFEField:: - InternalData &>(mapping_data); + MappingFEField::InternalData &>( + mapping_data); switch (mapping_kind) { @@ -1922,11 +1859,7 @@ namespace internal } - template + template void transform_differential_forms( const ArrayView> &input, @@ -1935,19 +1868,17 @@ namespace internal const ArrayView> & output) { AssertDimension(input.size(), output.size()); - Assert((dynamic_cast< - const typename dealii:: - MappingFEField:: - InternalData *>(&mapping_data) != nullptr), - ExcInternalError()); - const typename dealii::MappingFEField::InternalData - &data = static_cast< + Assert( + (dynamic_cast< + const typename dealii:: + MappingFEField::InternalData *>( + &mapping_data) != nullptr), + ExcInternalError()); + const typename dealii::MappingFEField:: + InternalData &data = static_cast< const typename dealii:: - MappingFEField:: - InternalData &>(mapping_data); + MappingFEField::InternalData &>( + mapping_data); switch (mapping_kind) { @@ -1973,9 +1904,9 @@ namespace internal -template +template void -MappingFEField::transform( +MappingFEField::transform( const ArrayView> & input, const MappingKind mapping_kind, const typename Mapping::InternalDataBase &mapping_data, @@ -1984,17 +1915,17 @@ MappingFEField::transform( AssertDimension(input.size(), output.size()); internal::MappingFEFieldImplementation:: - transform_fields(input, - mapping_kind, - mapping_data, - output); + transform_fields(input, + mapping_kind, + mapping_data, + output); } -template +template void -MappingFEField::transform( +MappingFEField::transform( const ArrayView> &input, const MappingKind mapping_kind, const typename Mapping::InternalDataBase &mapping_data, @@ -2003,15 +1934,17 @@ MappingFEField::transform( AssertDimension(input.size(), output.size()); internal::MappingFEFieldImplementation:: - transform_differential_forms( - input, mapping_kind, mapping_data, output); + transform_differential_forms(input, + mapping_kind, + mapping_data, + output); } -template +template void -MappingFEField::transform( +MappingFEField::transform( const ArrayView> &input, const MappingKind, const typename Mapping::InternalDataBase &mapping_data, @@ -2027,9 +1960,9 @@ MappingFEField::transform( -template +template void -MappingFEField::transform( +MappingFEField::transform( const ArrayView> &input, const MappingKind mapping_kind, const typename Mapping::InternalDataBase &mapping_data, @@ -2075,9 +2008,9 @@ MappingFEField::transform( -template +template void -MappingFEField::transform( +MappingFEField::transform( const ArrayView> &input, const MappingKind /*mapping_kind*/, const typename Mapping::InternalDataBase &mapping_data, @@ -2093,12 +2026,11 @@ MappingFEField::transform( -template +template Point -MappingFEField:: - transform_unit_to_real_cell( - const typename Triangulation::cell_iterator &cell, - const Point & p) const +MappingFEField::transform_unit_to_real_cell( + const typename Triangulation::cell_iterator &cell, + const Point & p) const { // Use the get_data function to create an InternalData with data vectors of // the right size and transformation shape values already computed at point @@ -2115,10 +2047,10 @@ MappingFEField:: } -template +template Point -MappingFEField:: - do_transform_unit_to_real_cell(const InternalData &data) const +MappingFEField::do_transform_unit_to_real_cell( + const InternalData &data) const { Point p_real; @@ -2136,12 +2068,11 @@ MappingFEField:: -template +template Point -MappingFEField:: - transform_real_to_unit_cell( - const typename Triangulation::cell_iterator &cell, - const Point & p) const +MappingFEField::transform_real_to_unit_cell( + const typename Triangulation::cell_iterator &cell, + const Point & p) const { // first a Newton iteration based on the real mapping. It uses the center // point of the cell as a starting point @@ -2185,14 +2116,13 @@ MappingFEField:: } -template +template Point -MappingFEField:: - do_transform_real_to_unit_cell( - const typename Triangulation::cell_iterator &cell, - const Point & p, - const Point & initial_p_unit, - InternalData & mdata) const +MappingFEField::do_transform_real_to_unit_cell( + const typename Triangulation::cell_iterator &cell, + const Point & p, + const Point & initial_p_unit, + InternalData & mdata) const { const unsigned int n_shapes = mdata.shape_values.size(); (void)n_shapes; @@ -2297,43 +2227,43 @@ failure: } -template +template unsigned int -MappingFEField::get_degree() const +MappingFEField::get_degree() const { return euler_dof_handler->get_fe().degree; } -template +template ComponentMask -MappingFEField::get_component_mask() - const +MappingFEField::get_component_mask() const { return this->fe_mask; } -template +template std::unique_ptr> -MappingFEField::clone() const +MappingFEField::clone() const { - return std::make_unique< - MappingFEField>(*this); + return std::make_unique>( + *this); } -template +template void -MappingFEField::update_internal_dofs( +MappingFEField::update_internal_dofs( const typename Triangulation::cell_iterator &cell, - const typename MappingFEField:: - InternalData &data) const + const typename MappingFEField::InternalData + &data) const { Assert(euler_dof_handler != nullptr, ExcMessage("euler_dof_handler is empty")); - typename DoFHandlerType::cell_iterator dof_cell(*cell, euler_dof_handler); + typename DoFHandler::cell_iterator dof_cell(*cell, + euler_dof_handler); Assert(uses_level_dofs || dof_cell->is_active() == true, ExcInactiveCell()); if (uses_level_dofs) { diff --git a/source/fe/mapping_fe_field.inst.in b/source/fe/mapping_fe_field.inst.in index b5820b794b..da2dd79212 100644 --- a/source/fe/mapping_fe_field.inst.in +++ b/source/fe/mapping_fe_field.inst.in @@ -19,10 +19,8 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS; VEC : REAL_VECTOR_TYPES) { #if deal_II_dimension <= deal_II_space_dimension - template class MappingFEField< - deal_II_dimension, - deal_II_space_dimension, - VEC, - dealii::DoFHandler>; + template class MappingFEField; #endif }