From: Wolfgang Bangerth Date: Sat, 3 Sep 2011 00:23:29 +0000 (+0000) Subject: Implement GridTools::volume. X-Git-Tag: v8.0.0~3535 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=9d1e85861fd0da20c542cf207826eddb2ac7c24b;p=dealii.git Implement GridTools::volume. git-svn-id: https://svn.dealii.org/trunk@24245 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index 7e18a31931..c234605355 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -283,6 +283,11 @@ and DoF handlers embedded in higher dimensional space have been added.

Specific improvements

    +
  1. New: There is now a function GridTools::volume computing the volume +of a triangulation. +
    +(Wolfgang Bangerth, 2011/09/02) +
  2. New: Code like cell-@>face(1)-@>set_boundary_indicator(42); now also works in 1d.
    diff --git a/deal.II/include/deal.II/grid/grid_tools.h b/deal.II/include/deal.II/grid/grid_tools.h index c57c608c9c..3e32d1f132 100644 --- a/deal.II/include/deal.II/grid/grid_tools.h +++ b/deal.II/include/deal.II/grid/grid_tools.h @@ -18,6 +18,7 @@ #include #include #include +#include #include @@ -42,7 +43,6 @@ class SparsityPattern; * individual functions for more information. * * @ingroup grid - * @author Wolfgang Bangerth, 2001, 2003, 2004, Ralf Hartmann, 2005 */ class GridTools { @@ -62,6 +62,21 @@ class GridTools template static double diameter (const Triangulation &tria); + + /** + * Compute the volume (i.e. the dim-dimensional measure) of the + * triangulation. We compute the measure using the integral + * $\int 1 \; dx$. The integral approximated is approximated + * via quadrature for which we need the mapping argument. + * + * This function also works for objects of type + * parallel::distributed::Triangulation, in which case the + * function is a collective operation. + */ + template + static + double volume (const Triangulation &tria, + const Mapping &mapping = StaticMappingQ1::mapping); /** * Given a list of vertices (typically diff --git a/deal.II/source/grid/grid_tools.cc b/deal.II/source/grid/grid_tools.cc index 08131de96a..79ea68df19 100644 --- a/deal.II/source/grid/grid_tools.cc +++ b/deal.II/source/grid/grid_tools.cc @@ -12,7 +12,7 @@ //--------------------------------------------------------------------------- - +#include #include #include #include @@ -27,6 +27,8 @@ #include #include #include +#include +#include #include #include #include @@ -78,6 +80,10 @@ template double GridTools::diameter (const Triangulation &tria) { + Assert ((dynamic_cast*>(&tria) + == 0), + ExcNotImplemented()); + // the algorithm used simply // traverses all cells and picks // out the boundary vertices. it @@ -126,6 +132,58 @@ GridTools::diameter (const Triangulation &tria) +template +double +GridTools::volume (const Triangulation &triangulation, + const Mapping &mapping) +{ + // get the degree of the mapping if possible. if not, just assume 1 + const unsigned int mapping_degree + = (dynamic_cast*>(&mapping) != 0 ? + dynamic_cast*>(&mapping)->get_degree() : + 1); + + // then initialize an appropriate quadrature formula + const QGauss quadrature_formula (mapping_degree + 1); + const unsigned int n_q_points = quadrature_formula.size(); + + FE_DGQ dummy_fe(0); + FEValues fe_values (mapping, dummy_fe, quadrature_formula, + update_JxW_values); + + typename Triangulation::active_cell_iterator + cell = triangulation.begin_active(), + endc = triangulation.end(); + + double local_volume = 0; + + // compute the integral quantities by quadrature + for (; cell!=endc; ++cell) + if (!cell->is_ghost() && !cell->is_artificial()) + { + fe_values.reinit (cell); + for (unsigned int q=0; q* p_tria + = dynamic_cast*>(&triangulation)) + MPI_Allreduce (&local_volume, &global_volume, 1, MPI_DOUBLE, + MPI_SUM, + p_tria->get_communicator()); + else + global_volume = local_volume; +#else + global_volume = local_volume; +#endif + + return global_volume; +} + + template <> double GridTools::cell_measure<3>(const std::vector > &all_vertices, diff --git a/deal.II/source/grid/grid_tools.inst.in b/deal.II/source/grid/grid_tools.inst.in index 9d0533cd67..a6dcc68690 100644 --- a/deal.II/source/grid/grid_tools.inst.in +++ b/deal.II/source/grid/grid_tools.inst.in @@ -62,6 +62,18 @@ for (deal_II_dimension : DIMENSIONS) GridTools::diameter (const Triangulation &); #endif + template + double + GridTools::volume (const Triangulation &, + const Mapping &); + +#if deal_II_dimension < 3 + template + double + GridTools::volume (const Triangulation &, + const Mapping &); +#endif + template void GridTools::delete_unused_vertices (std::vector > &, std::vector > &, diff --git a/tests/deal.II/grid_tools_03.cc b/tests/deal.II/grid_tools_03.cc new file mode 100644 index 0000000000..0659ae4235 --- /dev/null +++ b/tests/deal.II/grid_tools_03.cc @@ -0,0 +1,94 @@ +//---------------------------- grid_tools.cc --------------------------- +// $Id$ +// Version: $Name$ +// +// Copyright (C) 2001, 2002, 2003, 2004, 2005, 2006, 2008, 2011 by the deal.II authors +// +// This file is subject to QPL 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. +// +//---------------------------- grid_tools.cc --------------------------- + +// check GridTools::volume + + +#include "../tests.h" +#include +#include +#include +#include +#include +#include + +#include +#include + +std::ofstream logfile("grid_tools_03/output"); + + + +template +void test1 () +{ + // test 1: hypercube + if (true) + { + Triangulation tria; + GridGenerator::hyper_cube(tria); + + for (unsigned int i=0; i<2; ++i) + { + tria.refine_global(2); + deallog << dim << "d, " + << "hypercube volume, " + << i*2 + << " refinements: " + << GridTools::volume (tria) + << std::endl; + }; + }; + + // test 2: hyperball + if (dim >= 2) + { + Triangulation tria; + GridGenerator::hyper_ball(tria, Point(), 1); + + static const HyperBallBoundary boundary; + tria.set_boundary (0, boundary); + + for (unsigned int i=0; i<4; ++i) + { + tria.refine_global(1); + deallog << dim << "d, " + << "hyperball volume, " + << i + << " refinements: " + << GridTools::volume (tria) + << std::endl; + } + deallog << "exact value=" + << (dim==2 ? numbers::PI : + 4./3.*numbers::PI) + << std::endl; + } +} + + +int main () +{ + deallog << std::setprecision(4); + logfile << std::setprecision(4); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + test1<1> (); + test1<2> (); + test1<3> (); + + return 0; +} + diff --git a/tests/deal.II/grid_tools_03/cmp/generic b/tests/deal.II/grid_tools_03/cmp/generic new file mode 100644 index 0000000000..f785507852 --- /dev/null +++ b/tests/deal.II/grid_tools_03/cmp/generic @@ -0,0 +1,17 @@ + +DEAL::1d, hypercube volume, 0 refinements: 1.000 +DEAL::1d, hypercube volume, 2 refinements: 1.000 +DEAL::2d, hypercube volume, 0 refinements: 1.000 +DEAL::2d, hypercube volume, 2 refinements: 1.000 +DEAL::2d, hyperball volume, 0 refinements: 2.828 +DEAL::2d, hyperball volume, 1 refinements: 3.061 +DEAL::2d, hyperball volume, 2 refinements: 3.121 +DEAL::2d, hyperball volume, 3 refinements: 3.137 +DEAL::exact value=3.142 +DEAL::3d, hypercube volume, 0 refinements: 1.000 +DEAL::3d, hypercube volume, 2 refinements: 1.000 +DEAL::3d, hyperball volume, 0 refinements: 3.210 +DEAL::3d, hyperball volume, 1 refinements: 3.917 +DEAL::3d, hyperball volume, 2 refinements: 4.119 +DEAL::3d, hyperball volume, 3 refinements: 4.171 +DEAL::exact value=4.189 diff --git a/tests/deal.II/grid_tools_04.cc b/tests/deal.II/grid_tools_04.cc new file mode 100644 index 0000000000..3ea6d4685a --- /dev/null +++ b/tests/deal.II/grid_tools_04.cc @@ -0,0 +1,68 @@ +//---------------------------- grid_tools.cc --------------------------- +// $Id$ +// Version: $Name$ +// +// Copyright (C) 2001, 2002, 2003, 2004, 2005, 2006, 2008, 2011 by the deal.II authors +// +// This file is subject to QPL 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. +// +//---------------------------- grid_tools.cc --------------------------- + +// check GridTools::volume for codim-one + + +#include "../tests.h" +#include +#include +#include +#include +#include +#include + +#include +#include + +std::ofstream logfile("grid_tools_04/output"); + + + +template +void test1 () +{ + // test 1: hypercube + if (true) + { + Triangulation tria; + GridGenerator::hyper_cube(tria); + + for (unsigned int i=0; i<2; ++i) + { + tria.refine_global(2); + deallog << dim << "d, " + << "hypercube volume, " + << i*2 + << " refinements: " + << GridTools::volume (tria) + << std::endl; + }; + }; +} + + +int main () +{ + deallog << std::setprecision(4); + logfile << std::setprecision(4); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + test1<1> (); + test1<2> (); + + return 0; +} + diff --git a/tests/deal.II/grid_tools_04/cmp/generic b/tests/deal.II/grid_tools_04/cmp/generic new file mode 100644 index 0000000000..3df88e4356 --- /dev/null +++ b/tests/deal.II/grid_tools_04/cmp/generic @@ -0,0 +1,5 @@ + +DEAL::1d, hypercube volume, 0 refinements: 1.000 +DEAL::1d, hypercube volume, 2 refinements: 1.000 +DEAL::2d, hypercube volume, 0 refinements: 1.000 +DEAL::2d, hypercube volume, 2 refinements: 1.000