From 38325cd7b7c65fdc678dd4c15d7b06060d0950c6 Mon Sep 17 00:00:00 2001 From: Peter Munch Date: Tue, 23 Jun 2020 22:50:44 +0200 Subject: [PATCH] Remove DoFHandlerType from base/parsed_convergence_table.h --- .../deal.II/base/parsed_convergence_table.h | 122 ++++++++---------- 1 file changed, 53 insertions(+), 69 deletions(-) diff --git a/include/deal.II/base/parsed_convergence_table.h b/include/deal.II/base/parsed_convergence_table.h index 9851303783..5681e237c9 100644 --- a/include/deal.II/base/parsed_convergence_table.h +++ b/include/deal.II/base/parsed_convergence_table.h @@ -270,26 +270,23 @@ public: * components does not coincide with the number of components of the * underlying finite element space. */ - template + template void - error_from_exact( - const DoFHandlerType & vspace, - const VectorType & solution, - const Function &exact, - const Function *weight = nullptr); + error_from_exact(const DoFHandler &vspace, + const VectorType & solution, + const Function & exact, + const Function * weight = nullptr); /** * Same as above, with a different mapping. */ - template + template void - error_from_exact( - const Mapping - & mapping, - const DoFHandlerType & vspace, - const VectorType & solution, - const Function &exact, - const Function *weight = nullptr); + error_from_exact(const Mapping & mapping, + const DoFHandler &vspace, + const VectorType & solution, + const Function & exact, + const Function * weight = nullptr); /** * Add an additional column (with name @p column_name) to the table, by invoking @@ -365,24 +362,23 @@ public: /** * Difference between two solutions in the same vector space. */ - template + template void - difference(const DoFHandlerType &, + difference(const DoFHandler &, const VectorType &, const VectorType &, - const Function *weight = nullptr); + const Function *weight = nullptr); /** * Same as above, with a non default mapping. */ - template + template void - difference(const Mapping &mapping, - const DoFHandlerType &, + difference(const Mapping &mapping, + const DoFHandler &, const VectorType &, const VectorType &, - const Function *weight = nullptr); + const Function *weight = nullptr); /** * Write the error table to the @p out stream (in text format), and @@ -482,81 +478,69 @@ private: // ============================================================ // Template functions // ============================================================ -template +template void -ParsedConvergenceTable::difference( - const DoFHandlerType & dh, - const VectorType & solution1, - const VectorType & solution2, - const Function *weight) +ParsedConvergenceTable::difference(const DoFHandler &dh, + const VectorType & solution1, + const VectorType & solution2, + const Function * weight) { AssertThrow(solution1.size() == solution2.size(), ExcDimensionMismatch(solution1.size(), solution2.size())); VectorType solution(solution1); solution -= solution2; - error_from_exact(dh, - solution, - Functions::ConstantFunction( - 0, component_names.size()), - weight); + error_from_exact( + dh, + solution, + Functions::ConstantFunction(0, component_names.size()), + weight); } -template +template void -ParsedConvergenceTable::difference( - const Mapping - & mapping, - const DoFHandlerType & dh, - const VectorType & solution1, - const VectorType & solution2, - const Function *weight) +ParsedConvergenceTable::difference(const Mapping & mapping, + const DoFHandler &dh, + const VectorType & solution1, + const VectorType & solution2, + const Function * weight) { AssertThrow(solution1.size() == solution2.size(), ExcDimensionMismatch(solution1.size(), solution2.size())); VectorType solution(solution1); solution -= solution2; - error_from_exact(mapping, - dh, - solution, - Functions::ConstantFunction( - 0, component_names.size()), - weight); + error_from_exact( + mapping, + dh, + solution, + Functions::ConstantFunction(0, component_names.size()), + weight); } -template +template void -ParsedConvergenceTable::error_from_exact( - const DoFHandlerType & dh, - const VectorType & solution, - const Function &exact, - const Function *weight) +ParsedConvergenceTable::error_from_exact(const DoFHandler &dh, + const VectorType & solution, + const Function &exact, + const Function *weight) { - error_from_exact(StaticMappingQ1::mapping, - dh, - solution, - exact, - weight); + error_from_exact( + StaticMappingQ1::mapping, dh, solution, exact, weight); } -template +template void -ParsedConvergenceTable::error_from_exact( - const Mapping - & mapping, - const DoFHandlerType & dh, - const VectorType & solution, - const Function &exact, - const Function *weight) +ParsedConvergenceTable::error_from_exact(const Mapping &mapping, + const DoFHandler &dh, + const VectorType & solution, + const Function &exact, + const Function *weight) { - const int dim = DoFHandlerType::dimension; - const int spacedim = DoFHandlerType::space_dimension; const auto n_components = component_names.size(); if (compute_error) -- 2.39.5