]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Fixed a couple of bugs and completed the implementation of the convergence table.
authorheltai <heltai@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 16 Apr 2008 07:56:21 +0000 (07:56 +0000)
committerheltai <heltai@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 16 Apr 2008 07:56:21 +0000 (07:56 +0000)
git-svn-id: https://svn.dealii.org/trunk@15999 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/base/include/base/convergence_table.h
deal.II/base/source/convergence_table.cc
deal.II/base/source/utilities.cc
deal.II/deal.II/include/numerics/fe_field_function.h
deal.II/doc/news/changes.h

index 35202a88e657f0f6db58d30adf89e094d5d533df..401ff0f1d1ac277274a15df605a8724499b88b9f 100644 (file)
@@ -104,7 +104,25 @@ class ConvergenceTable: public TableHandler
                                      * 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,
index c165e73f1fa2169b095ecc8bfc5c263849de2646..445320028c1a748636d1e34de0431c2d5dc2ce52 100644 (file)
@@ -34,11 +34,16 @@ void ConvergenceTable::evaluate_convergence_rates(const std::string &data_column
     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)
@@ -51,6 +56,18 @@ void ConvergenceTable::evaluate_convergence_rates(const std::string &data_column
        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)
@@ -58,7 +75,25 @@ void ConvergenceTable::evaluate_convergence_rates(const std::string &data_column
       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());  
     }
index be1ae32579b1acfbe60becb84605f29ae54fa7bb..d6aef3d73396f1a2580aa38891d67d1df50aef9b 100644 (file)
@@ -66,14 +66,19 @@ namespace Utilities
                     (digits==2 && i>=100)  ||
                     (digits==3 && i>=1000) ||
                     (digits==4 && i>=10000)||
-                    (i>=100000)),
+                    (digits==5 && i>=100000)||
+                    (i>=1000000)),
                  ExcInvalidNumber2StringConversersion(i, digits));
   
     std::string s;
     switch (digits) 
       {
+       case 6:
+             s += '0' + i/100000;
+       case 5:
+             s += '0' + (i%100000)/10000;
        case 4:
-             s += '0' + i/1000;
+             s += '0' + (i%10000)/1000;
        case 3:
              s += '0' + (i%1000)/100;
        case 2:
index 24bb7e48a340940623991f9e7756e054bfb285e7..54786fe11e37d564b4582f222eab48fc47b6fbf7 100644 (file)
@@ -153,7 +153,7 @@ namespace Functions
                                        * function. This will speed
                                        * things up a little.
                                        */
-      inline void set_active_cell(typename DH::active_cell_iterator &newcell);
+      void set_active_cell(typename DH::active_cell_iterator &newcell);
   
                                       /**
                                        * Get ONE vector value at the
index 3df64c89438684ca1cc14b195f1d0655474b46aa..2b0ce30567053f2abee2b69786e2cb15d1ffa121 100644 (file)
@@ -139,6 +139,21 @@ inconvenience this causes.
 
 <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>
@@ -290,6 +305,12 @@ constraints individually.
 
 <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)

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.