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 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232
|
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
/**
* @class vtkCellGridSidesQuery
* @brief Compute external faces of a cell-grid.
*/
#ifndef vtkCellGridSidesQuery_h
#define vtkCellGridSidesQuery_h
#include "vtkCellGridQuery.h"
#include "vtkStringToken.h" // For API.
#include <functional>
#include <map>
#include <set>
#include <unordered_map>
#include <vector>
VTK_ABI_NAMESPACE_BEGIN
/**\brief A cell-grid query for enumerating sides of cells.
*
* As responders invoke the AddSides() method on this query,
* entries are added to storage indicating the number of times
* a given shape + connectivity appear. After the query is
* finalized, calling GetSides() will return a map holding
* sets of sides organized by cell ID for any shape + connectivity
* repeated an odd number of times.
*/
class VTKCOMMONDATAMODEL_EXPORT vtkCellGridSidesQuery : public vtkCellGridQuery
{
public:
static vtkCellGridSidesQuery* New();
vtkTypeMacro(vtkCellGridSidesQuery, vtkCellGridQuery);
void PrintSelf(ostream& os, vtkIndent indent) override;
void Initialize() override;
void Finalize() override;
std::map<vtkStringToken,
std::unordered_map<vtkStringToken, std::unordered_map<vtkIdType, std::set<int>>>>&
GetSides()
{
return this->Sides;
}
// std::map<vtkStringToken, std::unordered_map<vtkIdType, std::set<int>>>& GetSides() { return
// this->Sides; }
/// Add a \a side with the given \a shape and connectivity to the request's state.
///
/// The \a shape, \a conn size, and \a conn entries are hashed together into a key
/// which is mapped to a set of all the matching sides.
/// The \a cellType and \a cell ID are also stored with each matching side; these
/// are used during Finalize() to generate the output map-of-maps returned
/// by GetSides() so that the sides are reported by cell type, cell ID, and then side ID.
///
/// Note that the \a conn entries are hashed in a particular canonical order so
/// that the same hash is generated for sides with point IDs that have been shifted
/// and/or reversed.
/// The hash always starts at the smallest entry of \a conn and goes
/// in the direction that has the largest next entry.
/// Examples:
/// (3, 2, 0, 1) → starts at index 2 (0) and hashes backwards: (0, 2, 3, 1)
/// (4, 5, 6, 7) → starts at index 0 (4) and hashes backwards: (4, 7, 6, 5)
/// (7, 3, 6, 2) → starts at index 3 (2) and hashes forwards: (2, 7, 3, 6)
///
/// By storing the \a cellType, we avoid requiring a global-to-local cell numbering
/// in vtkCellGrid instances (as vtkPolyData incurs) which may hold multiple types of cells.
//@{
template <typename T, std::size_t NN>
void AddSide(vtkStringToken cellType, vtkIdType cell, int side, vtkStringToken shape,
const std::array<T, NN>& conn)
{
// Find where we should start hashing and in what direction.
std::size_t ss = 0;
T smin = conn[0];
for (std::size_t jj = 1; jj < NN; ++jj)
{
if (conn[jj] < smin)
{
smin = conn[jj];
ss = jj;
}
}
bool forward = conn[(ss + 1) % NN] > conn[(ss + NN - 1) % NN];
std::size_t hashedValue = std::hash<std::size_t>{}(NN);
vtkCellGridSidesQuery::HashCombiner()(hashedValue, shape.GetId());
if (forward)
{
for (std::size_t ii = 0; ii < NN; ++ii)
{
std::size_t hashedToken = std::hash<T>{}(conn[(ss + ii) % NN]);
vtkCellGridSidesQuery::HashCombiner()(hashedValue, hashedToken);
}
}
else // backward
{
for (std::size_t ii = 0; ii < NN; ++ii)
{
std::size_t hashedToken = std::hash<T>{}(conn[(ss + NN - ii) % NN]);
// hashedValue = hashedValue ^ (hashedToken << (ii + 1));
vtkCellGridSidesQuery::HashCombiner()(hashedValue, hashedToken);
}
}
this->Hashes[hashedValue].Sides.insert(Side{ cellType, shape, cell, side });
}
template <typename T>
void AddSide(vtkStringToken cellType, vtkIdType cell, int side, vtkStringToken shape,
const std::vector<T>& conn)
{
std::size_t ss = 0;
std::size_t NN = conn.size();
if (NN == 0)
{
return;
}
T smin = conn[0];
for (std::size_t jj = 1; jj < NN; ++jj)
{
if (conn[jj] < smin)
{
smin = conn[jj];
ss = jj;
}
}
bool forward = conn[(ss + 1) % NN] > conn[(ss + NN - 1) % NN];
std::size_t hashedValue = std::hash<std::size_t>{}(NN);
vtkCellGridSidesQuery::HashCombiner()(hashedValue, shape.GetId());
// std::cout << "Hash(" << (forward ? "F" : "R") << ")";
if (forward)
{
for (std::size_t ii = 0; ii < NN; ++ii)
{
std::size_t hashedToken = std::hash<T>{}(conn[(ss + ii) % NN]);
vtkCellGridSidesQuery::HashCombiner()(hashedValue, hashedToken);
// std::cout << " " << conn[(ss + ii) % NN];
}
}
else // backward
{
for (std::size_t ii = 0; ii < NN; ++ii)
{
std::size_t hashedToken = std::hash<T>{}(conn[(ss + NN - ii) % NN]);
// hashedValue = hashedValue ^ (hashedToken << (ii + 1));
vtkCellGridSidesQuery::HashCombiner()(hashedValue, hashedToken);
// std::cout << " " << conn[(ss + NN - ii) % NN];
}
}
// std::cout << " = " << std::hex << hashedValue << std::dec << "\n";
this->Hashes[hashedValue].Sides.insert(Side{ cellType, shape, cell, side });
}
//@}
protected:
vtkCellGridSidesQuery() = default;
~vtkCellGridSidesQuery() override = default;
// Hash combiner adapted from boost::hash_combine
class HashCombiner
{
public:
template <typename T>
void operator()(T& h, typename std::enable_if<sizeof(T) == 8, std::size_t>::type k)
{
constexpr T m = 0xc6a4a7935bd1e995ull;
constexpr int r = 47;
k *= m;
k ^= k >> r;
k *= m;
h ^= k;
h *= m;
// Completely arbitrary number, to prevent 0's
// from hashing to 0.
h += 0xe6546b64;
}
template <typename T>
void operator()(T& h, typename std::enable_if<sizeof(T) == 4, std::size_t>::type k)
{
constexpr std::uint32_t c1 = 0xcc9e2d51;
constexpr std::uint32_t c2 = 0x1b873593;
constexpr std::uint32_t r1 = 15;
constexpr std::uint32_t r2 = 13;
k *= c1;
k = (k << r1) | (k >> (32 - r1));
k *= c2;
h ^= k;
h = (h << r2) | (h >> (32 - r2));
h = h * 5 + 0xe6546b64;
}
};
struct Side
{
vtkStringToken CellType;
vtkStringToken SideShape;
vtkIdType DOF;
int SideId;
bool operator<(const Side& other) const
{
return (this->CellType < other.CellType) ||
(this->CellType == other.CellType &&
((this->DOF < other.DOF) || (this->DOF == other.DOF && this->SideId < other.SideId)));
}
};
struct Entry
{
std::set<Side> Sides;
};
std::unordered_map<std::size_t, Entry> Hashes;
std::map<vtkStringToken,
std::unordered_map<vtkStringToken, std::unordered_map<vtkIdType, std::set<int>>>>
Sides;
private:
vtkCellGridSidesQuery(const vtkCellGridSidesQuery&) = delete;
void operator=(const vtkCellGridSidesQuery&) = delete;
};
VTK_ABI_NAMESPACE_END
#endif // vtkCellGridSidesQuery_h
|