]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added tangent vector for FunctionManifold
authorLuca Heltai <luca.heltai@sissa.it>
Fri, 18 Mar 2016 19:10:02 +0000 (20:10 +0100)
committerLuca Heltai <luca.heltai@sissa.it>
Sun, 20 Mar 2016 11:33:46 +0000 (12:33 +0100)
doc/news/changes.h
include/deal.II/base/function_parser.h
include/deal.II/grid/manifold_lib.h
source/base/function_parser.cc
source/grid/manifold_lib.cc
tests/manifold/function_manifold_03.cc [new file with mode: 0644]
tests/manifold/function_manifold_03.with_muparser=true.output [new file with mode: 0644]

index e103c7af5ae93f1e13e03b08f1e8888f394a258c..503c06e1997914fdf823dad3c48c49f0a0b1270c 100644 (file)
@@ -98,7 +98,6 @@ inconvenience this causes.
 <a name="specific"></a>
 <h3>Specific improvements</h3>
 
-
 <ol>
  <li> Improved: The distribution of degrees of freedom on multigrid levels,
  DoFHandler::distribute_mg_dofs(), contained a few steps that scaled
index 5ae72919dfc04c3cd63d35df93364d7858f766ee..3fe91a0c4e161df8ee22b5f8f59195f71c9856c1 100644 (file)
@@ -19,7 +19,7 @@
 
 #include <deal.II/base/config.h>
 #include <deal.II/base/exceptions.h>
-#include <deal.II/base/function.h>
+#include <deal.II/base/auto_derivative_function.h>
 #include <deal.II/base/tensor.h>
 #include <deal.II/base/point.h>
 #include <deal.II/base/thread_local_storage.h>
@@ -177,7 +177,7 @@ template <typename> class Vector;
  * @author Luca Heltai, Timo Heister 2005, 2014
  */
 template <int dim>
