/usr/include/dolfin/generation/BoxMesh.h is in libdolfin-dev 2017.2.0.post0-2.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 | // Copyright (C) 2005-2017 Anders Logg and Garth N. Wells
//
// This file is part of DOLFIN.
//
// DOLFIN is free software: you can 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 3 of the License, or
// (at your option) any later version.
//
// DOLFIN is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
// GNU Lesser General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public License
// along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
#ifndef __BOX_H
#define __BOX_H
#include <array>
#include <cstddef>
#include <dolfin/log/log.h>
#include <dolfin/common/MPI.h>
#include <dolfin/mesh/Mesh.h>
namespace dolfin
{
/// Tetrahedral mesh of the 3D rectangular prism spanned by two
/// points p0 and p1. Given the number of cells (nx, ny, nz) in each
/// direction, the total number of tetrahedra will be 6*nx*ny*nz and
/// the total number of vertices will be (nx + 1)*(ny + 1)*(nz + 1).
class BoxMesh : public Mesh
{
public:
/// Create a uniform finite element _Mesh_ over the rectangular
/// prism spanned by the two _Point_s p0 and p1. The order of the
/// two points is not important in terms of minimum and maximum
/// coordinates.
///
/// @param p (std::array<_Point_, 2>)
/// Points of box.
/// @param n (std::array<double, 3> )
/// Number of cells in each direction.
/// @param cell_type
/// Tetrahedron or Hexahedron
///
/// @code{.cpp}
/// // Mesh with 8 cells in each direction on the
/// // set [-1,2] x [-1,2] x [-1,2].
/// Point p0(-1, -1, -1);
/// Point p1(2, 2, 2);
/// auto mesh = BoxMesh::create({p0, p1}, {8, 8, 8});
/// @endcode
static Mesh create(const std::array<Point,2 >& p, std::array<std::size_t, 3> n, CellType::Type cell_type)
{ return create(MPI_COMM_WORLD, p, n, cell_type); }
/// Create a uniform finite element _Mesh_ over the rectangular
/// prism spanned by the two _Point_s p0 and p1. The order of the
/// two points is not important in terms of minimum and maximum
/// coordinates.
///
/// @param comm (MPI_Comm)
/// MPI communicator
/// @param p (std::array<_Point_, 2>)
/// Points of box.
/// @param n (std::array<double, 3> )
/// Number of cells in each direction.
/// @param cell_type
/// Tetrahedron or hexahedron
///
/// @code{.cpp}
/// // Mesh with 8 cells in each direction on the
/// // set [-1,2] x [-1,2] x [-1,2].
/// Point p0(-1, -1, -1);
/// Point p1(2, 2, 2);
/// auto mesh = BoxMesh::create({p0, p1}, {8, 8, 8});
/// @endcode
static Mesh create(MPI_Comm comm, const std::array<Point, 2>& p,
std::array<std::size_t, 3> n, CellType::Type cell_type)
{
Mesh mesh(comm);
if (cell_type == CellType::Type::tetrahedron)
build_tet(mesh, p, n);
else if (cell_type == CellType::Type::hexahedron)
build_hex(mesh, n);
else
{
dolfin_error("BoxMesh.h",
"generate box mesh",
"Wrong cell type '%d'", cell_type);
}
return mesh;
}
/// Deprecated
/// Create a uniform finite element _Mesh_ over the rectangular
/// prism spanned by the two _Point_s p0 and p1. The order of the
/// two points is not important in terms of minimum and maximum
/// coordinates.
///
/// @param p0 (_Point_)
/// First point.
/// @param p1 (_Point_)
/// Second point.
/// @param nx (double)
/// Number of cells in x-direction.
/// @param ny (double)
/// Number of cells in y-direction.
/// @param nz (double)
/// Number of cells in z-direction.
///
/// @code{.cpp}
/// // Mesh with 8 cells in each direction on the
/// // set [-1,2] x [-1,2] x [-1,2].
/// Point p0(-1, -1, -1);
/// Point p1(2, 2, 2);
/// BoxMesh mesh(p0, p1, 8, 8, 8);
/// @endcode
BoxMesh(const Point& p0, const Point& p1,
std::size_t nx, std::size_t ny, std::size_t nz);
/// Deprecated
/// Create a uniform finite element _Mesh_ over the rectangular
/// prism spanned by the two _Point_s p0 and p1. The order of the
/// two points is not important in terms of minimum and maximum
/// coordinates.
///
/// @param comm (MPI_Comm)
/// MPI communicator
/// @param p0 (_Point_)
/// First point.
/// @param p1 (_Point_)
/// Second point.
/// @param nx (double)
/// Number of cells in x-direction.
/// @param ny (double)
/// Number of cells in y-direction.
/// @param nz (double)
/// Number of cells in z-direction.
///
/// @code{.cpp}
/// // Mesh with 8 cells in each direction on the
/// // set [-1,2] x [-1,2] x [-1,2].
/// Point p0(-1, -1, -1);
/// Point p1(2, 2, 2);
/// BoxMesh mesh(MPI_COMM_WORLD, p0, p1, 8, 8, 8);
/// @endcode
BoxMesh(MPI_Comm comm, const Point& p0, const Point& p1,
std::size_t nx, std::size_t ny, std::size_t nz);
private:
// Build mesh
static void build_tet(Mesh& mesh, const std::array<Point, 2>& p,
std::array<std::size_t, 3> n);
static void build_hex(Mesh& mesh, std::array<std::size_t, 3> n);
};
}
#endif
|