From: Luca Heltai Date: Wed, 1 May 2019 13:11:52 +0000 (+0200) Subject: exclude_from_rates -> compute_rate X-Git-Tag: v9.1.0-rc1~143^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=167e9c83862077bf810f33fd9340f56d934786dc;p=dealii.git exclude_from_rates -> compute_rate --- diff --git a/include/deal.II/base/parsed_convergence_table.h b/include/deal.II/base/parsed_convergence_table.h index 2132019ab6..32b14fe04a 100644 --- a/include/deal.II/base/parsed_convergence_table.h +++ b/include/deal.II/base/parsed_convergence_table.h @@ -222,7 +222,7 @@ public: * # Number of digits to use when printing the error. * set Error precision = 3 * - * # Extra columns to add to the table. Available options are dofs and cells + * # Extra columns to add to the table. Available options are dofs and cells. * set Extra columns = dofs, cells * * # The exponent to use when computing p-norms. @@ -306,7 +306,7 @@ public: * error_from_exact() or difference(). * * Make sure you add all extra columns before the first call to - * error_from_exact() or difference(), since adding additional columns to the + * error_from_exact() or difference(). Adding additional columns to the * convergence table after you already started filling the table will trigger * an exception. * @@ -321,7 +321,7 @@ public: * return dt; * }; * - * table.add_extra_column("dt", dt_function, true); + * table.add_extra_column("dt", dt_function, false); * * for (unsigned int i = 0; i < n_cycles; ++i) * { @@ -343,9 +343,9 @@ public: * provided that you use the following parameter file (only non default * entries are shown here): * @code - * set Extra columns = - * set List of error norms to compute = L2_norm - * set Rate key = dt + * set Extra columns = + * set List of error norms to compute = L2_norm + * set Rate key = dt * @endcode * * @@ -353,13 +353,16 @@ public: * @param custom_function Function that will be called to fill the given * entry. You need to make sure that the scope of this function is valid * up to the call to error_from_exact() or difference(); - * @param exclude_from_rates If set to true, then no error rates will be - * computed for this column. + * @param compute_rate If set to true, then this column will be included in + * the list of columns for which error rates are computed. You may want to set + * this to false if you want to compute error rates with respect to this + * column. In this case, you should also specify @p column_name as the rate + * key in the parameter file. */ void add_extra_column(const std::string & column_name, const std::function &custom_function, - const bool & exclude_from_rates = false); + const bool & compute_rate = true); /** * Difference between two solutions in the same vector space. @@ -434,7 +437,7 @@ private: std::vector> norms_per_unique_component; /** - * Exponent to use in the p-norm types. + * Exponent to use in p-norm types. */ double exponent; @@ -469,8 +472,8 @@ private: std::string error_file_name; /** - * Compute the error. If this is false, all functions regarding - * errors are disabled and don't do anything. + * Compute the error. If this is false, all methods that perform the + * computation of the error are disabled and don't do anything. */ bool compute_error; }; diff --git a/source/base/parsed_convergence_table.cc b/source/base/parsed_convergence_table.cc index c4e1445788..72761b3065 100644 --- a/source/base/parsed_convergence_table.cc +++ b/source/base/parsed_convergence_table.cc @@ -171,12 +171,11 @@ ParsedConvergenceTable::prepare_table_for_output() has_key = true; if (col != "") - table.omit_column_from_convergence_rate_evaluation( - Patterns::Tools::to_string(col)); + table.omit_column_from_convergence_rate_evaluation(col); } for (const auto &extra_col : extra_column_functions) - if (extra_col.second.second) + if (extra_col.second.second == false) { if (rate_key == extra_col.first) has_key = true; @@ -188,12 +187,10 @@ ParsedConvergenceTable::prepare_table_for_output() { if (rate_mode == "reduction_rate_log2") table.evaluate_all_convergence_rates( - Patterns::Tools::to_string(rate_key), - ConvergenceTable::reduction_rate_log2); + rate_key, ConvergenceTable::reduction_rate_log2); else if (rate_mode == "reduction_rate") table.evaluate_all_convergence_rates( - Patterns::Tools::to_string(rate_key), - ConvergenceTable::reduction_rate); + rate_key, ConvergenceTable::reduction_rate); else { Assert(rate_mode != "none", ExcInternalError()); @@ -264,9 +261,9 @@ void ParsedConvergenceTable::add_extra_column( const std::string & column_name, const std::function &custom_function, - const bool & exclude_from_rates) + const bool & compute_rate) { - extra_column_functions[column_name] = {custom_function, exclude_from_rates}; + extra_column_functions[column_name] = {custom_function, compute_rate}; } DEAL_II_NAMESPACE_CLOSE diff --git a/tests/base/parsed_convergence_table_10.cc b/tests/base/parsed_convergence_table_10.cc index 834f995923..d217a138df 100644 --- a/tests/base/parsed_convergence_table_10.cc +++ b/tests/base/parsed_convergence_table_10.cc @@ -66,7 +66,7 @@ main() auto dt = [&]() { return i + 1.0; }; table.add_extra_column("cycle", cycle); - table.add_extra_column("dt", dt, true); + table.add_extra_column("dt", dt, false); table.error_from_exact(dh, sol, exact); }