-class FunctionParser : public Function<dim>
+class FunctionParser : public AutoDerivativeFunction<dim>
 {
 public:
   /**
index 21de7b05552a1eefca372bd3609ec65b09f8826a..15b3d7e4dc67258e9979a07d660228ee7cbc608e 100644 (file)
@@ -265,6 +265,22 @@ public:
   virtual Point<spacedim>
   push_forward(const Point<chartdim> &chart_point) const;
 
+  /**
+   * Given a point in the chartdim dimensional Euclidean space, this
+   * method returns the derivatives of the function $F$ that maps from
+   * the sub_manifold coordinate system to the Euclidean coordinate
+   * system. In other words, it is a matrix of size
+   * $\text{spacedim}\times\text{chartdim}$.
+   *
+   * This function is used in the computations required by the
+   * get_tangent_vector() function.
+   *
+   * Refer to the general documentation of this class for more information.
+   */
+  virtual
+  DerivativeForm<1,chartdim,spacedim>
+  push_forward_gradient(const Point<chartdim> &chart_point) const;
+
   /**
    * Given a point in the spacedim coordinate system, uses the
    * pull_back_function to compute the pull_back of points in @p spacedim
index 50703b27a67ff805426ce68b46959299210ac52a..094869b40fae57d81e28e8240ad83d4b5f3d3230 100644 (file)
@@ -48,7 +48,7 @@ template <int dim>
 FunctionParser<dim>::FunctionParser(const unsigned int n_components,
                                     const double       initial_time)
   :
-  Function<dim>(n_components, initial_time)
+  AutoDerivativeFunction<dim>(1e-8, n_components, initial_time)
 {}
 
 
index 4740f4e65c6f77ec8926e31061f8db0d3a57c07e..0d039f8c8683c1dc24c2bd940289895d2d4115fd 100644 (file)
@@ -320,6 +320,20 @@ FunctionManifold<dim,spacedim,chartdim>::push_forward(const Point<chartdim> &cha
 }
 
 
+template <int dim, int spacedim, int chartdim>
+DerivativeForm<1,chartdim, spacedim>
+FunctionManifold<dim,spacedim,chartdim>::push_forward_gradient(const Point<chartdim> &chart_point) const
+{
+  DerivativeForm<1, chartdim, spacedim> DF;
+  std::vector<Tensor<1, chartdim> > gradients(spacedim);
+  push_forward_function->vector_gradient(chart_point, gradients);
+  for (unsigned int i=0; i<spacedim; ++i)
+    for (unsigned int j=0; j<chartdim; ++j)
+      DF[i][j] = gradients[i][j];
+  return DF;
+}
+
+
 template <int dim, int spacedim, int chartdim>
 Point<chartdim>
 FunctionManifold<dim,spacedim,chartdim>::pull_back(const Point<spacedim> &space_point) const
diff --git a/tests/manifold/function_manifold_03.cc b/tests/manifold/function_manifold_03.cc
new file mode 100644 (file)
index 0000000..a698f7c
--- /dev/null
@@ -0,0 +1,96 @@
+//----------------------------  function_manifold_chart ---------------------------
+//    Copyright (C) 2011 - 2015 by the mathLab team.
+//
+//    This file is subject to LGPL 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.
+//
+//---------------------------- function_manifold_chart ---------------------------
+
+
+// Test a simple parabolic manifold, including gradients and tangent vector
+
+#include "../tests.h"
+#include <fstream>
+#include <deal.II/base/logstream.h>
+
+
+// all include files you need here
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria_boundary_lib.h>
+#include <deal.II/grid/manifold_lib.h>
+#include <deal.II/grid/grid_out.h>
+
+// Helper function
+template <int dim, int spacedim>
+void test(unsigned int ref=1)
+{
+  deallog << "Testing dim " << dim
+          << ", spacedim " << spacedim << std::endl;
+
+  // Here the only allowed axis is z. In cylinder the default is x.
+  std::string push_forward_expression;
+  std::string pull_back_expression;
+
+  switch (spacedim)
+    {
+    case 2:
+      push_forward_expression = "x; x^2";
+      pull_back_expression = "x";
+      break;
+    case 3:
+      push_forward_expression = "x; x^2; 0";
+      pull_back_expression = "x";
+      break;
+    default:
+      Assert(false, ExcInternalError());
+    }
+
+  FunctionManifold<dim,spacedim,1> manifold(push_forward_expression,
+                                            pull_back_expression);
+
+  // Two points and two weights
+  std::vector<Point<spacedim> > p(2);
+  p[1][0] = 1.0;
+  p[1][1] = 1.0;
+  std::vector<double> w(2);
+
+  unsigned int n_intermediates = 16;
+
+  deallog << "P0: " << p[0]
+          << ", P1: " << p[1] << std::endl;
+
+  for (unsigned int i=0; i<n_intermediates+1; ++i)
+    {
+      w[0] = 1.0-(double)i/((double)n_intermediates);
+      w[1] = 1.0 - w[0];
+
+      Point<spacedim> ip = manifold.get_new_point(Quadrature<spacedim>(p, w));
+      Tensor<1,spacedim> t1 = manifold.get_tangent_vector(ip, p[0]);
+      Tensor<1,spacedim> t2 = manifold.get_tangent_vector(ip, p[1]);
+
+      deallog << "P: " << ip
+              << ", T(P, P0): " << t1
+              << ", T(P, P1): " << t2 << std::endl;
+
+    }
+}
+
+int main ()
+{
+  std::ofstream logfile("output");
+  deallog.attach(logfile);
+  deallog.threshold_double(1.e-10);
+
+
+  test<2,2>();
+  test<2,3>();
+  test<3,3>();
+
+  return 0;
+}
+
diff --git a/tests/manifold/function_manifold_03.with_muparser=true.output b/tests/manifold/function_manifold_03.with_muparser=true.output
new file mode 100644 (file)
index 0000000..c901d84
--- /dev/null
@@ -0,0 +1,58 @@
+
+DEAL::Testing dim 2, spacedim 2
+DEAL::P0: 0.00000 0.00000, P1: 1.00000 1.00000
+DEAL::P: 0.00000 0.00000, T(P, P0): 0.00000 0.00000, T(P, P1): 1.00000 0.00000
+DEAL::P: 0.0625000 0.00390625, T(P, P0): -0.0625000 -0.00781250, T(P, P1): 0.937500 0.117187
+DEAL::P: 0.125000 0.0156250, T(P, P0): -0.125000 -0.0312500, T(P, P1): 0.875000 0.218750
+DEAL::P: 0.187500 0.0351562, T(P, P0): -0.187500 -0.0703125, T(P, P1): 0.812500 0.304687
+DEAL::P: 0.250000 0.0625000, T(P, P0): -0.250000 -0.125000, T(P, P1): 0.750000 0.375000
+DEAL::P: 0.312500 0.0976562, T(P, P0): -0.312500 -0.195312, T(P, P1): 0.687500 0.429687
+DEAL::P: 0.375000 0.140625, T(P, P0): -0.375000 -0.281250, T(P, P1): 0.625000 0.468750
+DEAL::P: 0.437500 0.191406, T(P, P0): -0.437500 -0.382812, T(P, P1): 0.562500 0.492187
+DEAL::P: 0.500000 0.250000, T(P, P0): -0.500000 -0.500000, T(P, P1): 0.500000 0.500000
+DEAL::P: 0.562500 0.316406, T(P, P0): -0.562500 -0.632813, T(P, P1): 0.437500 0.492188
+DEAL::P: 0.625000 0.390625, T(P, P0): -0.625000 -0.781250, T(P, P1): 0.375000 0.468750
+DEAL::P: 0.687500 0.472656, T(P, P0): -0.687500 -0.945313, T(P, P1): 0.312500 0.429688
+DEAL::P: 0.750000 0.562500, T(P, P0): -0.750000 -1.12500, T(P, P1): 0.250000 0.375000
+DEAL::P: 0.812500 0.660156, T(P, P0): -0.812500 -1.32031, T(P, P1): 0.187500 0.304688
+DEAL::P: 0.875000 0.765625, T(P, P0): -0.875000 -1.53125, T(P, P1): 0.125000 0.218750
+DEAL::P: 0.937500 0.878906, T(P, P0): -0.937500 -1.75781, T(P, P1): 0.0625000 0.117188
+DEAL::P: 1.00000 1.00000, T(P, P0): -1.00000 -2.00000, T(P, P1): 0.00000 0.00000
+DEAL::Testing dim 2, spacedim 3
+DEAL::P0: 0.00000 0.00000 0.00000, P1: 1.00000 1.00000 0.00000
+DEAL::P: 0.00000 0.00000 0.00000, T(P, P0): 0.00000 0.00000 0.00000, T(P, P1): 1.00000 0.00000 0.00000
+DEAL::P: 0.0625000 0.00390625 0.00000, T(P, P0): -0.0625000 -0.00781250 0.00000, T(P, P1): 0.937500 0.117187 0.00000
+DEAL::P: 0.125000 0.0156250 0.00000, T(P, P0): -0.125000 -0.0312500 0.00000, T(P, P1): 0.875000 0.218750 0.00000
+DEAL::P: 0.187500 0.0351562 0.00000, T(P, P0): -0.187500 -0.0703125 0.00000, T(P, P1): 0.812500 0.304687 0.00000
+DEAL::P: 0.250000 0.0625000 0.00000, T(P, P0): -0.250000 -0.125000 0.00000, T(P, P1): 0.750000 0.375000 0.00000
+DEAL::P: 0.312500 0.0976562 0.00000, T(P, P0): -0.312500 -0.195312 0.00000, T(P, P1): 0.687500 0.429687 0.00000
+DEAL::P: 0.375000 0.140625 0.00000, T(P, P0): -0.375000 -0.281250 0.00000, T(P, P1): 0.625000 0.468750 0.00000
+DEAL::P: 0.437500 0.191406 0.00000, T(P, P0): -0.437500 -0.382812 0.00000, T(P, P1): 0.562500 0.492187 0.00000
+DEAL::P: 0.500000 0.250000 0.00000, T(P, P0): -0.500000 -0.500000 0.00000, T(P, P1): 0.500000 0.500000 0.00000
+DEAL::P: 0.562500 0.316406 0.00000, T(P, P0): -0.562500 -0.632813 0.00000, T(P, P1): 0.437500 0.492188 0.00000
+DEAL::P: 0.625000 0.390625 0.00000, T(P, P0): -0.625000 -0.781250 0.00000, T(P, P1): 0.375000 0.468750 0.00000
+DEAL::P: 0.687500 0.472656 0.00000, T(P, P0): -0.687500 -0.945313 0.00000, T(P, P1): 0.312500 0.429688 0.00000
+DEAL::P: 0.750000 0.562500 0.00000, T(P, P0): -0.750000 -1.12500 0.00000, T(P, P1): 0.250000 0.375000 0.00000
+DEAL::P: 0.812500 0.660156 0.00000, T(P, P0): -0.812500 -1.32031 0.00000, T(P, P1): 0.187500 0.304688 0.00000
+DEAL::P: 0.875000 0.765625 0.00000, T(P, P0): -0.875000 -1.53125 0.00000, T(P, P1): 0.125000 0.218750 0.00000
+DEAL::P: 0.937500 0.878906 0.00000, T(P, P0): -0.937500 -1.75781 0.00000, T(P, P1): 0.0625000 0.117188 0.00000
+DEAL::P: 1.00000 1.00000 0.00000, T(P, P0): -1.00000 -2.00000 0.00000, T(P, P1): 0.00000 0.00000 0.00000
+DEAL::Testing dim 3, spacedim 3
+DEAL::P0: 0.00000 0.00000 0.00000, P1: 1.00000 1.00000 0.00000
+DEAL::P: 0.00000 0.00000 0.00000, T(P, P0): 0.00000 0.00000 0.00000, T(P, P1): 1.00000 0.00000 0.00000
+DEAL::P: 0.0625000 0.00390625 0.00000, T(P, P0): -0.0625000 -0.00781250 0.00000, T(P, P1): 0.937500 0.117187 0.00000
+DEAL::P: 0.125000 0.0156250 0.00000, T(P, P0): -0.125000 -0.0312500 0.00000, T(P, P1): 0.875000 0.218750 0.00000
+DEAL::P: 0.187500 0.0351562 0.00000, T(P, P0): -0.187500 -0.0703125 0.00000, T(P, P1): 0.812500 0.304687 0.00000
+DEAL::P: 0.250000 0.0625000 0.00000, T(P, P0): -0.250000 -0.125000 0.00000, T(P, P1): 0.750000 0.375000 0.00000
+DEAL::P: 0.312500 0.0976562 0.00000, T(P, P0): -0.312500 -0.195312 0.00000, T(P, P1): 0.687500 0.429687 0.00000
+DEAL::P: 0.375000 0.140625 0.00000, T(P, P0): -0.375000 -0.281250 0.00000, T(P, P1): 0.625000 0.468750 0.00000
+DEAL::P: 0.437500 0.191406 0.00000, T(P, P0): -0.437500 -0.382812 0.00000, T(P, P1): 0.562500 0.492187 0.00000
+DEAL::P: 0.500000 0.250000 0.00000, T(P, P0): -0.500000 -0.500000 0.00000, T(P, P1): 0.500000 0.500000 0.00000
+DEAL::P: 0.562500 0.316406 0.00000, T(P, P0): -0.562500 -0.632813 0.00000, T(P, P1): 0.437500 0.492188 0.00000
+DEAL::P: 0.625000 0.390625 0.00000, T(P, P0): -0.625000 -0.781250 0.00000, T(P, P1): 0.375000 0.468750 0.00000
+DEAL::P: 0.687500 0.472656 0.00000, T(P, P0): -0.687500 -0.945313 0.00000, T(P, P1): 0.312500 0.429688 0.00000
+DEAL::P: 0.750000 0.562500 0.00000, T(P, P0): -0.750000 -1.12500 0.00000, T(P, P1): 0.250000 0.375000 0.00000
+DEAL::P: 0.812500 0.660156 0.00000, T(P, P0): -0.812500 -1.32031 0.00000, T(P, P1): 0.187500 0.304688 0.00000
+DEAL::P: 0.875000 0.765625 0.00000, T(P, P0): -0.875000 -1.53125 0.00000, T(P, P1): 0.125000 0.218750 0.00000
+DEAL::P: 0.937500 0.878906 0.00000, T(P, P0): -0.937500 -1.75781 0.00000, T(P, P1): 0.0625000 0.117188 0.00000
+DEAL::P: 1.00000 1.00000 0.00000, T(P, P0): -1.00000 -2.00000 0.00000, T(P, P1): 0.00000 0.00000 0.00000

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.