-
Notifications
You must be signed in to change notification settings - Fork 707
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Add ReferenceCell::volume() and ::barycenter(). #13678
Changes from 3 commits
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,5 @@ | ||
New: There are now functions ReferenceCell::volume() and | ||
ReferenceCells::barycenter() that compute what their names suggest for | ||
the various kinds of reference cells deal.II supports. | ||
<br> | ||
(Wolfgang Bangerth, 2022/05/04) |
Original file line number | Diff line number | Diff line change | ||||||||||||||||||||||||||||||||||||||||||||
---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|
@@ -0,0 +1,118 @@ | ||||||||||||||||||||||||||||||||||||||||||||||
// --------------------------------------------------------------------- | ||||||||||||||||||||||||||||||||||||||||||||||
// | ||||||||||||||||||||||||||||||||||||||||||||||
// Copyright (C) 2020 - 2021 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.md at | ||||||||||||||||||||||||||||||||||||||||||||||
// the top level directory of deal.II. | ||||||||||||||||||||||||||||||||||||||||||||||
// | ||||||||||||||||||||||||||||||||||||||||||||||
// --------------------------------------------------------------------- | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
// Test ReferenceCell::volume() | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
#include <deal.II/base/quadrature_lib.h> | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_pyramid_p.h> | ||||||||||||||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_q.h> | ||||||||||||||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_simplex_p.h> | ||||||||||||||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_simplex_p_bubbles.h> | ||||||||||||||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_system.h> | ||||||||||||||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_values.h> | ||||||||||||||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_wedge_p.h> | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
#include <deal.II/grid/grid_generator.h> | ||||||||||||||||||||||||||||||||||||||||||||||
#include <deal.II/grid/tria.h> | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
#include "../tests.h" | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
using namespace dealii; | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
template <int dim> | ||||||||||||||||||||||||||||||||||||||||||||||
void | ||||||||||||||||||||||||||||||||||||||||||||||
test(const ReferenceCell &reference_cell) | ||||||||||||||||||||||||||||||||||||||||||||||
{ | ||||||||||||||||||||||||||||||||||||||||||||||
Triangulation<dim> triangulation; | ||||||||||||||||||||||||||||||||||||||||||||||
GridGenerator::reference_cell(triangulation, reference_cell); | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
std::unique_ptr<Quadrature<dim>> q; | ||||||||||||||||||||||||||||||||||||||||||||||
if (reference_cell == ReferenceCells::Line) | ||||||||||||||||||||||||||||||||||||||||||||||
q = std::make_unique<QGauss<dim>>(2); | ||||||||||||||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Triangle) | ||||||||||||||||||||||||||||||||||||||||||||||
q = std::make_unique<QGaussSimplex<dim>>(2); | ||||||||||||||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Quadrilateral) | ||||||||||||||||||||||||||||||||||||||||||||||
q = std::make_unique<QGauss<dim>>(2); | ||||||||||||||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Tetrahedron) | ||||||||||||||||||||||||||||||||||||||||||||||
q = std::make_unique<QGaussSimplex<dim>>(2); | ||||||||||||||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Wedge) | ||||||||||||||||||||||||||||||||||||||||||||||
q = std::make_unique<QGaussWedge<dim>>(2); | ||||||||||||||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Pyramid) | ||||||||||||||||||||||||||||||||||||||||||||||
q = std::make_unique<QGaussPyramid<dim>>(2); | ||||||||||||||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Hexahedron) | ||||||||||||||||||||||||||||||||||||||||||||||
q = std::make_unique<QGauss<dim>>(2); | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
std::unique_ptr<FiniteElement<dim>> fe; | ||||||||||||||||||||||||||||||||||||||||||||||
if (reference_cell == ReferenceCells::Line) | ||||||||||||||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_Q<dim>>(1); | ||||||||||||||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Triangle) | ||||||||||||||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_SimplexP<dim>>(1); | ||||||||||||||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Quadrilateral) | ||||||||||||||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_Q<dim>>(1); | ||||||||||||||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Tetrahedron) | ||||||||||||||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_SimplexP<dim>>(1); | ||||||||||||||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Wedge) | ||||||||||||||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_WedgeP<dim>>(1); | ||||||||||||||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Pyramid) | ||||||||||||||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_PyramidP<dim>>(1); | ||||||||||||||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Hexahedron) | ||||||||||||||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_Q<dim>>(1); | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
// Set up the objects to compute an integral on the reference cell | ||||||||||||||||||||||||||||||||||||||||||||||
FEValues<dim> fe_values(*fe, *q, update_JxW_values); | ||||||||||||||||||||||||||||||||||||||||||||||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Sorry for nit-picking, but shouldn't this work:
Suggested change
Since you never evaluate any of the shape functions? There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. I would hope this does not work ;) But there is a |
||||||||||||||||||||||||||||||||||||||||||||||
fe_values.reinit(triangulation.begin_active()); | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
double volume = 0; | ||||||||||||||||||||||||||||||||||||||||||||||
for (unsigned int i = 0; i < q->size(); ++i) | ||||||||||||||||||||||||||||||||||||||||||||||
volume += fe_values.JxW(i); | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
deallog << "ReferenceCell: " << reference_cell.to_string() << std::endl; | ||||||||||||||||||||||||||||||||||||||||||||||
deallog << " computed volume = " << volume << std::endl; | ||||||||||||||||||||||||||||||||||||||||||||||
deallog << " self-reported volume = " << reference_cell.volume() | ||||||||||||||||||||||||||||||||||||||||||||||
<< std::endl; | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
Assert(std::fabs(volume - reference_cell.volume()) < 1e-12, | ||||||||||||||||||||||||||||||||||||||||||||||
ExcInternalError()); | ||||||||||||||||||||||||||||||||||||||||||||||
} | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
int | ||||||||||||||||||||||||||||||||||||||||||||||
main() | ||||||||||||||||||||||||||||||||||||||||||||||
{ | ||||||||||||||||||||||||||||||||||||||||||||||
initlog(); | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
{ | ||||||||||||||||||||||||||||||||||||||||||||||
deallog.push("1D"); | ||||||||||||||||||||||||||||||||||||||||||||||
test<1>(ReferenceCells::Line); | ||||||||||||||||||||||||||||||||||||||||||||||
deallog.pop(); | ||||||||||||||||||||||||||||||||||||||||||||||
} | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
{ | ||||||||||||||||||||||||||||||||||||||||||||||
deallog.push("2D"); | ||||||||||||||||||||||||||||||||||||||||||||||
test<2>(ReferenceCells::Quadrilateral); | ||||||||||||||||||||||||||||||||||||||||||||||
test<2>(ReferenceCells::Triangle); | ||||||||||||||||||||||||||||||||||||||||||||||
deallog.pop(); | ||||||||||||||||||||||||||||||||||||||||||||||
} | ||||||||||||||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||||||||||||||
{ | ||||||||||||||||||||||||||||||||||||||||||||||
deallog.push("3D"); | ||||||||||||||||||||||||||||||||||||||||||||||
test<3>(ReferenceCells::Tetrahedron); | ||||||||||||||||||||||||||||||||||||||||||||||
test<3>(ReferenceCells::Pyramid); | ||||||||||||||||||||||||||||||||||||||||||||||
test<3>(ReferenceCells::Wedge); | ||||||||||||||||||||||||||||||||||||||||||||||
test<3>(ReferenceCells::Hexahedron); | ||||||||||||||||||||||||||||||||||||||||||||||
deallog.pop(); | ||||||||||||||||||||||||||||||||||||||||||||||
} | ||||||||||||||||||||||||||||||||||||||||||||||
} |
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,22 @@ | ||
|
||
DEAL:1D::ReferenceCell: Line | ||
DEAL:1D:: computed volume = 1.00000 | ||
DEAL:1D:: self-reported volume = 1.00000 | ||
DEAL:2D::ReferenceCell: Quad | ||
DEAL:2D:: computed volume = 1.00000 | ||
DEAL:2D:: self-reported volume = 1.00000 | ||
DEAL:2D::ReferenceCell: Tri | ||
DEAL:2D:: computed volume = 0.500000 | ||
DEAL:2D:: self-reported volume = 0.500000 | ||
DEAL:3D::ReferenceCell: Tet | ||
DEAL:3D:: computed volume = 0.166667 | ||
DEAL:3D:: self-reported volume = 0.166667 | ||
DEAL:3D::ReferenceCell: Pyramid | ||
DEAL:3D:: computed volume = 1.33333 | ||
DEAL:3D:: self-reported volume = 1.33333 | ||
DEAL:3D::ReferenceCell: Wedge | ||
DEAL:3D:: computed volume = 0.500000 | ||
DEAL:3D:: self-reported volume = 0.500000 | ||
DEAL:3D::ReferenceCell: Hex | ||
DEAL:3D:: computed volume = 1.00000 | ||
DEAL:3D:: self-reported volume = 1.00000 |
Original file line number | Diff line number | Diff line change | ||||||||||||||||||||||||||||||||
---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|
@@ -0,0 +1,125 @@ | ||||||||||||||||||||||||||||||||||
// --------------------------------------------------------------------- | ||||||||||||||||||||||||||||||||||
// | ||||||||||||||||||||||||||||||||||
// Copyright (C) 2020 - 2021 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.md at | ||||||||||||||||||||||||||||||||||
// the top level directory of deal.II. | ||||||||||||||||||||||||||||||||||
// | ||||||||||||||||||||||||||||||||||
// --------------------------------------------------------------------- | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
// Test ReferenceCell::barycenter() | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
#include <deal.II/base/quadrature_lib.h> | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_pyramid_p.h> | ||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_q.h> | ||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_simplex_p.h> | ||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_simplex_p_bubbles.h> | ||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_system.h> | ||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_values.h> | ||||||||||||||||||||||||||||||||||
#include <deal.II/fe/fe_wedge_p.h> | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
#include <deal.II/grid/grid_generator.h> | ||||||||||||||||||||||||||||||||||
#include <deal.II/grid/tria.h> | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
#include "../tests.h" | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
using namespace dealii; | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
template <int dim> | ||||||||||||||||||||||||||||||||||
void | ||||||||||||||||||||||||||||||||||
test(const ReferenceCell &reference_cell) | ||||||||||||||||||||||||||||||||||
{ | ||||||||||||||||||||||||||||||||||
Triangulation<dim> triangulation; | ||||||||||||||||||||||||||||||||||
GridGenerator::reference_cell(triangulation, reference_cell); | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
std::unique_ptr<Quadrature<dim>> q; | ||||||||||||||||||||||||||||||||||
if (reference_cell == ReferenceCells::Line) | ||||||||||||||||||||||||||||||||||
q = std::make_unique<QGauss<dim>>(2); | ||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Triangle) | ||||||||||||||||||||||||||||||||||
q = std::make_unique<QGaussSimplex<dim>>(2); | ||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Quadrilateral) | ||||||||||||||||||||||||||||||||||
q = std::make_unique<QGauss<dim>>(2); | ||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Tetrahedron) | ||||||||||||||||||||||||||||||||||
q = std::make_unique<QGaussSimplex<dim>>(2); | ||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Wedge) | ||||||||||||||||||||||||||||||||||
q = std::make_unique<QGaussWedge<dim>>(2); | ||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Pyramid) | ||||||||||||||||||||||||||||||||||
q = std::make_unique<QGaussPyramid<dim>>(2); | ||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Hexahedron) | ||||||||||||||||||||||||||||||||||
q = std::make_unique<QGauss<dim>>(2); | ||||||||||||||||||||||||||||||||||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. and here |
||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
std::unique_ptr<FiniteElement<dim>> fe; | ||||||||||||||||||||||||||||||||||
if (reference_cell == ReferenceCells::Line) | ||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_Q<dim>>(1); | ||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Triangle) | ||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_SimplexP<dim>>(1); | ||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Quadrilateral) | ||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_Q<dim>>(1); | ||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Tetrahedron) | ||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_SimplexP<dim>>(1); | ||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Wedge) | ||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_WedgeP<dim>>(1); | ||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Pyramid) | ||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_PyramidP<dim>>(1); | ||||||||||||||||||||||||||||||||||
else if (reference_cell == ReferenceCells::Hexahedron) | ||||||||||||||||||||||||||||||||||
fe = std::make_unique<FE_Q<dim>>(1); | ||||||||||||||||||||||||||||||||||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Same here
Suggested change
|
||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
// Set up the objects to compute an integral on the reference cell | ||||||||||||||||||||||||||||||||||
FEValues<dim> fe_values(*fe, | ||||||||||||||||||||||||||||||||||
*q, | ||||||||||||||||||||||||||||||||||
update_JxW_values | update_quadrature_points); | ||||||||||||||||||||||||||||||||||
fe_values.reinit(triangulation.begin_active()); | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
double volume = 0; | ||||||||||||||||||||||||||||||||||
Point<dim> barycenter; | ||||||||||||||||||||||||||||||||||
for (unsigned int i = 0; i < q->size(); ++i) | ||||||||||||||||||||||||||||||||||
{ | ||||||||||||||||||||||||||||||||||
volume += fe_values.JxW(i); | ||||||||||||||||||||||||||||||||||
barycenter += fe_values.quadrature_point(i) * fe_values.JxW(i); | ||||||||||||||||||||||||||||||||||
} | ||||||||||||||||||||||||||||||||||
barycenter /= volume; | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
deallog << "ReferenceCell: " << reference_cell.to_string() << std::endl; | ||||||||||||||||||||||||||||||||||
deallog << " computed barycenter = " << barycenter << std::endl; | ||||||||||||||||||||||||||||||||||
deallog << " self-reported barycenter = " << reference_cell.barycenter<dim>() | ||||||||||||||||||||||||||||||||||
<< std::endl; | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
Assert((barycenter - reference_cell.barycenter<dim>()).norm() <= 1e-12, | ||||||||||||||||||||||||||||||||||
ExcInternalError()); | ||||||||||||||||||||||||||||||||||
} | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
int | ||||||||||||||||||||||||||||||||||
main() | ||||||||||||||||||||||||||||||||||
{ | ||||||||||||||||||||||||||||||||||
initlog(); | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
{ | ||||||||||||||||||||||||||||||||||
deallog.push("1D"); | ||||||||||||||||||||||||||||||||||
test<1>(ReferenceCells::Line); | ||||||||||||||||||||||||||||||||||
deallog.pop(); | ||||||||||||||||||||||||||||||||||
} | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
{ | ||||||||||||||||||||||||||||||||||
deallog.push("2D"); | ||||||||||||||||||||||||||||||||||
test<2>(ReferenceCells::Quadrilateral); | ||||||||||||||||||||||||||||||||||
test<2>(ReferenceCells::Triangle); | ||||||||||||||||||||||||||||||||||
deallog.pop(); | ||||||||||||||||||||||||||||||||||
} | ||||||||||||||||||||||||||||||||||
|
||||||||||||||||||||||||||||||||||
{ | ||||||||||||||||||||||||||||||||||
deallog.push("3D"); | ||||||||||||||||||||||||||||||||||
test<3>(ReferenceCells::Tetrahedron); | ||||||||||||||||||||||||||||||||||
test<3>(ReferenceCells::Pyramid); | ||||||||||||||||||||||||||||||||||
test<3>(ReferenceCells::Wedge); | ||||||||||||||||||||||||||||||||||
test<3>(ReferenceCells::Hexahedron); | ||||||||||||||||||||||||||||||||||
deallog.pop(); | ||||||||||||||||||||||||||||||||||
} | ||||||||||||||||||||||||||||||||||
} |
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,22 @@ | ||
|
||
DEAL:1D::ReferenceCell: Line | ||
DEAL:1D:: computed barycenter = 0.500000 | ||
DEAL:1D:: self-reported barycenter = 0.500000 | ||
DEAL:2D::ReferenceCell: Quad | ||
DEAL:2D:: computed barycenter = 0.500000 0.500000 | ||
DEAL:2D:: self-reported barycenter = 0.500000 0.500000 | ||
DEAL:2D::ReferenceCell: Tri | ||
DEAL:2D:: computed barycenter = 0.333333 0.333333 | ||
DEAL:2D:: self-reported barycenter = 0.333333 0.333333 | ||
DEAL:3D::ReferenceCell: Tet | ||
DEAL:3D:: computed barycenter = 0.250000 0.250000 0.250000 | ||
DEAL:3D:: self-reported barycenter = 0.250000 0.250000 0.250000 | ||
DEAL:3D::ReferenceCell: Pyramid | ||
DEAL:3D:: computed barycenter = 2.08167e-17 6.24500e-17 0.250000 | ||
bangerth marked this conversation as resolved.
Show resolved
Hide resolved
|
||
DEAL:3D:: self-reported barycenter = 0.00000 0.00000 0.250000 | ||
DEAL:3D::ReferenceCell: Wedge | ||
DEAL:3D:: computed barycenter = 0.333333 0.333333 0.500000 | ||
DEAL:3D:: self-reported barycenter = 0.333333 0.333333 0.500000 | ||
DEAL:3D::ReferenceCell: Hex | ||
DEAL:3D:: computed barycenter = 0.500000 0.500000 0.500000 | ||
DEAL:3D:: self-reported barycenter = 0.500000 0.500000 0.500000 |
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Also picky, but we could use
ReferenceCell::get_gauss_type_quadrature(2))
here for improved genericity.There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Would you allow me to do this as a follow-up for cloned tests? I think there is value in testing the equivalence of the output since what constitutes a "Gauss" formula is not unique for anything other than tensor-product formulas (and the function you propose may or may not return the formulas I use above).