--- /dev/null
+//---------------------------- 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;
+}
+
--- /dev/null
+
+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