From: Timo Heister Date: Fri, 18 Mar 2016 14:20:13 +0000 (+0100) Subject: add GridTools::rotate in 3d X-Git-Tag: v8.5.0-rc1~1190^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=df01d79bb47609a5820721c817633dccaacb12ce;p=dealii.git add GridTools::rotate in 3d --- diff --git a/doc/news/changes.h b/doc/news/changes.h index e103c7af5a..7aeaec7e88 100644 --- a/doc/news/changes.h +++ b/doc/news/changes.h @@ -74,6 +74,11 @@ inconvenience this causes.

General

    +
  1. New: Added GridTools::rotate() in three space dimensions. +
    + (Timo Heister, 2016/03/18) +
  2. +
  3. New: Added unit tests for complex-valued PETSc and SLEPc.
    (Toby D. Young, Denis Davydov, 2016/03/11) diff --git a/include/deal.II/grid/grid_tools.h b/include/deal.II/grid/grid_tools.h index d42576ee2f..548f8b5d81 100644 --- a/include/deal.II/grid/grid_tools.h +++ b/include/deal.II/grid/grid_tools.h @@ -263,6 +263,24 @@ namespace GridTools 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 + void + rotate (const double angle, + const unsigned int axis, + Triangulation &triangulation); + /** * Transform the given triangulation smoothly to a different domain where, * typically, each of the vertices at the boundary of the triangulation is diff --git a/source/grid/grid_tools.cc b/source/grid/grid_tools.cc index e18cf24bfc..61653fd841 100644 --- a/source/grid/grid_tools.cc +++ b/source/grid/grid_tools.cc @@ -523,6 +523,36 @@ namespace GridTools 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 class Scale @@ -559,7 +589,16 @@ namespace GridTools transform (Rotate2d(angle), triangulation); } + template + void + rotate (const double angle, + const unsigned int axis, + Triangulation &triangulation) + { + Assert(axis<3, ExcMessage("Invalid axis given!")); + transform (Rotate3d(angle, axis), triangulation); + } template void diff --git a/source/grid/grid_tools.inst.in b/source/grid/grid_tools.inst.in index cce84f3517..542e3a5611 100644 --- a/source/grid/grid_tools.inst.in +++ b/source/grid/grid_tools.inst.in @@ -133,6 +133,14 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS void scale (const double, Triangulation &); +#if deal_II_space_dimension==3 + template + void + rotate (const double, + const unsigned int , + Triangulation &); +#endif + template void distort_random (const double, Triangulation &, diff --git a/tests/grid/rotate_01.cc b/tests/grid/rotate_01.cc new file mode 100644 index 0000000000..533726c4e2 --- /dev/null +++ b/tests/grid/rotate_01.cc @@ -0,0 +1,88 @@ +// --------------------------------------------------------------------- +// +// 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 +#include +#include +#include +#include + +#include +#include + + +template +void test (); + +template <> +void test<2,2> () +{ + const int dim = 2; + const int spacedim = 2; + Triangulation tria; + Point origin; + std_cxx11::array,dim> edges; + edges[0] = Point(2.0, 0.0)-Point(); + edges[1] = Point(0.2, 0.8)-Point(); + GridGenerator::subdivided_parallelepiped (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 +void test () +{ + Triangulation tria; + Point origin(0.1, 0.2, 0.3); + std_cxx11::array,dim> edges; + edges[0] = Point(2.0, 0.0, 0.1)-Point(); + if (dim>=2) + edges[1] = Point(0.2, 0.8, 0.15)-Point(); + if (dim>=3) + edges[2] = Point(0.0, 0.0, 0.1)-Point(); + + GridGenerator::subdivided_parallelepiped (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; +} diff --git a/tests/grid/rotate_01.output b/tests/grid/rotate_01.output new file mode 100644 index 0000000000..cc6fd4fb74 --- /dev/null +++ b/tests/grid/rotate_01.output @@ -0,0 +1,43 @@ + +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 +