]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Changed function name and added test 1285/head
authorJason Sheldon <Jason.Sheldon@psu.edu>
Fri, 7 Aug 2015 19:46:52 +0000 (15:46 -0400)
committerJason Sheldon <Jason.Sheldon@psu.edu>
Fri, 7 Aug 2015 19:46:52 +0000 (15:46 -0400)
include/deal.II/fe/mapping.h
source/fe/mapping.cc
tests/fe/mapping_project_01.cc [new file with mode: 0644]
tests/fe/mapping_project_01.output [new file with mode: 0644]

index 4cb023ce56d20069bb130feb3b26076f59a42eb2..dce14e066685f005748d4779c077721478f2ec67 100644 (file)
@@ -285,7 +285,7 @@ public:
    * so it throws an exception in this case.
    */
   virtual Point<dim-1>
-  transform_real_to_unit_projected_to_face (
+  project_real_point_to_unit_point_on_face (
     const typename Triangulation<dim,spacedim>::cell_iterator &cell,
     const unsigned int &face_no,
     const Point<spacedim> &p) const;
index f3e7432afe062eafcbd166625c38ddf3f0f77fe1..46b5378c2aba0616158b4bccf0e3b6b1da88a310 100644 (file)
@@ -42,7 +42,7 @@ Mapping<dim, spacedim>::get_vertices (
 template<int dim, int spacedim>
 Point<dim-1>
 Mapping<dim,spacedim>::
-transform_real_to_unit_projected_to_face (
+project_real_point_to_unit_point_on_face (
   const typename Triangulation<dim,spacedim>::cell_iterator &cell,
   const unsigned int &face_no,
   const Point<spacedim> &p) const
diff --git a/tests/fe/mapping_project_01.cc b/tests/fe/mapping_project_01.cc
new file mode 100644 (file)
index 0000000..dd85fa4
--- /dev/null
@@ -0,0 +1,127 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2015 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.
+//
+// ---------------------------------------------------------------------
+
+// Tests for project_real_point_to_unit_point_on_face in 2D and 3D on a cube,
+// and also in 3D with a parallelepiped, checks against expected values.
+
+#include "../tests.h"
+#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/fe/mapping_q1.h>
+#include <iostream>
+#include <fstream>
+#include <cmath>
+
+using namespace dealii;
+void dim2_grid ()
+{
+  Triangulation<2> triangulation;
+
+  const Point<2> p_ll (-1, -1);//lower left point
+  const Point<2> p_ur (1, 1);//upper right point
+  GridGenerator::hyper_rectangle (triangulation, p_ll, p_ur, false);
+
+  const Point<2> testp (.5, -.5);  //test point
+
+  MappingQ1<2> mapping;
+
+  deallog<<"Check project for 2D cube from (-1,-1), to (1,1)."<<std::endl;
+
+  Triangulation<2>::active_cell_iterator
+  cell = triangulation.begin_active(),
+  endc = triangulation.end();
+  deallog<<"The test point has real coordinates: "<<testp
+        <<", and unit coordinates: "
+        <<mapping.transform_real_to_unit_cell(cell, testp)<<std::endl;
+  for (; cell!=endc; ++cell)
+    for (unsigned int face=0; face<GeometryInfo<2>::faces_per_cell; ++face)
+    {
+      deallog<<"For face: "<<face<<", the test point projects to: "
+            <<mapping.project_real_point_to_unit_point_on_face(cell, face, testp)<<std::endl;
+    }
+}
+void dim3_grid ()
+{
+  Triangulation<3> triangulation;
+
+  const Point<3> p_ll (-1, -1, -1);//lower left point
+  const Point<3> p_ur (1, 1, 1);//upper right point
+  GridGenerator::hyper_rectangle (triangulation, p_ll, p_ur, false);
+
+  const Point<3> testp (.5, -.5, 0);  //test point
+
+  MappingQ1<3> mapping;
+
+  deallog<<"Check project for 3D cube from (-1,-1,-1) to (1,1,1)."<<std::endl;
+
+  Triangulation<3>::active_cell_iterator
+  cell = triangulation.begin_active(),
+  endc = triangulation.end();
+  deallog<<"The test point has real coordinates: "<<testp
+        <<", and unit coordinates: "
+        <<mapping.transform_real_to_unit_cell(cell, testp)<<std::endl;
+  for (; cell!=endc; ++cell)
+    for (unsigned int face=0; face<GeometryInfo<3>::faces_per_cell; ++face)
+    {
+      deallog<<"For face: "<<face<<", the test point projects to: "
+            <<mapping.project_real_point_to_unit_point_on_face(cell, face, testp)<<std::endl;
+    }
+}
+void dim3_parallelepiped_grid ()
+{
+  Triangulation<3> triangulation;
+
+  Point<3> (corners) [3];
+  corners[0] = Point<3> (2, 0, 0);
+  corners[1] = Point<3> (0, 2, 0);
+  corners[2] = Point<3> (0, 1, 2);
+
+  GridGenerator::parallelepiped (triangulation, corners);
+
+  const Point<3> testp (1, 1, 1);  //test point
+
+  MappingQ1<3> mapping;
+
+  deallog<<"Check project for 3D parallelepiped with vectors (2, 0, 0), (0, 2, 0), and (0, 1, 2)."<<std::endl;
+  Triangulation<3>::active_cell_iterator
+  cell = triangulation.begin_active(),
+  endc = triangulation.end();
+  deallog<<"The test point has real coordinates: "<<testp
+        <<", and unit coordinates: "
+        <<mapping.transform_real_to_unit_cell(cell, testp)<<std::endl;
+
+  for (; cell!=endc; ++cell)
+    for (unsigned int face=0; face<GeometryInfo<3>::faces_per_cell; ++face)
+    {
+      deallog<<"For face: "<<face<<", the test point projects to: "
+            <<mapping.project_real_point_to_unit_point_on_face(cell, face, testp)<<std::endl;
+    }
+}
+int main ()
+{
+
+  std::ofstream logfile("output");
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.precision(3);
+
+  dim2_grid ();
+  dim3_grid ();
+  dim3_parallelepiped_grid ();
+
+  return 0;
+}
diff --git a/tests/fe/mapping_project_01.output b/tests/fe/mapping_project_01.output
new file mode 100644 (file)
index 0000000..e8e989d
--- /dev/null
@@ -0,0 +1,23 @@
+
+DEAL::Check project for 2D cube from (-1,-1), to (1,1).
+DEAL::The test point has real coordinates: 0.500 -0.500, and unit coordinates: 0.750 0.250
+DEAL::For face: 0, the test point projects to: 0.250
+DEAL::For face: 1, the test point projects to: 0.250
+DEAL::For face: 2, the test point projects to: 0.750
+DEAL::For face: 3, the test point projects to: 0.750
+DEAL::Check project for 3D cube from (-1,-1,-1) to (1,1,1).
+DEAL::The test point has real coordinates: 0.500 -0.500 0.00, and unit coordinates: 0.750 0.250 0.500
+DEAL::For face: 0, the test point projects to: 0.250 0.500
+DEAL::For face: 1, the test point projects to: 0.250 0.500
+DEAL::For face: 2, the test point projects to: 0.750 0.500
+DEAL::For face: 3, the test point projects to: 0.750 0.500
+DEAL::For face: 4, the test point projects to: 0.750 0.250
+DEAL::For face: 5, the test point projects to: 0.750 0.250
+DEAL::Check project for 3D parallelepiped with vectors (2, 0, 0), (0, 2, 0), and (0, 1, 2).
+DEAL::The test point has real coordinates: 1.00 1.00 1.00, and unit coordinates: 0.500 0.250 0.500
+DEAL::For face: 0, the test point projects to: 0.250 0.500
+DEAL::For face: 1, the test point projects to: 0.250 0.500
+DEAL::For face: 2, the test point projects to: 0.500 0.500
+DEAL::For face: 3, the test point projects to: 0.500 0.500
+DEAL::For face: 4, the test point projects to: 0.500 0.250
+DEAL::For face: 5, the test point projects to: 0.500 0.250

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.