]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
slit domain and add postprocessing
authorkanschat <kanschat@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 7 Apr 2010 16:15:24 +0000 (16:15 +0000)
committerkanschat <kanschat@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 7 Apr 2010 16:15:24 +0000 (16:15 +0000)
git-svn-id: https://svn.dealii.org/trunk@20961 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/examples/step-39/Makefile
deal.II/examples/step-39/postprocess.pl [new file with mode: 0644]
deal.II/examples/step-39/step-39.cc

index db1de210668737e0880ac6fef22054e307a823d6..4a33ab61d87d26dc3e7bc32f34c1bf116e03de96 100644 (file)
@@ -97,7 +97,8 @@ $(target) : $(libraries)
 run: $(target)
        @echo ============================ Running $<
        @./$(target)$(EXEEXT)
-
+       @echo ============================ Postprocessing to step-39.dat
+       @$(PERL) postprocess.pl step-39.log > step-39.dat
 
 # As a last rule to the `make' program, we define what to do when
 # cleaning up a directory. This usually involves deleting object files
diff --git a/deal.II/examples/step-39/postprocess.pl b/deal.II/examples/step-39/postprocess.pl
new file mode 100644 (file)
index 0000000..c830dc6
--- /dev/null
@@ -0,0 +1,25 @@
+######################################################################
+# $Id$
+######################################################################
+# Postprocess logstream output and create a data file for gnuplot
+######################################################################
+
+use strict;
+
+my $step;     # The iteration step in the adaptive loop
+my @dofs;     # The number of degrees of freedom in each step
+my @error;    # The error of the solution
+my @estimate; # The a posteriori error estimate
+
+while(<>)
+{
+    $step = $1 if m/DEAL::Step\s*(\d+)/;
+    $dofs[$step] = $1 if m/DEAL::DoFHandler\s*(\d+)/;
+    $error[$step] = $1 if m/DEAL::Error\s*(\S+)/;
+    $estimate[$step] = $1 if m/DEAL::Estimate\s*(\S+)/;
+}
+
+for (my $i=0;$i<=$step;++$i)
+{
+    printf "%-3d\t%-7d\t%g\t%g\n", $i, $dofs[$i], $error[$i], $estimate[$i];
+}
index 1c86d1e32de7ba213e326d729ec23398a172a367..04dfa1bc484f6739d8af7be7da104238a9b5032c 100644 (file)
@@ -69,7 +69,7 @@ using namespace dealii;
                                 // This is the function we use to set
                                 // the boundary values and also the
                                 // exact solution we compare to.
-Functions::LSingularityFunction exact_solution;
+Functions::SlitSingularityFunction<2> exact_solution;
 
                                 // @sect3{The local integrators}
 
@@ -273,6 +273,7 @@ void Estimator<dim>::cell(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::
       const double t = dinfo.cell->diameter() * trace(DDuh[k]);
       dinfo.value(0) +=  t*t * fe.JxW(k);
     }
+  dinfo.value(0) = std::sqrt(dinfo.value(0));
 }
 
 
@@ -291,7 +292,8 @@ void Estimator<dim>::bdry(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::
   
   for (unsigned k=0;k<fe.n_quadrature_points;++k)
     dinfo.value(0) += penalty * (boundary_values[k] - uh[k]) * (boundary_values[k] - uh[k])
-                * fe.JxW(k);
+                     * fe.JxW(k);
+  dinfo.value(0) = std::sqrt(dinfo.value(0));  
 }
 
 
@@ -320,6 +322,7 @@ void Estimator<dim>::face(MeshWorker::DoFInfo<dim>& dinfo1,
       dinfo1.value(0) += (penalty * diff1*diff1 + h * diff2*diff2)
                    * fe.JxW(k);
     }
+  dinfo1.value(0) = std::sqrt(dinfo1.value(0));
   dinfo2.value(0) = dinfo1.value(0);
 }
 
@@ -375,7 +378,7 @@ Step39<dim>::Step39(const FiniteElement<dim>& fe)
                dof_handler(mg_dof_handler),
                estimates(1)
 {
-  GridGenerator::hyper_L(triangulation, -1, 1);
+  GridGenerator::hyper_cube_slit(triangulation, -1, 1);
 }
 
 
@@ -582,7 +585,7 @@ Step39<dim>::error()
   
   QGauss<dim> quadrature(n_gauss_points);
   VectorTools::integrate_difference(mapping, dof_handler, solution, exact_solution,
-                                   cell_errors, quadrature, VectorTools::L2_norm);
+                                   cell_errors, quadrature, VectorTools::H1_norm);
   deallog << "Error    " << cell_errors.l2_norm() << std::endl;
 }
 
@@ -645,8 +648,8 @@ void Step39<dim>::output_results (const unsigned int cycle) const
 
   std::string filename(fn);
   filename += ".gnuplot";
-  std::cout << "Writing solution to <" << filename << ">..."
-           << std::endl << std::endl;
+  deallog << "Writing solution to <" << filename << ">..."
+         << std::endl << std::endl;
   std::ofstream gnuplot_output (filename.c_str());
   
   DataOut<dim> data_out;
@@ -708,7 +711,8 @@ Step39<dim>::run(unsigned int n_steps)
 
 int main()
 {
-  deallog.log_execution_time(true);
+  std::ofstream logfile("step-39.log");
+  deallog.attach(logfile);
   FE_DGQ<2> fe1(3);
   Step39<2> test1(fe1);
   test1.run(20);

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.