Tensor<1,spacedim> tangent = ((p-face->vertex(0)).norm_square() > (p-face->vertex(1)).norm_square() ?
get_tangent_vector(p, face->vertex(0)) :
-get_tangent_vector(p, face->vertex(1)));
- return cross_product_2d(tangent);
+ Tensor<1,spacedim> normal = cross_product_2d(tangent);
+ return normal/normal.norm();
}
template<>
Tensor<1,spacedim> t1,t2;
// Take the difference between p and all four vertices
- Tensor<1,spacedim> dp[4];
- double dpns[4];
- int min_index=-1;
- double min_distance = 0;
- for (unsigned int i=0; i<4; ++i)
+ int min_index=0;
+ Tensor<1,spacedim> dp = p-face->vertex(0);
+ double min_distance = dp.norm_square();
+
+ for (unsigned int i=1; i<4; ++i)
{
- dp[i] = p-face->vertex(i);
- dpns[i] = dp[i].norm_square();
- min_index = (min_index == -1 ? (int)i : dpns[i] < min_distance ? i : min_index);
- min_distance = dpns[min_index];
+ dp = p-face->vertex(i);
+ double distance = dp.norm_square();
+ if (distance < min_distance)
+ {
+ min_index = i;
+ min_distance = distance;
+ }
}
// Verify we have a valid vertex index
AssertIndexRange(min_index, 4);
break;
}
}
- return cross_product_3d(t1,t2);
+ Tensor<1,spacedim> normal = cross_product_3d(t1,t2);
+ return normal/normal.norm();
}
}
+
+template <>
+void
+Manifold<2, 2>::
+get_normals_at_vertices (const typename Triangulation<2, 2>::face_iterator &face,
+ FaceVertexNormals &n) const
+{
+ n[0] = cross_product_2d(get_tangent_vector(face->vertex(0), face->vertex(1)));
+ n[1] = -cross_product_2d(get_tangent_vector(face->vertex(1), face->vertex(0)));
+
+ n[0] /= n[0].norm();
+ n[1] /= n[1].norm();
+}
+
+
+template <>
+void
+Manifold<3, 3>::
+get_normals_at_vertices (const typename Triangulation<3, 3>::face_iterator &face,
+ FaceVertexNormals &n) const
+{
+ n[0] = cross_product_3d
+ (get_tangent_vector(face->vertex(0), face->vertex(1)),
+ get_tangent_vector(face->vertex(0), face->vertex(2)));
+
+ n[1] = cross_product_3d
+ (get_tangent_vector(face->vertex(1), face->vertex(3)),
+ get_tangent_vector(face->vertex(1), face->vertex(0)));
+
+ n[2] = cross_product_3d
+ (get_tangent_vector(face->vertex(2), face->vertex(0)),
+ get_tangent_vector(face->vertex(2), face->vertex(3)));
+
+ n[3] = cross_product_3d
+ (get_tangent_vector(face->vertex(3), face->vertex(2)),
+ get_tangent_vector(face->vertex(3), face->vertex(1)));
+
+ for (unsigned int i=0; i<4; ++i)
+ n[i] /=n[i].norm();
+}
+
+
template <int dim, int spacedim>
void
Manifold<dim, spacedim>::
//
// ---------------------------------------------------------------------
-// test GridTools::torus() and TorusManifold, output visually checked
+// test get_normals_at_vertices for a SphericalManifold.
#include "../tests.h"
#include <deal.II/base/logstream.h>
<< "set view equal xyz" << std::endl
<< "set size ratio -1" << std::endl
<< "set multiplot" << std::endl
- << "splot '-' with vectors " << std::endl;
+ << (dim == 3 ? "s" : "")
+ << "plot '-' with vectors " << std::endl;
for (typename Triangulation<dim,spacedim>::active_cell_iterator
cell = triangulation.begin_active();
set view equal xyz
set size ratio -1
set multiplot
-splot '-' with vectors
--0.707107 -0.707107 0.0707107 0.0707107
--1.83697e-16 -1.00000 1.83697e-17 0.100000
--1.83697e-16 -1.00000 1.83697e-17 0.100000
-0.707107 -0.707107 -0.0707107 0.0707107
+plot '-' with vectors
-0.707107 -0.707107 -0.0707107 -0.0707107
--1.00000 1.22465e-16 -0.100000 1.22465e-17
--1.00000 1.22465e-16 -0.100000 1.22465e-17
--0.707107 0.707107 -0.0707107 0.0707107
-0.707107 -0.707107 -0.0707107 0.0707107
-1.00000 0.00000 -0.100000 0.00000
-1.00000 0.00000 -0.100000 0.00000
-0.707107 0.707107 -0.0707107 -0.0707107
--0.707107 0.707107 -0.0707107 0.0707107
-6.12323e-17 1.00000 6.12323e-18 0.100000
-6.12323e-17 1.00000 6.12323e-18 0.100000
+-1.83697e-16 -1.00000 -1.83697e-17 -0.100000
+-1.83697e-16 -1.00000 -1.83697e-17 -0.100000
+0.707107 -0.707107 0.0707107 -0.0707107
+-0.707107 -0.707107 0.0707107 0.0707107
+-1.00000 1.22465e-16 0.100000 -1.22465e-17
+-1.00000 1.22465e-16 0.100000 -1.22465e-17
+-0.707107 0.707107 0.0707107 -0.0707107
+0.707107 -0.707107 0.0707107 -0.0707107
+1.00000 0.00000 0.100000 0.00000
+1.00000 0.00000 0.100000 0.00000
0.707107 0.707107 0.0707107 0.0707107
+-0.707107 0.707107 0.0707107 -0.0707107
+6.12323e-17 1.00000 -6.12323e-18 -0.100000
+6.12323e-17 1.00000 -6.12323e-18 -0.100000
+0.707107 0.707107 -0.0707107 -0.0707107
e
set view equal xyz
set size ratio -1