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 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187
|
/*=========================================================================
*
* Copyright NumFOCUS
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* https://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
#ifndef itkHexahedronCell_h
#define itkHexahedronCell_h
#include "itkQuadrilateralCell.h"
#include "itkHexahedronCellTopology.h"
#include "itkMakeFilled.h"
#include <array>
namespace itk
{
/** \class HexahedronCell
* \brief Represents a hexahedron (cuboid) for a Mesh.
*
* HexahedronCell represents a hexahedron, more precisely, a cuboid, for a Mesh.
*
* \tparam TPixelType The type associated with a point, cell, or boundary
* for use in storing its data.
*
* \tparam TCellTraits Type information of mesh containing cell.
*
* \todo When reviewing this class, the documentation of the template
* parameters MUST be fixed.
*
* NOTE: ONLY 3D implementations are instrumented. All other dimensions
* result in incorrect processing.
*
* \ingroup MeshObjects
* \ingroup ITKCommon
*/
template <typename TCellInterface>
class ITK_TEMPLATE_EXPORT HexahedronCell
: public TCellInterface
, private HexahedronCellTopology
{
public:
ITK_DISALLOW_COPY_AND_MOVE(HexahedronCell);
/** Standard class type aliases. */
itkCellCommonTypedefs(HexahedronCell);
itkCellInheritedTypedefs(TCellInterface);
/** \see LightObject::GetNameOfClass() */
itkOverrideGetNameOfClassMacro(HexahedronCell);
/** The type of boundary for this triangle's vertices. */
using VertexType = VertexCell<TCellInterface>;
using VertexAutoPointer = typename VertexType::SelfAutoPointer;
/** The type of boundary for this triangle's edges. */
using EdgeType = LineCell<TCellInterface>;
using EdgeAutoPointer = typename EdgeType::SelfAutoPointer;
/** The type of boundary for this hexahedron's faces. */
using FaceType = QuadrilateralCell<TCellInterface>;
using FaceAutoPointer = typename FaceType::SelfAutoPointer;
/** Hexahedron-specific topology numbers. */
static constexpr unsigned int NumberOfPoints = 8;
static constexpr unsigned int NumberOfVertices = 8;
static constexpr unsigned int NumberOfEdges = 12;
static constexpr unsigned int NumberOfFaces = 6;
static constexpr unsigned int CellDimension = 3;
/** HARDCODE Implementation requirements, while
* allowing general interface. The General interface
* is needed to facilitate the general SpatialObject
* loader.
*/
static constexpr unsigned int CellDimension3D = 3;
static constexpr unsigned int PointDimension3D = 3;
/** Implement the standard CellInterface. */
CellGeometryEnum
GetType() const override
{
return CellGeometryEnum::HEXAHEDRON_CELL;
}
void
MakeCopy(CellAutoPointer &) const override;
unsigned int
GetDimension() const override;
unsigned int
GetNumberOfPoints() const override;
CellFeatureCount
GetNumberOfBoundaryFeatures(int dimension) const override;
bool
GetBoundaryFeature(int dimension, CellFeatureIdentifier, CellAutoPointer &) override;
void
SetPointIds(PointIdConstIterator first) override;
void
SetPointIds(PointIdConstIterator first, PointIdConstIterator last) override;
void
SetPointId(int localId, PointIdentifier) override;
PointIdIterator
PointIdsBegin() override;
PointIdConstIterator
PointIdsBegin() const override;
PointIdIterator
PointIdsEnd() override;
PointIdConstIterator
PointIdsEnd() const override;
/** Hexahedron-specific interface. */
virtual CellFeatureCount
GetNumberOfVertices() const;
virtual CellFeatureCount
GetNumberOfEdges() const;
virtual CellFeatureCount
GetNumberOfFaces() const;
virtual bool
GetVertex(CellFeatureIdentifier, VertexAutoPointer &);
virtual bool
GetEdge(CellFeatureIdentifier, EdgeAutoPointer &);
virtual bool
GetFace(CellFeatureIdentifier, FaceAutoPointer &);
/** Evaluate the position inside the cell */
bool
EvaluatePosition(CoordRepType *,
PointsContainer *,
CoordRepType *,
CoordRepType[],
double *,
InterpolationWeightType *) override;
/** Visitor interface */
itkCellVisitMacro(CellGeometryEnum::HEXAHEDRON_CELL);
protected:
/** Store the number of points needed for a hexahedron. */
std::array<PointIdentifier, NumberOfPoints> m_PointIds{ MakeFilled<std::array<PointIdentifier, NumberOfPoints>>(
NumericTraits<PointIdentifier>::max()) };
void
InterpolationDerivs(CoordRepType pcoords[Self::CellDimension],
CoordRepType derivs[Self::CellDimension * Self::NumberOfPoints]);
void
InterpolationFunctions(CoordRepType pcoords[Self::CellDimension], InterpolationWeightType sf[Self::NumberOfPoints]);
void
EvaluateLocation(int & itkNotUsed(subId),
PointsContainer * points,
CoordRepType pcoords[Self::CellDimension],
CoordRepType x[Self::CellDimension],
InterpolationWeightType * weights);
public:
HexahedronCell() = default;
~HexahedronCell() override = default;
};
} // end namespace itk
#ifndef ITK_MANUAL_INSTANTIATION
# include "itkHexahedronCell.hxx"
#endif
#endif
|