From: Marc Fehling Date: Tue, 4 Dec 2018 22:18:31 +0000 (+0100) Subject: FESeries::calculate accepts different VectorType for local dof values. Minor cleanups. X-Git-Tag: v9.1.0-rc1~512^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F7503%2Fhead;p=dealii.git FESeries::calculate accepts different VectorType for local dof values. Minor cleanups. --- diff --git a/doc/news/changes/minor/20181204Fehling b/doc/news/changes/minor/20181204Fehling new file mode 100644 index 0000000000..6fbb83551a --- /dev/null +++ b/doc/news/changes/minor/20181204Fehling @@ -0,0 +1,6 @@ +Changed: FESeries::Fourier::calculate and FESeries::Legendre::calculate +accept both Vector and Vector as input parameters for +local dof values. +
+(Marc Fehling, 2018/12/04) + diff --git a/include/deal.II/fe/fe_series.h b/include/deal.II/fe/fe_series.h index 1b76e67e54..dff8b2d173 100644 --- a/include/deal.II/fe/fe_series.h +++ b/include/deal.II/fe/fe_series.h @@ -91,6 +91,8 @@ namespace FESeries class Fourier : public Subscriptor { public: + using CoefficientType = typename std::complex; + /** * A non-default constructor. The @p size_in_each_direction defines the number * of modes in each direction, @p fe_collection is the hp::FECollection @@ -107,10 +109,11 @@ namespace FESeries * @p local_dof_values corresponding to FiniteElement with * @p cell_active_fe_index . */ + template void - calculate(const dealii::Vector & local_dof_values, - const unsigned int cell_active_fe_index, - Table> &fourier_coefficients); + calculate(const dealii::Vector &local_dof_values, + const unsigned int cell_active_fe_index, + Table & fourier_coefficients); private: /** @@ -131,12 +134,12 @@ namespace FESeries /** * Transformation matrices for each FiniteElement. */ - std::vector>> fourier_transform_matrices; + std::vector> fourier_transform_matrices; /** * Auxiliary vector to store unrolled coefficients. */ - std::vector> unrolled_coefficients; + std::vector unrolled_coefficients; }; /** @@ -186,6 +189,8 @@ namespace FESeries class Legendre : public Subscriptor { public: + using CoefficientType = double; + /** * A non-default constructor. The @p size_in_each_direction defines the number * of coefficients in each direction, @p fe_collection is the hp::FECollection @@ -202,10 +207,11 @@ namespace FESeries * @p local_dof_values corresponding to FiniteElement with * @p cell_active_fe_index . */ + template void - calculate(const dealii::Vector &local_dof_values, + calculate(const dealii::Vector &local_dof_values, const unsigned int cell_active_fe_index, - Table & legendre_coefficients); + Table & legendre_coefficients); private: /** @@ -226,12 +232,12 @@ namespace FESeries /** * Transformation matrices for each FiniteElement. */ - std::vector> legendre_transform_matrices; + std::vector> legendre_transform_matrices; /** * Auxiliary vector to store unrolled coefficients. */ - std::vector unrolled_coefficients; + std::vector unrolled_coefficients; }; @@ -249,12 +255,12 @@ namespace FESeries * in this case: mean, L1_norm, L2_norm, Linfty_norm. The mean norm can only * be applied to real valued coefficients. */ - template + template std::pair, std::vector> - process_coefficients(const Table & coefficients, + process_coefficients(const Table &coefficients, const std::function( - const TableIndices &)> &predicate, - const VectorTools::NormType norm); + const TableIndices &)> & predicate, + const VectorTools::NormType norm); @@ -278,32 +284,34 @@ namespace internal { namespace FESeriesImplementation { - template + template void - fill_map_index(const Table & coefficients, - const TableIndices & ind, - const std::function( - const TableIndices &)> & predicate, - std::map> &pred_to_values) + fill_map_index( + const Table &coefficients, + const TableIndices & ind, + const std::function< + std::pair(const TableIndices &)> &predicate, + std::map> &pred_to_values) { const std::pair pred_pair = predicate(ind); // don't add a value if predicate is false if (pred_pair.first == false) return; - const unsigned int pred_value = pred_pair.second; - const T & coeff_value = coefficients(ind); + const unsigned int pred_value = pred_pair.second; + const CoefficientType &coeff_value = coefficients(ind); // If pred_value is not in the pred_to_values map, the element will be // created. Otherwise a reference to the existing element is returned. pred_to_values[pred_value].push_back(coeff_value); } - template + template void - fill_map(const Table<1, T> & coefficients, - const std::function( - const TableIndices<1> &)> & predicate, - std::map> &pred_to_values) + fill_map( + const Table<1, CoefficientType> &coefficients, + const std::function< + std::pair(const TableIndices<1> &)> &predicate, + std::map> & pred_to_values) { for (unsigned int i = 0; i < coefficients.size(0); i++) { @@ -312,12 +320,13 @@ namespace internal } } - template + template void - fill_map(const Table<2, T> & coefficients, - const std::function( - const TableIndices<2> &)> & predicate, - std::map> &pred_to_values) + fill_map( + const Table<2, CoefficientType> &coefficients, + const std::function< + std::pair(const TableIndices<2> &)> &predicate, + std::map> & pred_to_values) { for (unsigned int i = 0; i < coefficients.size(0); i++) for (unsigned int j = 0; j < coefficients.size(1); j++) @@ -327,12 +336,13 @@ namespace internal } } - template + template void - fill_map(const Table<3, T> & coefficients, - const std::function( - const TableIndices<3> &)> & predicate, - std::map> &pred_to_values) + fill_map( + const Table<3, CoefficientType> &coefficients, + const std::function< + std::pair(const TableIndices<3> &)> &predicate, + std::map> & pred_to_values) { for (unsigned int i = 0; i < coefficients.size(0); i++) for (unsigned int j = 0; j < coefficients.size(1); j++) @@ -344,16 +354,16 @@ namespace internal } - template + template double - complex_mean_value(const T &value) + complex_mean_value(const Number &value) { return value; } - template + template double - complex_mean_value(const std::complex &value) + complex_mean_value(const std::complex &value) { AssertThrow(false, ExcMessage( @@ -365,10 +375,10 @@ namespace internal } // namespace internal -template +template std::pair, std::vector> FESeries::process_coefficients( - const Table &coefficients, + const Table &coefficients, const std::function(const TableIndices &)> & predicate, const VectorTools::NormType norm) @@ -379,19 +389,17 @@ FESeries::process_coefficients( // first, parse all table elements into a map of predicate values and // coefficients. We could have stored (predicate values ->TableIndicies) map, // but its processing would have been much harder later on. - std::map> pred_to_values; + std::map> pred_to_values; internal::FESeriesImplementation::fill_map(coefficients, predicate, pred_to_values); // now go through the map and populate the @p norm_values based on @p norm: - for (typename std::map>::const_iterator it = - pred_to_values.begin(); - it != pred_to_values.end(); - ++it) + for (const auto &pred_to_value : pred_to_values) { - predicate_values.push_back(it->first); - Vector values(it->second.begin(), it->second.end()); + predicate_values.push_back(pred_to_value.first); + Vector values(pred_to_value.second.cbegin(), + pred_to_value.second.cend()); switch (norm) { diff --git a/source/fe/fe_series_fourier.cc b/source/fe/fe_series_fourier.cc index 7763754986..260a221aa2 100644 --- a/source/fe/fe_series_fourier.cc +++ b/source/fe/fe_series_fourier.cc @@ -182,11 +182,12 @@ namespace FESeries template + template void Fourier::calculate( - const Vector & local_dof_values, - const unsigned int cell_active_fe_index, - Table> &fourier_coefficients) + const Vector & local_dof_values, + const unsigned int cell_active_fe_index, + Table &fourier_coefficients) { ensure_existence(*fe_collection, *q_collection, @@ -194,12 +195,12 @@ namespace FESeries cell_active_fe_index, fourier_transform_matrices); - const FullMatrix> &matrix = + const FullMatrix &matrix = fourier_transform_matrices[cell_active_fe_index]; std::fill(unrolled_coefficients.begin(), unrolled_coefficients.end(), - std::complex(0.)); + CoefficientType(0.)); Assert(unrolled_coefficients.size() == matrix.m(), ExcInternalError()); diff --git a/source/fe/fe_series_fourier.inst.in b/source/fe/fe_series_fourier.inst.in index d07d7c1b96..379c956e98 100644 --- a/source/fe/fe_series_fourier.inst.in +++ b/source/fe/fe_series_fourier.inst.in @@ -22,3 +22,17 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS) deal_II_space_dimension>; #endif } + +for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS; + SCALAR : REAL_SCALARS) + { +#if deal_II_dimension <= deal_II_space_dimension + template void + FESeries::Fourier::calculate( + const Vector &, + const unsigned int, + Table::CoefficientType> &); +#endif + } diff --git a/source/fe/fe_series_legendre.cc b/source/fe/fe_series_legendre.cc index 566c104863..a9171989ca 100644 --- a/source/fe/fe_series_legendre.cc +++ b/source/fe/fe_series_legendre.cc @@ -202,11 +202,12 @@ namespace FESeries template + template void Legendre::calculate( - const dealii::Vector &local_dof_values, + const dealii::Vector &local_dof_values, const unsigned int cell_active_fe_index, - Table & legendre_coefficients) + Table & legendre_coefficients) { ensure_existence(*fe_collection, *q_collection, @@ -214,10 +215,12 @@ namespace FESeries cell_active_fe_index, legendre_transform_matrices); - const FullMatrix &matrix = + const FullMatrix &matrix = legendre_transform_matrices[cell_active_fe_index]; - std::fill(unrolled_coefficients.begin(), unrolled_coefficients.end(), 0.); + std::fill(unrolled_coefficients.begin(), + unrolled_coefficients.end(), + CoefficientType(0.)); Assert(unrolled_coefficients.size() == matrix.m(), ExcInternalError()); diff --git a/source/fe/fe_series_legendre.inst.in b/source/fe/fe_series_legendre.inst.in index 587120ed8d..ba54311eb4 100644 --- a/source/fe/fe_series_legendre.inst.in +++ b/source/fe/fe_series_legendre.inst.in @@ -22,3 +22,17 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS) deal_II_space_dimension>; #endif } + +for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS; + SCALAR : REAL_SCALARS) + { +#if deal_II_dimension <= deal_II_space_dimension + template void + FESeries::Legendre::calculate( + const Vector &, + const unsigned int, + Table::CoefficientType> &); +#endif + }