]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Proper output file.
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 11 Jan 2007 03:52:46 +0000 (03:52 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 11 Jan 2007 03:52:46 +0000 (03:52 +0000)
git-svn-id: https://svn.dealii.org/trunk@14303 0785d39b-7218-0410-832d-ea1e28bc413d

tests/hp/step-11.cc [new file with mode: 0644]
tests/hp/step-11/cmp/generic [new file with mode: 0644]

diff --git a/tests/hp/step-11.cc b/tests/hp/step-11.cc
new file mode 100644 (file)
index 0000000..258ea26
--- /dev/null
@@ -0,0 +1,251 @@
+//----------------------------  step-11.cc  ---------------------------
+//    $Id$
+//    Version: $Name$ 
+//
+//    Copyright (C) 2005, 2006, 2007 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//----------------------------  step-11.cc  ---------------------------
+
+
+// a hp-ified version of step-11
+
+
+#include "../tests.h"
+
+#include <base/logstream.h>
+#include <fstream>
+std::ofstream logfile("step-11/output");
+
+#include <base/quadrature_lib.h>
+#include <base/function.h>
+#include <base/logstream.h>
+#include <base/table_handler.h>
+#include <lac/vector.h>
+#include <lac/sparse_matrix.h>
+#include <lac/solver_cg.h>
+#include <lac/precondition.h>
+#include <grid/tria.h>
+#include <grid/grid_generator.h>
+#include <grid/tria_boundary_lib.h>
+#include <grid/tria_accessor.h>
+#include <grid/tria_iterator.h>
+#include <dofs/dof_handler.h>
+#include <dofs/dof_constraints.h>
+#include <dofs/dof_accessor.h>
+#include <dofs/dof_tools.h>
+#include <fe/fe_q.h>
+#include <fe/hp_fe_values.h>
+#include <fe/mapping_q.h>
+#include <numerics/vectors.h>
+#include <numerics/matrices.h>
+
+#include <lac/compressed_sparsity_pattern.h>
+
+#include <algorithm>
+#include <iostream>
+#include <iomanip>
+#include <cmath>
+
+
+template <int dim>
+class LaplaceProblem 
+{
+  public:
+    LaplaceProblem (const unsigned int mapping_degree);
+    void run ();
+    
+  private:
+    void setup_system ();
+    void assemble_and_solve ();
+    void solve ();
+
+    Triangulation<dim>   triangulation;
+    hp::FECollection<dim>            fe;
+    hp::DoFHandler<dim>      dof_handler;
+    hp::MappingCollection<dim>        mapping;
+
+    SparsityPattern      sparsity_pattern;
+    SparseMatrix<double> system_matrix;
+    ConstraintMatrix     mean_value_constraints;
+
+    Vector<double>       solution;
+    Vector<double>       system_rhs;
+
+    TableHandler         output_table;
+};
+
+
+
+template <int dim>
+LaplaceProblem<dim>::LaplaceProblem (const unsigned int mapping_degree) :
+                fe (FE_Q<dim>(1)),
+               dof_handler (triangulation),
+               mapping (MappingQ<dim>(mapping_degree))
+{
+  deallog << "Using mapping with degree " << mapping_degree << ":"
+           << std::endl
+           << "============================"
+           << std::endl;
+}
+
+
+
+template <int dim>
+void LaplaceProblem<dim>::setup_system ()
+{
+  dof_handler.distribute_dofs (fe);
+  solution.reinit (dof_handler.n_dofs());
+  system_rhs.reinit (dof_handler.n_dofs());
+
+  std::vector<bool> boundary_dofs (dof_handler.n_dofs(), false);
+  DoFTools::extract_boundary_dofs (dof_handler, std::vector<bool>(1,true),
+                                  boundary_dofs);
+
+  const unsigned int first_boundary_dof
+    = std::distance (boundary_dofs.begin(),
+                    std::find (boundary_dofs.begin(),
+                               boundary_dofs.end(),
+                               true));
+
+  mean_value_constraints.clear ();
+  mean_value_constraints.add_line (first_boundary_dof);
+  for (unsigned int i=first_boundary_dof+1; i<dof_handler.n_dofs(); ++i)
+    if (boundary_dofs[i] == true)
+      mean_value_constraints.add_entry (first_boundary_dof,
+                                       i, -1);
+  mean_value_constraints.close ();
+
+  CompressedSparsityPattern csp (dof_handler.n_dofs(),
+                                dof_handler.n_dofs());
+  DoFTools::make_sparsity_pattern (dof_handler, csp);
+  mean_value_constraints.condense (csp);
+
+  sparsity_pattern.copy_from (csp);
+  system_matrix.reinit (sparsity_pattern);
+}
+
+
+
+template <int dim>
+void LaplaceProblem<dim>::assemble_and_solve () 
+{
+
+  const unsigned int gauss_degree
+    = std::max (static_cast<unsigned int>(std::ceil(1.*(static_cast<const MappingQ<dim>&>(mapping[0]).get_degree()+1)/2)),
+               2U);
+  MatrixTools::create_laplace_matrix (mapping, dof_handler,
+                                     hp::QCollection<dim>(QGauss<dim>(gauss_degree)),
+                                     system_matrix);
+  VectorTools::create_right_hand_side (mapping, dof_handler,
+                                      hp::QCollection<dim>(QGauss<dim>(gauss_degree)),
+                                      ConstantFunction<dim>(-2),
+                                      system_rhs);
+  Vector<double> tmp (system_rhs.size());
+  VectorTools::create_boundary_right_hand_side (mapping, dof_handler,
+                                               hp::QCollection<dim-1>(QGauss<dim-1>(gauss_degree)),
+                                               ConstantFunction<dim>(1),
+                                               tmp);
+  system_rhs += tmp;
+
+  mean_value_constraints.condense (system_matrix);
+  mean_value_constraints.condense (system_rhs);  
+
+  solve ();
+  mean_value_constraints.distribute (solution);
+
+  Vector<float> norm_per_cell (triangulation.n_active_cells());
+  VectorTools::integrate_difference (mapping, dof_handler,
+                                    solution,
+                                    ZeroFunction<dim>(),
+                                    norm_per_cell,
+                                    hp::QCollection<dim>(QGauss<dim>(gauss_degree+1)),
+                                    VectorTools::H1_seminorm);
+  const double norm = norm_per_cell.l2_norm();
+
+  output_table.add_value ("cells", triangulation.n_active_cells());
+  output_table.add_value ("|u|_1", norm);
+  output_table.add_value ("error", std::fabs(norm-std::sqrt(3.14159265358/2)));
+}
+
+
+
+template <int dim>
+void LaplaceProblem<dim>::solve () 
+{
+  SolverControl           solver_control (1000, 1e-12);
+  SolverCG<>              cg (solver_control);
+
+  PreconditionSSOR<> preconditioner;
+  preconditioner.initialize(system_matrix, 1.2);
+
+  cg.solve (system_matrix, solution, system_rhs,
+           preconditioner);
+}
+
+
+
+template <int dim>
+void LaplaceProblem<dim>::run () 
+{
+  GridGenerator::hyper_ball (triangulation);
+  static const HyperBallBoundary<dim> boundary;
+  triangulation.set_boundary (0, boundary);
+  
+  for (unsigned int cycle=0; cycle<6; ++cycle, triangulation.refine_global(1))
+    {
+      setup_system ();
+      assemble_and_solve ();
+    };
+
+  output_table.set_precision("|u|_1", 6);
+  output_table.set_precision("error", 6);
+  output_table.write_text (deallog.get_file_stream());
+  deallog << std::endl;
+}
+
+    
+
+int main () 
+{
+  try
+    {
+      logfile.precision(2);
+  
+      deallog.attach(logfile);
+      deallog.depth_console(0);
+      deallog.threshold_double(1.e-10);  
+
+      for (unsigned int mapping_degree=1; mapping_degree<=3; ++mapping_degree)
+       LaplaceProblem<2>(mapping_degree).run ();
+    }
+  catch (std::exception &exc)
+    {
+      std::cerr << std::endl << std::endl
+               << "----------------------------------------------------"
+               << std::endl;
+      std::cerr << "Exception on processing: " << std::endl
+               << exc.what() << std::endl
+               << "Aborting!" << std::endl
+               << "----------------------------------------------------"
+               << std::endl;
+      return 1;
+    }
+  catch (...) 
+    {
+      std::cerr << std::endl << std::endl
+               << "----------------------------------------------------"
+               << std::endl;
+      std::cerr << "Unknown exception!" << std::endl
+               << "Aborting!" << std::endl
+               << "----------------------------------------------------"
+               << std::endl;
+      return 1;
+    };
+
+  return 0;
+}
diff --git a/tests/hp/step-11/cmp/generic b/tests/hp/step-11/cmp/generic
new file mode 100644 (file)
index 0000000..4a28f6c
--- /dev/null
@@ -0,0 +1,67 @@
+
+DEAL::Using mapping with degree 1:
+DEAL::============================
+DEAL:cg::Starting value 1.1
+DEAL:cg::Convergence step 7 value 0
+DEAL:cg::Starting value 0.98
+DEAL:cg::Convergence step 18 value 0
+DEAL:cg::Starting value 0.65
+DEAL:cg::Convergence step 29 value 0
+DEAL:cg::Starting value 0.37
+DEAL:cg::Convergence step 53 value 0
+DEAL:cg::Starting value 0.19
+DEAL:cg::Convergence step 108 value 0
+DEAL:cg::Starting value 0.10
+DEAL:cg::Convergence step 238 value 0
+cells  |u|_1    error   
+5     0.680402 0.572912 
+20    1.085518 0.167796 
+80    1.208981 0.044334 
+320   1.242041 0.011273 
+1280  1.250482 0.002832 
+5120  1.252605 0.000709 
+DEAL::
+DEAL::Using mapping with degree 2:
+DEAL::============================
+DEAL:cg::Starting value 1.338917
+DEAL:cg::Convergence step 7 value 0
+DEAL:cg::Starting value 1.028493
+DEAL:cg::Convergence step 17 value 0
+DEAL:cg::Starting value 0.654097
+DEAL:cg::Convergence step 29 value 0
+DEAL:cg::Starting value 0.368099
+DEAL:cg::Convergence step 53 value 0
+DEAL:cg::Starting value 0.194917
+DEAL:cg::Convergence step 108 value 0
+DEAL:cg::Starting value 0.100213
+DEAL:cg::Convergence step 238 value 0
+cells  |u|_1    error   
+5     1.050963 0.202351 
+20    1.199642 0.053672 
+80    1.239913 0.013401 
+320   1.249987 0.003327 
+1280  1.252486 0.000828 
+5120  1.253108 0.000206 
+DEAL::
+DEAL::Using mapping with degree 3:
+DEAL::============================
+DEAL:cg::Starting value 1.358459
+DEAL:cg::Convergence step 7 value 0
+DEAL:cg::Starting value 1.030204
+DEAL:cg::Convergence step 17 value 0
+DEAL:cg::Starting value 0.654234
+DEAL:cg::Convergence step 29 value 0
+DEAL:cg::Starting value 0.368109
+DEAL:cg::Convergence step 53 value 0
+DEAL:cg::Starting value 0.194918
+DEAL:cg::Convergence step 108 value 0
+DEAL:cg::Starting value 0.100213
+DEAL:cg::Convergence step 238 value 0
+cells  |u|_1    error   
+5     1.086161 0.167153 
+20    1.204349 0.048965 
+80    1.240502 0.012812 
+320   1.250059 0.003255 
+1280  1.252495 0.000819 
+5120  1.253109 0.000205 
+DEAL::

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.