* the <tt>set_tex_supercaption (...)</tt> function
* of the base class TableHandler.
*
- * Still not implemented.
+ * This method behaves in the following way:
+ *
+ * If RateMode is reduction_rate, then the computed
+ * output is
+ * \f[ \frac{e_{n-1}/k_{n-1}}{e_n/k_n} \f].
+ *
+ * Where \f$k\f$ is the reference column.
+ *
+ * If RateMode is reduction_rate_log2, then the
+ * computed output is:
+ * \f[
+ * 2\frac{\log |e_{n-1}/e_{n}|}{\log |k_n/k_{n-1}|}
+ * \f].
+ *
+ * This is useful, for example, if we use as
+ * reference key the number of degrees of freedom.
+ * Assuming that the error is proportional to
+ * \f$ C (1/\sqrt{k})^r \f$, then this method will
+ * produce the rate \f$ r \f$ as a result.
*/
void
evaluate_convergence_rates (const std::string &data_column_key,
return;
std::vector<TableEntryBase *> &entries=columns[data_column_key].entries;
+ std::vector<TableEntryBase *> &ref_entries=columns[reference_column_key].entries;
std::string rate_key=data_column_key;
const unsigned int n=entries.size();
+ const unsigned int n_ref=ref_entries.size();
+ Assert(n == n_ref, ExcDimensionMismatch(n, n_ref));
std::vector<double> values(n);
+ std::vector<double> ref_values(n_ref);
+
for (unsigned int i=0; i<n; ++i)
{
if (dynamic_cast<TableEntry<double>*>(entries[i]) != 0)
values[i]=dynamic_cast<TableEntry<unsigned int>*>(entries[i])->value();
else
Assert(false, ExcWrongValueType());
+
+ // And now the reference values.
+ if (dynamic_cast<TableEntry<double>*>(ref_entries[i]) != 0)
+ ref_values[i]=dynamic_cast<TableEntry<double>*>(ref_entries[i])->value();
+ else if (dynamic_cast<TableEntry<float>*>(ref_entries[i]) != 0)
+ ref_values[i]=dynamic_cast<TableEntry<float>*>(ref_entries[i])->value();
+ else if (dynamic_cast<TableEntry<int>*>(ref_entries[i]) != 0)
+ ref_values[i]=dynamic_cast<TableEntry<int>*>(ref_entries[i])->value();
+ else if (dynamic_cast<TableEntry<unsigned int>*>(ref_entries[i]) != 0)
+ ref_values[i]=dynamic_cast<TableEntry<unsigned int>*>(ref_entries[i])->value();
+ else
+ Assert(false, ExcWrongValueType());
}
switch (rate_mode)
case none:
break;
case reduction_rate:
+ rate_key+="red.rate";
+ Assert(columns.count(rate_key)==0, ExcRateColumnAlreadyExists(rate_key));
+ // no value available for the
+ // first row
+ add_value(rate_key, std::string("-"));
+ for (unsigned int i=1; i<n; ++i)
+ add_value(rate_key, values[i-1]/values[i] *
+ ref_values[i]/ref_values[i-1]);
+ break;
case reduction_rate_log2:
+ rate_key+="red.rate";
+ Assert(columns.count(rate_key)==0, ExcRateColumnAlreadyExists(rate_key));
+ // no value available for the
+ // first row
+ add_value(rate_key, std::string("-"));
+ for (unsigned int i=1; i<n; ++i)
+ add_value(rate_key, 2*std::log(std::fabs(values[i-1]/values[i])) /
+ std::log(std::fabs(ref_values[i]/ref_values[i-1])));
+ break;
default:
Assert(false, ExcNotImplemented());
}
<ol>
+<li> <p>New: ConvergenceTable.evaluate_convergence_rates has been incomplete
+for a long time. I added an implementation for reduction_rate and
+reduction_rate_lo2 in the case where one specifies a reference column
+that allows to compute error rates also for local refinement.
+<br>
+(Luca Heltai 2008/04/16)
+</p></li>
+
+<li> <p> Fixed: The Utilities::int_to_string function did not work for 5 or 6
+digits, while its counterpart Utilities::needed_digits does. This is now fixed,
+and the two are coherent.
+<br>
+(Luca Heltai 2008/04/16)
+</p> </li>
+
<li> <p> New: The SymmetricTensor class now has a constructor that creates
an object from an array consisting its independent components.
<br>
<ol>
+ <li> <p>Fixed: Functions::FEFieldFunction<dim>.set_active_cell was declared
+ inline in the declaration but not in the implementation. This is now fixed.
+ <br>
+ (Luca Heltai 2008/04/16)
+ </p></li>
+
<li> <p>New: The FE_RaviartThomasNodal now supports hp functionality.
<br>
(Zhu Liang, 2008/04/04)