<h3>General</h3>
<ol>
+ <li> New: Added GridTools::rotate() in three space dimensions.
+ <br>
+ (Timo Heister, 2016/03/18)
+ </li>
+
<li> New: Added unit tests for complex-valued PETSc and SLEPc.
<br>
(Toby D. Young, Denis Davydov, 2016/03/11)
void rotate (const double angle,
Triangulation<2> &triangulation);
+ /**
+ * Rotate all vertices of the given @p triangulation in counter-clockwise
+ * direction around the axis with the given index. Otherwise like the
+ * function above.
+ *
+ * @param[in] angle Angle in radians to rotate the Triangulation by.
+ * @param[in] axis Index of the coordinate axis to rotate around, keeping
+ * that coordinate fixed (0=x axis, 1=y axis, 2=z axis).
+ * @param[in,out] triangulation The Triangulation object to rotate.
+ *
+ * @note Implemented for dim=1, 2, and 3.
+ */
+ template<int dim>
+ void
+ rotate (const double angle,
+ const unsigned int axis,
+ Triangulation<dim,3> &triangulation);
+
/**
* Transform the given triangulation smoothly to a different domain where,
* typically, each of the vertices at the boundary of the triangulation is
const double angle;
};
+ // Transformation to rotate around one of the cartesian axes.
+ class Rotate3d
+ {
+ public:
+ Rotate3d (const double angle,
+ const unsigned int axis)
+ :
+ angle(angle),
+ axis(axis)
+ {}
+
+ Point<3> operator() (const Point<3> &p) const
+ {
+ if (axis==0)
+ return Point<3> (p(0),
+ std::cos(angle)*p(1) - std::sin(angle) * p(2),
+ std::sin(angle)*p(1) + std::cos(angle) * p(2));
+ else if (axis==1)
+ return Point<3> (std::cos(angle)*p(0) + std::sin(angle) * p(2),
+ p(1),
+ -std::sin(angle)*p(0) + std::cos(angle) * p(2));
+ else
+ return Point<3> (std::cos(angle)*p(0) - std::sin(angle) * p(1),
+ std::sin(angle)*p(0) + std::cos(angle) * p(1),
+ p(2));
+ }
+ private:
+ const double angle;
+ const unsigned int axis;
+ };
template <int spacedim>
class Scale
transform (Rotate2d(angle), triangulation);
}
+ template<int dim>
+ void
+ rotate (const double angle,
+ const unsigned int axis,
+ Triangulation<dim,3> &triangulation)
+ {
+ Assert(axis<3, ExcMessage("Invalid axis given!"));
+ transform (Rotate3d(angle, axis), triangulation);
+ }
template <int dim, int spacedim>
void
void scale<deal_II_dimension> (const double,
Triangulation<deal_II_dimension, deal_II_space_dimension> &);
+#if deal_II_space_dimension==3
+ template
+ void
+ rotate<deal_II_dimension> (const double,
+ const unsigned int ,
+ Triangulation<deal_II_dimension, deal_II_space_dimension> &);
+#endif
+
template
void distort_random<deal_II_dimension> (const double,
Triangulation<deal_II_dimension, deal_II_space_dimension> &,
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+// test GridTools::rotate(), output is checked for correctness visually
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/grid_out.h>
+
+#include <fstream>
+#include <iomanip>
+
+
+template <int dim, int spacedim>
+void test ();
+
+template <>
+void test<2,2> ()
+{
+ const int dim = 2;
+ const int spacedim = 2;
+ Triangulation<dim, spacedim> tria;
+ Point<spacedim> origin;
+ std_cxx11::array<Tensor<1,spacedim>,dim> edges;
+ edges[0] = Point<spacedim>(2.0, 0.0)-Point<spacedim>();
+ edges[1] = Point<spacedim>(0.2, 0.8)-Point<spacedim>();
+ GridGenerator::subdivided_parallelepiped<dim, spacedim> (tria,
+ origin,
+ edges);
+
+ //GridOut().write_gnuplot (tria, deallog.get_file_stream());
+ GridTools::rotate(numbers::PI/3.0, tria);
+ GridOut().write_gnuplot (tria, deallog.get_file_stream());
+}
+
+// version for <1,3>, <2,3>, and <3,3>
+template <int dim, int spacedim>
+void test ()
+{
+ Triangulation<dim, spacedim> tria;
+ Point<spacedim> origin(0.1, 0.2, 0.3);
+ std_cxx11::array<Tensor<1,spacedim>,dim> edges;
+ edges[0] = Point<spacedim>(2.0, 0.0, 0.1)-Point<spacedim>();
+ if (dim>=2)
+ edges[1] = Point<spacedim>(0.2, 0.8, 0.15)-Point<spacedim>();
+ if (dim>=3)
+ edges[2] = Point<spacedim>(0.0, 0.0, 0.1)-Point<spacedim>();
+
+ GridGenerator::subdivided_parallelepiped<dim, spacedim> (tria,
+ origin,
+ edges);
+
+ //GridOut().write_gnuplot (tria, deallog.get_file_stream());
+ GridTools::rotate(numbers::PI/3.0, 0, tria);
+ //GridOut().write_gnuplot (tria, deallog.get_file_stream());
+ GridTools::rotate(numbers::PI/10.0, 1, tria);
+ //GridOut().write_gnuplot (tria, deallog.get_file_stream());
+ GridTools::rotate(-numbers::PI/5.0, 2, tria);
+ GridOut().write_gnuplot (tria, deallog.get_file_stream());
+}
+
+
+int main ()
+{
+ initlog ();
+
+ test<2,2> ();
+ test<1,3> ();
+ test<2,3> ();
+ test<3,3> ();
+
+ return 0;
+}
--- /dev/null
+
+0.00000 0.00000 0 0
+1.00000 1.73205 0 0
+0.407180 2.30526 0 0
+-0.592820 0.573205 0 0
+0.00000 0.00000 0 0
+
+
+0.0638108 -0.243894 0.276485 0 0
+1.56425 -1.44107 -0.293997 0 0
+
+
+0.0638108 -0.243894 0.276485 0 0
+1.56425 -1.44107 -0.293997 0 0
+2.06885 -1.47383 0.374441 0 0
+0.568409 -0.276649 0.944922 0 0
+0.0638108 -0.243894 0.276485 0 0
+
+
+0.0638108 -0.243894 0.276485 0 0
+1.56425 -1.44107 -0.293997 0 0
+1.52585 -1.52022 -0.246444 0 0
+0.0254071 -0.323039 0.324037 0 0
+0.0638108 -0.243894 0.276485 0 0
+
+0.568409 -0.276649 0.944922 0 0
+2.06885 -1.47383 0.374441 0 0
+2.03044 -1.55297 0.421993 0 0
+0.530005 -0.355794 0.992475 0 0
+0.568409 -0.276649 0.944922 0 0
+
+0.0638108 -0.243894 0.276485 0 0
+0.568409 -0.276649 0.944922 0 0
+
+1.56425 -1.44107 -0.293997 0 0
+2.06885 -1.47383 0.374441 0 0
+
+1.52585 -1.52022 -0.246444 0 0
+2.03044 -1.55297 0.421993 0 0
+
+0.0254071 -0.323039 0.324037 0 0
+0.530005 -0.355794 0.992475 0 0
+