]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Implement VectorTools::compute_mean_value for dim!=spacedim.
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 12 May 2011 03:10:56 +0000 (03:10 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 12 May 2011 03:10:56 +0000 (03:10 +0000)
git-svn-id: https://svn.dealii.org/trunk@23691 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/numerics/vectors.templates.h
deal.II/source/numerics/vectors.inst.in
tests/codim_one/mean_value.cc [new file with mode: 0644]
tests/codim_one/mean_value/cmp/generic [new file with mode: 0644]

index 6e9293591d57625dee18175497d2c1ff8a477aec..ad435b58ebd37a609881a8d6df89c4eb8a14e076 100644 (file)
@@ -5243,9 +5243,9 @@ VectorTools::compute_mean_value (const Mapping<dim, spacedim>    &mapping,
   Assert (component < dof.get_fe().n_components(),
          ExcIndexRange(component, 0, dof.get_fe().n_components()));
 
-  FEValues<dim> fe(mapping, dof.get_fe(), quadrature,
-                  UpdateFlags(update_JxW_values
-                              | update_values));
+  FEValues<dim,spacedim> fe(mapping, dof.get_fe(), quadrature,
+                           UpdateFlags(update_JxW_values
+                                       | update_values));
 
   typename DoFHandler<dim,spacedim>::active_cell_iterator cell;
   std::vector<Vector<double> > values(quadrature.size(),
@@ -5298,7 +5298,7 @@ VectorTools::compute_mean_value (const DoFHandler<dim,spacedim> &dof,
                                 const unsigned int     component)
 {
   Assert (DEAL_II_COMPAT_MAPPING, ExcCompatibility("mapping"));
-  return compute_mean_value(StaticMappingQ1<dim>::mapping, dof, quadrature, v, component);
+  return compute_mean_value(StaticMappingQ1<dim,spacedim>::mapping, dof, quadrature, v, component);
 }
 
 DEAL_II_NAMESPACE_CLOSE
index b921e567d976dc8b26890e645ae45e50f11d6630..2d3b5e8b9521839ddfe06193cd9bcd2583d369e7 100644 (file)
@@ -2,7 +2,7 @@
 //    $Id$
 //    Version: $Name$
 //
-//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008 by the deal.II authors
+//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2011 by the deal.II authors
 //
 //    This file is subject to QPL and may not be  distributed
 //    without copyright and license information. Please refer
@@ -24,7 +24,7 @@ for (VEC : SERIAL_VECTORS ; deal_II_dimension : DIMENSIONS)
       (const DoFHandler<deal_II_dimension>&,
       const Function<deal_II_dimension>&,
       VEC&);
-      
+
   template 
     void VectorTools::interpolate
     (const Mapping<deal_II_dimension>&,
@@ -252,6 +252,22 @@ for (VEC : SERIAL_VECTORS ; deal_II_dimension : DIMENSIONS)
      const VEC&,
      const unsigned int);
 
+#if deal_II_dimension < 3
+  template
+    double VectorTools::compute_mean_value<deal_II_dimension>
+    (const Mapping<deal_II_dimension,deal_II_dimension+1>&,
+     const DoFHandler<deal_II_dimension,deal_II_dimension+1>&,
+     const Quadrature<deal_II_dimension>&,
+     const VEC&,
+     const unsigned int);
+  template
+    double VectorTools::compute_mean_value<deal_II_dimension>
+    (const DoFHandler<deal_II_dimension,deal_II_dimension+1>&,
+     const Quadrature<deal_II_dimension>&,
+     const VEC&,
+     const unsigned int);
+#endif
+
   template
     void VectorTools::project
     (const Mapping<deal_II_dimension>      &,
diff --git a/tests/codim_one/mean_value.cc b/tests/codim_one/mean_value.cc
new file mode 100644 (file)
index 0000000..233d328
--- /dev/null
@@ -0,0 +1,76 @@
+//----------------------------  template.cc  ---------------------------
+//    $Id$
+//    Version: $Name$ 
+//
+//    Copyright (C) 2005, 2008, 2011 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.
+//
+//----------------------------  template.cc  ---------------------------
+
+
+// continous projection of a function on the surface of a hypersphere
+
+#include "../tests.h"
+#include <fstream>
+#include <base/logstream.h>
+#include <grid/tria.h>
+#include <grid/grid_generator.h>
+#include <fe/mapping.h>
+#include <fe/mapping_q1.h>
+#include <fe/fe_q.h>
+#include <fe/fe_values.h>
+#include <base/quadrature_lib.h>
+#include <dofs/dof_handler.h>
+#include <dofs/dof_accessor.h>
+#include <lac/constraint_matrix.h>
+#include <numerics/vectors.h>
+#include <numerics/data_out.h>
+#include <base/function.h>
+#include <base/function_lib.h>
+
+#include <base/quadrature_lib.h>
+
+#include <fstream>
+#include <string>
+
+
+std::ofstream logfile("mean_value/output");
+
+
+template <int dim, int spacedim>
+void test()
+{
+  Triangulation<dim, spacedim> triangulation;
+  GridGenerator::hyper_cube (triangulation, -1, 1);
+
+  FE_Q<dim,spacedim>     fe (1);
+  DoFHandler<dim,spacedim> dof_handler (triangulation);
+
+  dof_handler.distribute_dofs (fe);
+  Functions::CosineFunction<spacedim> cosine;
+  Vector<double> x(dof_handler.n_dofs());
+  VectorTools::interpolate(dof_handler, cosine, x);
+
+  const double mean = VectorTools::compute_mean_value (dof_handler,
+                                                      QGauss<dim>(2),
+                                                      x, 0);
+  // we have a symmetric domain and a symmetric function. the result
+  // should be close to zero
+  Assert (std::fabs(mean) < 1e-15, ExcInternalError());
+  deallog << "OK" << std::endl;
+}
+
+
+
+int main () 
+{
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  
+  test<1,2>();
+  test<2,3>();
+}
diff --git a/tests/codim_one/mean_value/cmp/generic b/tests/codim_one/mean_value/cmp/generic
new file mode 100644 (file)
index 0000000..8b3b075
--- /dev/null
@@ -0,0 +1,3 @@
+
+DEAL::OK
+DEAL::OK

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.