tk namespace

Toolkit declarations and definitions for general purpose utilities.

Classes

class Around
Helper class simplifying client code for iterating on entries surrounding entries via linked lists derived from unstructured mesh connectivity.
class ASCMeshReader
ASCMeshReader : tk::Reader.
class ConjugateGradients
ConjugateGradients Charm++ chare array used to perform a distributed linear solve with the conjugate gradients algorithm.
class CSR
Compressed sparse row (CSR) storage for a sparse matrix.
template<uint8_t Layout>
class Data
Zero-runtime-cost data-layout wrappers with type-based compile-time dispatch.
template<class List>
struct DeepTuplePrinter
class DiagWriter
DiagWriter : tk::Writer.
class Exception
Basic exception class for producing file:func:line info + call trace.
class ExodusIIMeshReader
class ExodusIIMeshWriter
class GmshMeshReader
class GmshMeshWriter
template<typename, typename = std::void_t<>>
struct HasFunction_expect_choices
Detect if a type defines function 'expect::choices()'.
template<typename, typename = std::void_t<>>
struct HasFunction_expect_description
Detect if a type defines function 'expect::description()'.
template<typename, typename = std::void_t<>>
struct HasTypedef_alias
Detect if a type defines type 'alias'.
template<typename, typename = std::void_t<>>
struct HasTypedef_code
Detect if a type defines type 'code'.
template<typename, typename = std::void_t<>>
struct HasTypedef_i_am_tagged_tuple
Detect if a type defines type 'i_am_tagged_tuple'.
template<typename, typename = std::void_t<>>
struct HasVar_expect_lower
Detect if a type defines variable 'expect::lower'.
template<typename, typename = std::void_t<>>
struct HasVar_expect_upper
Detect if a type defines variable 'expect::upper'.
template<bool Condition, typename Then, typename Else = void>
struct if_
Type selection: if_< Condition, Then, Else >::type.
class LBSwitch
class MeditMeshReader
MeditMeshReader : tk::Reader.
class MeshWriter
Charm++ group used to output particle data to file in parallel.
class NetgenMeshReader
NetgenMeshReader : tk::Reader.
class NetgenMeshWriter
class PDFWriter
PDFWriter : Writer.
class Print
template<std::size_t N>
class Progress
class QuietCerr
class RDGFLOMeshReader
RDGFLOMeshReader : tk::Reader.
class Reader
template<class List>
class TaggedTuple
Tagged tuple, allowing tag-based access.
class Timer
Simple class to do timing.
template<class List>
struct TuplePrinter
class UGRIDMeshReader
UGRIDMeshReader : tk::Reader.
class UniPDF
Univariate PDF estimator.
class UnsMesh
3D unstructured mesh class
class Writer

Enums

enum ErrCode { SUCCESS = EXIT_SUCCESS, FAILURE = EXIT_FAILURE }
Error codes for the OS (or whatever calls us)
enum class ExoElemType: int { TET = 0, TRI = 1 }
enum class ExoWriter { CREATE, OPEN }
ExodusII writer constructor modes.
enum GmshElemType { TRI = 2, TET = 4, PNT = 15 }
Identifiers of supported Gmsh elements.
enum class GmshFileType { UNDEFINED = -1, ASCII = 0, BINARY = 1 }
Supported Gmsh mesh file types.
enum class MeshReaderType: uint8_t { GMSH = 0, NETGEN, EXODUSII, ASC, UGRID, RDGFLO, MEDIT }
Supported mesh readers.
enum class MeshWriterType: uint8_t { GMSH = 0, NETGEN, EXODUSII }
Supported mesh writers.
enum class HeaderType: uint8_t { INCITER, UNITTEST, MESHCONV }
Executable types for which an ascii logo is available in tk::Print.
enum class Centering: char { NODE = 'n', ELEM = 'e' }

Typedefs

template<std::size_t N>
using Table = std::vector<std::array<real, N+1>>
using real = double
Real number type used throughout the whole code.

Functions

template<class Container>
void unique(Container& c)
template<class Container>
auto uniquecopy(const Container& src) -> Container
template<typename Container>
auto cref_find(const Container& map, const typename Container::key_type& key) -> auto noexcept(…)
Find and return a constant reference to value for key in container that provides a find() member function with error handling.
template<typename Container>
auto ref_find(const Container& map, const typename Container::key_type& key) -> auto noexcept(…)
Find and return a reference to value for key in a container that provides a find() member function with error handling.
template<typename T>
auto extents(const std::vector<T>& vec) -> std::array<T, 2>
Return minimum and maximum values of a vector.
template<typename Container>
auto extents(const Container& map) -> std::array< typename Container::mapped_type, 2 > -> auto
Find and return minimum and maximum values in associative container.
template<class T, std::size_t N>
auto operator+=(std::array<T, N>& dst, const std::array<T, N>& src) -> std::array<T, N>&
template<class T, class Allocator>
auto operator+=(std::vector<T, Allocator>& dst, const std::vector<T, Allocator>& src) -> std::vector<T, Allocator>&
template<class T, class Allocator>
auto operator-=(std::vector<T, Allocator>& dst, const std::vector<T, Allocator>& src) -> std::vector<T, Allocator>&
template<class T, class Allocator>
auto operator/=(std::vector<T, Allocator>& dst, const std::vector<T, Allocator>& src) -> std::vector<T, Allocator>&
template<class C1, class C2>
auto keyEqual(const C1& a, const C2& b) -> bool
template<class Container>
auto sumsize(const Container& c) -> std::size_t
template<class Container>
auto numunique(const Container& c) -> std::size_t
template<class Map>
auto sumvalsize(const Map& c) -> std::size_t
template<class Container>
void destroy(Container& c)
template<typename Container, typename Predicate>
void erase_if(Container& items, const Predicate& predicate)
template<class T>
void concat(std::vector<T>&& src, std::vector<T>& dst)
template<class T>
void concat(std::vector<std::pair<bool, T>>&& src, std::vector<std::pair<bool, T>>& dst)
template<class Key, class Hash = std::hash<Key>, class Eq = std::equal_to<Key>>
void concat(std::unordered_set<Key, Hash, Eq>&& src, std::unordered_set<Key, Hash, Eq>& dst)
template<class Key, class Value>
auto operator<<(std::ostream& os, const std::pair<const Key, Value>& v) -> std::ostream&
template<typename T>
auto parameter(const T& v) -> std::string
Convert and return value as string.
template<typename V>
auto parameters(const V& v) -> std::string
Convert and return values from container as string.
template<uint8_t Layout>
auto operator*(tk::real lhs, const Data<Layout>& rhs) -> Data<Layout>
template<uint8_t Layout>
auto min(const Data<Layout>& a, const Data<Layout>& b) -> Data<Layout>
template<uint8_t Layout>
auto max(const Data<Layout>& a, const Data<Layout>& b) -> Data<Layout>
template<uint8_t Layout>
auto operator==(const Data<Layout>& lhs, const Data<Layout>& rhs) -> bool
template<uint8_t Layout>
auto operator!=(const Data<Layout>& lhs, const Data<Layout>& rhs) -> bool
template<uint8_t Layout>
auto maxdiff(const Data<Layout>& lhs, const Data<Layout>& rhs) -> std::pair<std::size_t, tk::real>
template<typename A, typename B>
auto flip_pair(const std::pair<A, B>& p) -> std::pair<B, A>
template<typename A, typename B>
auto flip_map(const std::map<A, B>& src) -> std::multimap<B, A>
auto linearLoadDistributor(real virtualization, uint64_t load, int npe, uint64_t& chunksize, uint64_t& remainder) -> uint64_t
Compute linear load distribution for given total work and virtualization.
template<class List>
auto operator<<(std::ostream& os, const tk::TaggedTuple<List>& t) -> std::ostream&
template<class List>
void print(std::ostream& os, const tk::TaggedTuple<List>& t)
template<class Tuple>
void print(std::ostream& os, const Tuple& c)
template<typename Enum, typename Ch, typename Tr, typename std::enable_if_t<std::is_enum_v<Enum>, int> = 0>
auto operator<<(std::basic_ostream<Ch, Tr>& os, const Enum& e) -> std::basic_ostream<Ch, Tr>&
template<class T, typename Ch, typename Tr>
auto operator<<(std::basic_ostream<Ch, Tr>& os, const std::vector<T>& v) -> std::basic_ostream<Ch, Tr>&
template<typename Ch, typename Tr, class Key, class Value, class Compare = std::less<Key>>
auto operator<<(std::basic_ostream<Ch, Tr>& os, const std::map<Key, Value, Compare>& m) -> std::basic_ostream<Ch, Tr>&
template<typename T, typename Ch, typename Tr>
auto operator<<(std::basic_string<Ch, Tr>& lhs, const T& e) -> std::basic_string<Ch, Tr>
template<typename T, typename Ch, typename Tr>
auto operator<<(std::basic_string<Ch, Tr>&& lhs, const T& e) -> std::basic_string<Ch, Tr>
void signalHandler(int signum)
Signal handler for multiple signals, SIGABRT, SIGSEGV, etc.
auto setSignalHandlers() -> int
Set signal handlers for multiple signals, SIGABRT, SIGSEGV, etc.
void processExceptionCharm()
Process an exception from the Charm++ runtime system.
template<std::size_t N>
auto sample(real x, const Table<N>& table) -> std::array<real, N>
auto hms(tk::real stamp) -> Timer::Watch
Convert existing time stamp as a real to Watch (global scope)
void cross(real v1x, real v1y, real v1z, real v2x, real v2y, real v2z, real& rx, real& ry, real& rz)
auto cross(const std::array<real, 3>& v1, const std::array<real, 3>& v2) -> std::array<real, 3>
void crossdiv(real v1x, real v1y, real v1z, real v2x, real v2y, real v2z, real j, real& rx, real& ry, real& rz)
auto crossdiv(const std::array<real, 3>& v1, const std::array<real, 3>& v2, real j) -> std::array<real, 3>
auto dot(real v1x, real v1y, real v1z, real v2x, real v2y, real v2z) -> real
auto dot(const std::array<real, 3>& v1, const std::array<real, 3>& v2) -> real
auto length(real x, real y, real z) -> real
auto length(const std::array<real, 3>& v) -> real
void unit(std::array<real, 3>& v) noexcept(…)
auto triple(real v1x, real v1y, real v1z, real v2x, real v2y, real v2z, real v3x, real v3y, real v3z) -> tk::real
auto triple(const std::array<real, 3>& v1, const std::array<real, 3>& v2, const std::array<real, 3>& v3) -> real
auto rotateX(const std::array<real, 3>& v, real angle) -> std::array<real, 3>
auto rotateY(const std::array<real, 3>& v, real angle) -> std::array<real, 3>
auto rotateZ(const std::array<real, 3>& v, real angle) -> std::array<real, 3>
auto normal(real x1, real x2, real x3, real y1, real y2, real y3, real z1, real z2, real z3, real& nx, real& ny, real& nz) -> real
auto serialize(const std::unordered_map<std::size_t, std::vector<std::size_t>>& d) -> std::pair<int, std::unique_ptr<char[]>>
Serialize to raw memory stream.
auto mergeGraph(int nmsg, CkReductionMsg** msgs) -> CkReductionMsg*
Charm++ custom reducer for merging during reduction across PEs.
auto serialize(const std::unordered_map<std::size_t, std::size_t>& d) -> std::pair<int, std::unique_ptr<char[]>>
Serialize to raw memory stream.
auto mergeParts(int nmsg, CkReductionMsg** msgs) -> CkReductionMsg*
Charm++ custom reducer for merging during reduction across PEs.
auto detectInput(const std::string& filename) -> MeshReaderType
Detect input mesh file type.
auto pickOutput(const std::string& filename) -> MeshWriterType
Determine output mesh file type.
auto readUnsMesh(const tk::Print& print, const std::string& filename, std::pair<std::string, tk::real>& timestamp) -> UnsMesh
Read unstructured mesh from file.
auto writeUnsMesh(const tk::Print& print, const std::string& filename, UnsMesh& mesh, bool reorder) -> std::vector<std::pair<std::string, tk::real>>
Write unstructured mesh to file.
static auto workdir() -> std::string
auto curtime() -> std::string
Wrapper for the standard C library's gettimeofday() from.
void echoHeader(HeaderType header)
Echo program header.
void echoBuildEnv(const std::string& executable)
Echo build environment.
void echoRunEnv(int argc, char** argv, int quiescence)
Echo runtime environment.
void finalize(const std::vector<tk::Timer>& timer, std::vector<std::pair<std::string, tk::Timer::Watch>>& timestamp, bool clean)
Finalize function for different executables.
auto npoin_in_graph(const std::vector<std::size_t>& inpoel) -> std::size_t
Compute number of points (nodes) in mesh from connectivity.
auto genEsup(const std::vector<std::size_t>& inpoel, std::size_t nnpe) -> std::pair<std::vector<std::size_t>, std::vector<std::size_t>>
Generate derived data structure, elements surrounding points.
auto genPsup(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup) -> std::pair<std::vector<std::size_t>, std::vector<std::size_t>>
Generate derived data structure, points surrounding points.
auto genEdsup(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup) -> std::pair<std::vector<std::size_t>, std::vector<std::size_t>>
Generate derived data structure, edges surrounding points.
auto genInpoed(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup) -> std::vector<std::size_t>
Generate derived data structure, edge connectivity.
auto genEsupel(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup) -> std::pair<std::vector<std::size_t>, std::vector<std::size_t>>
Generate derived data structure, elements surrounding points of elements.
auto genEsuel(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup) -> std::pair<std::vector<std::size_t>, std::vector<std::size_t>>
Generate derived data structure, elements surrounding elements.
auto genInedel(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::vector<std::size_t>& inpoed) -> std::vector<std::size_t>
Generate derived data structure, edges of elements.
auto genEsued(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup) -> std::unordered_map<UnsMesh::Edge, std::vector<std::size_t>, UnsMesh::Hash<2>, UnsMesh::Eq<2>>
Generate derived data structure, elements surrounding edges.
auto genEdpas(int mvecl, std::size_t nnpe, std::size_t npoin, const std::vector<std::size_t>& inpoed) -> std::pair<std::vector<std::size_t>, std::vector<std::size_t>>
Generate vector-groups for edges.
auto genNbfacTet(std::size_t tnbfac, const std::vector<std::size_t>& inpoel, const std::vector<std::size_t>& triinpoel_complete, const std::map<int, std::vector<std::size_t>>& bface_complete, const std::unordered_map<std::size_t, std::size_t>& lid, std::vector<std::size_t>& triinpoel, std::map<int, std::vector<std::size_t>>& bface) -> std::size_t
Generate number of boundary-faces and triangle boundary-face connectivity.
auto genEsuelTet(const std::vector<std::size_t>& inpoel, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup) -> std::vector<int>
Generate derived data structure, elements surrounding elements as a fixed length data structure as a full vector, including boundary elements as -1.
auto genNipfac(std::size_t nfpe, std::size_t nbfac, const std::vector<int>& esuelTet) -> std::size_t
Generate number of internal and physical-boundary faces.
auto genEsuf(std::size_t nfpe, std::size_t nipfac, std::size_t nbfac, const std::vector<std::size_t>& belem, const std::vector<int>& esuelTet) -> std::vector<int>
Generate derived data structure, elements surrounding faces.
auto genInpofaTet(std::size_t nipfac, std::size_t nbfac, const std::vector<std::size_t>& inpoel, const std::vector<std::size_t>& triinpoel, const std::vector<int>& esuelTet) -> std::vector<std::size_t>
Generate derived data structure, points on faces for tetrahedra only.
auto genBelemTet(std::size_t nbfac, const std::vector<std::size_t>& inpofa, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup) -> std::vector<std::size_t>
Generate derived data structure, boundary elements.
auto leakyPartition(const std::vector<int>& esueltet, const std::vector<std::size_t>& inpoel, const std::array<std::vector<real>, 3>& coord) -> bool
Perform leak-test on mesh (partition)
auto conforming(const std::vector<std::size_t>& inpoel, const std::array<std::vector<real>, 3>& coord, bool cerr, const std::vector<std::size_t>& rid) -> bool
Check if mesh (partition) is conforming.
auto intet(const std::array<std::vector<real>, 3>& coord, const std::vector<std::size_t>& inpoel, const std::vector<real>& p, std::size_t e, std::array<real, 4>& N) -> bool
Determine if a point is in a tetrahedron.
auto nodegrad(std::size_t node, const std::array<std::vector<tk::real>, 3>& coord, const std::vector<std::size_t>& inpoel, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup, const tk::Fields& U, uint64_t c) -> std::array<tk::real, 3>
Compute gradient at a mesh node.
auto edgegrad(const std::array<std::vector<tk::real>, 3>& coord, const std::vector<std::size_t>& inpoel, const std::vector<std::size_t>& esued, const tk::Fields& U, uint64_t c) -> std::array<tk::real, 3>
Compute gradient at a mesh edge.
auto shiftToZero(std::vector<std::size_t>& inpoel) -> std::size_t
Shift node IDs to start with zero in element connectivity.
void remap(std::vector<std::size_t>& ids, const std::vector<std::size_t>& map)
Apply new mapping to vector of indices.
void remap(std::vector<tk::real>& r, const std::vector<std::size_t>& map)
Apply new mapping to vector of real numbers.
auto remap(const std::vector<std::size_t>& ids, const std::vector<std::size_t>& map) -> std::vector<std::size_t>
Create remapped vector of indices using a vector.
void remap(std::vector<std::size_t>& ids, const std::unordered_map<std::size_t, std::size_t>& map)
In-place remap vector of indices using a map.
auto remap(const std::vector<std::size_t>& ids, const std::unordered_map<std::size_t, std::size_t>& map) -> std::vector<std::size_t>
Create remapped vector of indices using a map.
auto remap(const std::map<int, std::vector<std::size_t>>& ids, const std::unordered_map<std::size_t, std::size_t>& map) -> std::map<int, std::vector<std::size_t>>
Create remapped map of vector of indices using a map.
auto renumber(const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& psup) -> std::vector<std::size_t>
Reorder mesh points with the advancing front technique.
auto assignLid(const std::vector<std::size_t>& gid) -> std::unordered_map<std::size_t, std::size_t>
Assign local ids to global ids.
auto global2local(const std::vector<std::size_t>& ginpoel) -> std::tuple<std::vector<std::size_t>, std::vector<std::size_t>, std::unordered_map<std::size_t, std::size_t>>
Generate element connectivity of local node IDs from connectivity of global node IDs also returning the mapping between local to global IDs.
auto positiveJacobians(const std::vector<std::size_t>& inpoel, const std::array<std::vector<real>, 3>& coord) -> bool
Test positivity of the Jacobian for all cells in a mesh.
auto bfacenodes(const std::map<int, std::vector<std::size_t>>& bface, const std::vector<std::size_t>& triinpoel) -> std::map<int, std::vector<std::size_t>>
Generate nodes of side set faces.
auto count(const std::unordered_map<int, std::unordered_set<std::size_t>>& map, std::size_t node) -> tk::real
Count the number of contributions to a node.
auto slave(const std::unordered_map<int, std::unordered_set<std::size_t>>& map, std::size_t node, int chare) -> bool
Decide if a node is not counted by a chare.
auto serialize(std::size_t meshid, const std::vector<tk::UniPDF>& u) -> std::pair<int, std::unique_ptr<char[]>>
Serialize univariate PDF to raw memory stream.
auto mergeUniPDFs(int nmsg, CkReductionMsg** msgs) -> CkReductionMsg*
Charm++ custom reducer for merging a univariate PDF during reduction across PEs.
static auto operator<<(std::ostream& os, const tk::UniPDF& p) -> std::ostream&

Variables

const uint8_t UnkEqComp
Tags for selecting data layout policies.
const std::array<std::size_t, 2> ExoNnpe
const std::array<std::array<std::size_t, 3>, 4> expofa
const std::array<UnsMesh::Face, 4> lpofa
const std::array<UnsMesh::Edge, 6> lpoed
Const array defining the node ordering convention for tetrahedron edges.
const std::array<UnsMesh::Edge, 3> lpoet
Const array defining the node ordering convention for triangle edges.
static highwayhash::HH_U64 hh_key constexpr
static std::stringstream cerr_quiet
std::tringstream used to quiet std::cerr's stream by redirecting to it
static std::streambuf* cerr_old
std::streambuf used to store state of std::cerr before redirecting it

Enum documentation

enum tk::ErrCode

Error codes for the OS (or whatever calls us)

Enumerators
SUCCESS

Everything went fine.

FAILURE

Exceptions occurred.

enum class tk::ExoElemType: int

Supported ExodusII mesh cell types

enum class tk::ExoWriter

ExodusII writer constructor modes.

enum tk::GmshElemType

Identifiers of supported Gmsh elements.

enum class tk::GmshFileType

Supported Gmsh mesh file types.

enum class tk::MeshReaderType: uint8_t

Supported mesh readers.

enum class tk::MeshWriterType: uint8_t

Supported mesh writers.

enum class tk::HeaderType: uint8_t

Executable types for which an ascii logo is available in tk::Print.

enum class tk::Centering: char

Mesh/scheme centering types

Typedef documentation

template<std::size_t N>
using tk::Table = std::vector<std::array<real, N+1>>

Template parameters
N Number of ordinates in the table

Type alias for storing a discrete (y1,y2,...,yN) = f(x) function

using tk::real = double

Real number type used throughout the whole code.

Function documentation

template<class Container>
void tk::unique(Container& c)

Parameters
in/out Container

Make elements of container unique (in-place, overwriting source container)

template<class Container>
Container tk::uniquecopy(const Container& src)

Parameters
src in Container
Returns Container containing only unique elements compared to src

Make elements of container unique (on a copy, leaving the source as is)

template<typename Container>
auto tk::cref_find(const Container& map, const typename Container::key_type& key) noexcept(…)

Find and return a constant reference to value for key in container that provides a find() member function with error handling.

Parameters
map in Map associating values to keys
key in Key to search for
Returns A constant reference to the value associated to the key in map

template<typename Container>
auto tk::ref_find(const Container& map, const typename Container::key_type& key) noexcept(…)

Find and return a reference to value for key in a container that provides a find() member function with error handling.

Parameters
map in Map associating values to keys
key in Key to search for
Returns A reference to the value associated to the key in map

template<typename T>
std::array<T, 2> tk::extents(const std::vector<T>& vec)

Return minimum and maximum values of a vector.

Parameters
vec in Vector whose extents to compute
Returns Array of two values with the minimum and maximum values

template<typename Container>
auto tk::extents(const Container& map) -> std::array< typename Container::mapped_type, 2 >

Find and return minimum and maximum values in associative container.

Parameters
map in Map whose extents of values to find
Returns Array of two values with the minimum and maximum values in the map

template<class T, std::size_t N>
std::array<T, N>& tk::operator+=(std::array<T, N>& dst, const std::array<T, N>& src)

Parameters
dst in/out Destination array, i.e., left-hand side of a1 += a2
src in Source array, i.e., righ-hand side of a1 += a2
Returns Destination containing a1[0] += a2[0], a1[1] += a2[1], ...

Add all elements of an array to those of another one

template<class T, class Allocator>
std::vector<T, Allocator>& tk::operator+=(std::vector<T, Allocator>& dst, const std::vector<T, Allocator>& src)

Parameters
dst in/out Destination vector, i.e., left-hand side of v1 += v2
src in Source vector, i.e., righ-hand side of v1 += v2
Returns Destination containing v1[0] += v2[0], v1[1] += v2[1], ...

Add all elements of a vector to those of another one If src.size() > dst.size() will grow dst to that of src.size() padding with zeros.

template<class T, class Allocator>
std::vector<T, Allocator>& tk::operator-=(std::vector<T, Allocator>& dst, const std::vector<T, Allocator>& src)

Parameters
dst in/out Destination vector, i.e., left-hand side of v1 -= v2
src in Source vector, i.e., righ-hand side of v1 -= v2
Returns Destination containing v1[0] -= v2[0], v1[1] -= v2[1], ...

Subtract all elements of a vector from those of another one If src.size() > dst.size() will grow dst to that of src.size() padding with zeros.

template<class T, class Allocator>
std::vector<T, Allocator>& tk::operator/=(std::vector<T, Allocator>& dst, const std::vector<T, Allocator>& src)

Parameters
dst in/out Destination vector, i.e., left-hand side of v1 /= v2
src in Source vector, i.e., righ-hand side of v1 /= v2
Returns Destination containing v1[0] /= v2[0], v1[1] /= v2[1], ...

Divide all elements of a vector with those of another one If src.size() > dst.size() will grow dst to that of src.size() padding with zeros.

template<class C1, class C2>
bool tk::keyEqual(const C1& a, const C2& b)

Parameters
in 1st container to compare
in 2nd container to compare
Returns True if the containers have the same size and all keys (and only the keys) of the two containers are equal

Test if all keys of two associative containers are equal

template<class Container>
std::size_t tk::sumsize(const Container& c)

Parameters
in Container of containers
Returns Sum of the sizes of the containers of the container

Compute the sum of the sizes of a container of containers

template<class Container>
std::size_t tk::numunique(const Container& c)

Parameters
in Container of containers
Returns Number of unique values in a container of containers

Compute the number of unique values in a container of containers

template<class Map>
std::size_t tk::sumvalsize(const Map& c)

Template parameters
Map Container of containers type
Parameters
in Container of containers
Returns Sum of the sizes of the values of the associative container

Compute the sum of the sizes of the values of an associative container

template<class Container>
void tk::destroy(Container& c)

Parameters
in Container defining a swap() member function

Free memory of a container See http://stackoverflow.com/a/10465032 as to why this is done with the swap() member function of the container.

template<typename Container, typename Predicate>
void tk::erase_if(Container& items, const Predicate& predicate)

Template parameters
Container Type of container to remove from
Predicate Type for functor defining the predicate
Parameters
items Container object to remove from
predicate Predicate object instance to use

Remove items from container based on predicate

template<class T>
void tk::concat(std::vector<T>&& src, std::vector<T>& dst)

Template parameters
T Vector value type
Parameters
src in/out Source vector (moved from)
dst in/out Destination vector

Concatenate vectors of T

template<class T>
void tk::concat(std::vector<std::pair<bool, T>>&& src, std::vector<std::pair<bool, T>>& dst)

Template parameters
T Vector value type
Parameters
src in/out Source vector (moved from)
dst in/out Destination vector

Overwrite vectors of pair< bool, tk::real >

template<class Key, class Hash = std::hash<Key>, class Eq = std::equal_to<Key>>
void tk::concat(std::unordered_set<Key, Hash, Eq>&& src, std::unordered_set<Key, Hash, Eq>& dst)

Template parameters
Key Set key
Hash Set hasher
Eq Set equality operator
Parameters
src in/out Source set (moved from)
dst in/out Destination set

Concatenate unordered sets

template<class Key, class Value>
std::ostream& tk::operator<<(std::ostream& os, const std::pair<const Key, Value>& v)

Parameters
os in/out Output stream to write to
in value_type entry of a map
Returns Updated output stream

Operator << for writing value_type of a standard map to output streams

template<typename T>
std::string tk::parameter(const T& v)

Convert and return value as string.

Template parameters
T Value type for input
Parameters
in Value for input to return as a string
Returns String for input value

template<typename V>
std::string tk::parameters(const V& v)

Convert and return values from container as string.

Template parameters
V Container range for works on
Parameters
in Container whose components to return as a string
Returns Concatenated string of values read from a container

template<uint8_t Layout>
Data<Layout> tk::operator*(tk::real lhs, const Data<Layout>& rhs)

Parameters
lhs in Scalar to multiply with
rhs in Date object to multiply
Returns New Data object with all items multipled with lhs

Operator * multiplying all items by a scalar from the left

template<uint8_t Layout>
Data<Layout> tk::min(const Data<Layout>& a, const Data<Layout>& b)

Parameters
in 1st Data object
in 2nd Data object
Returns New Data object containing the minimum of all values for each value in a and b

Operator min between two Data objects

template<uint8_t Layout>
Data<Layout> tk::max(const Data<Layout>& a, const Data<Layout>& b)

Parameters
in 1st Data object
in 2nd Data object
Returns New Data object containing the maximum of all values for each value in a and b

Operator max between two Data objects

template<uint8_t Layout>
bool tk::operator==(const Data<Layout>& lhs, const Data<Layout>& rhs)

Parameters
lhs in Data object to compare
rhs in Data object to compare
Returns True if all entries are equal up to epsilon

Operator == between two Data objects

template<uint8_t Layout>
bool tk::operator!=(const Data<Layout>& lhs, const Data<Layout>& rhs)

Parameters
lhs in Data object to compare
rhs in Data object to compare
Returns True if all entries are unequal up to epsilon

Operator != between two Data objects

template<uint8_t Layout>
std::pair<std::size_t, tk::real> tk::maxdiff(const Data<Layout>& lhs, const Data<Layout>& rhs)

Parameters
lhs in 1st Data object
rhs in 2nd Data object
Returns The index, i.e., the raw position, of and the largest absolute value of the difference between all corresponding elements of lhs and rhs.

Compute the maximum difference between the elements of two Data objects The position returned is the position in the underlying raw data structure, independent of components. If lhs == rhs with precision std::numeric_limits< tk::real >::epsilon(), a pair of (0,0.0) is returned.

template<typename A, typename B>
std::pair<B, A> tk::flip_pair(const std::pair<A, B>& p)

Parameters
in std::pair of arbitrary types, A and B
Returns std::pair of arbitrary types, B and A

Flip a std::pair of arbitrary types

template<typename A, typename B>
std::multimap<B, A> tk::flip_map(const std::map<A, B>& src)

Parameters
src in std::map of arbitrary key and value pairs of types A and B
Returns std::multimap of arbitrary key and value pairs of types B and A

Flip a std::map of arbitrary types, yielding a std::multimap sorted by std::map::value_type.

uint64_t tk::linearLoadDistributor(real virtualization, uint64_t load, int npe, uint64_t& chunksize, uint64_t& remainder)

Compute linear load distribution for given total work and virtualization.

Parameters
virtualization in Degree of virtualization [0.0...1.0]
load in Total load, e.g., number of particles, number of mesh cells
npe in Number of processing elements to distribute the load to
chunksize in/out Chunk size, see detailed description
remainder in/out Remainder, see detailed description
Returns Number of work units

Compute load distibution (number of chares and chunksize) based on total work (e.g., total number of particles) and virtualization

The virtualization parameter, specified by the user, is a real number between 0.0 and 1.0, inclusive, which controls the degree of virtualization or over-decomposition. Independent of the value of virtualization the work is approximately evenly distributed among the available processing elements, given by npe. For zero virtualization (no over-decomposition), the work is simply decomposed into total_work/numPEs, which yields the smallest number of Charm++ chares and the largest chunks of work units. The other extreme is unity virtualization, which decomposes the total work into the smallest size work units possible, yielding the largest number of Charm++ chares. Obviously, the optimum will be between 0.0 and 1.0, depending on the problem.

The formula implemented uses a linear relationship between the virtualization parameter and the number of work units with the extremes described above. The formula is given by

chunksize = (1 - n) * v + n;

where

  • v = degree of virtualization
  • n = load/npes
  • load = total work, e.g., number of particles, number of mesh cells
  • npes = number of hardware processing elements

template<class List>
std::ostream& tk::operator<<(std::ostream& os, const tk::TaggedTuple<List>& t)

Template parameters
List brigand::list of types in the tagged tuple
Parameters
os in/out Output stream to output to
in TaggedTuple to print
Returns Output stream

Simple (unformatted, one-line) output of a TaggedTuple to output stream

template<class List>
void tk::print(std::ostream& os, const tk::TaggedTuple<List>& t)

Template parameters
List brigand::list of types in the tagged tuple
Parameters
os in/out Output stream to output to
in TaggedTuple to print

Simple (unformatted, one-line) output of a TaggedTuple to output stream

template<class Tuple>
void tk::print(std::ostream& os, const Tuple& c)

Template parameters
Tuple Tuple object type
Parameters
os in/out Output stream to print to
in Command line object to output to file

Output command line object (a TaggedTuple) to file

template<typename Enum, typename Ch, typename Tr, typename std::enable_if_t<std::is_enum_v<Enum>, int> = 0>
std::basic_ostream<Ch, Tr>& tk::operator<<(std::basic_ostream<Ch, Tr>& os, const Enum& e)

Parameters
os in Output stream into to write to
in Value of enum-class type to write to stream
Returns Updated output stream for chain-use of the operator

Operator << for writing an enum class to an output stream

template<class T, typename Ch, typename Tr>
std::basic_ostream<Ch, Tr>& tk::operator<<(std::basic_ostream<Ch, Tr>& os, const std::vector<T>& v)

Parameters
os in Output stream to write to
in Vector to write to stream
Returns Updated output stream for chain-use of the operator

Operator << for writing a std::vector to an output stream

template<typename Ch, typename Tr, class Key, class Value, class Compare = std::less<Key>>
std::basic_ostream<Ch, Tr>& tk::operator<<(std::basic_ostream<Ch, Tr>& os, const std::map<Key, Value, Compare>& m)

Parameters
os in Output stream to write to
in Map to write to stream
Returns Updated output stream for chain-use of the operator

Operator << for writing an std::map to an output stream

template<typename T, typename Ch, typename Tr>
std::basic_string<Ch, Tr> tk::operator<<(std::basic_string<Ch, Tr>& lhs, const T& e)

Parameters
lhs in Output std::basic_string into which e is written
in Value of arbitrary type to write to string
Returns Updated string

Operator << for adding (concatenating) T to a std::basic_string for lvalues.

template<typename T, typename Ch, typename Tr>
std::basic_string<Ch, Tr> tk::operator<<(std::basic_string<Ch, Tr>&& lhs, const T& e)

Parameters
lhs in Output std::basic_string into which e is written
in Value of arbitrary type to write to string
Returns Updated string

Operator << for adding (concatenating) T to a std::basic_string for rvalues.

void tk::signalHandler(int signum)

Signal handler for multiple signals, SIGABRT, SIGSEGV, etc.

Parameters
signum in Signal number

Signals caught: SIGABRT is generated when the program calls the abort() function, such as when an assert() triggers SIGSEGV is generated when the program makes an illegal memory access, such as reading unaligned memory, dereferencing a null pointer, reading memory out of bounds etc. SIGILL is generated when the program tries to execute a malformed instruction. This happens when the execution pointer starts reading non-program data, or when a pointer to a function is corrupted. SIGFPE is generated when executing an illegal floating point instruction, most commonly division by zero or floating point overflow.

int tk::setSignalHandlers()

Set signal handlers for multiple signals, SIGABRT, SIGSEGV, etc.

Returns Ignore, used for calling in a constructor's initializer list

void tk::processExceptionCharm()

Process an exception from the Charm++ runtime system.

See Josuttis, The C++ Standard Library - A Tutorial and Reference, 2nd Edition, 2012.

template<std::size_t N>
std::array<real, N> tk::sample(real x, const Table<N>& table)

Template parameters
N Number of ordinates in the table
Parameters
in Abscissa to sample at
table in Table to sample
Returns Ordinates sampled

Sample a discrete (y1,y2,...,yN) = f(x) function at x

Timer::Watch tk::hms(tk::real stamp)

Convert existing time stamp as a real to Watch (global scope)

Parameters
stamp in Time stamp as a real number
Returns Time as hours, minutes, and seconds, as a Watch struct.

Convert existing time stamp as a real to Watch (global-scope)

void tk::cross(real v1x, real v1y, real v1z, real v2x, real v2y, real v2z, real& rx, real& ry, real& rz)

Parameters
v1x in x coordinate of the 1st vector
v1y in y coordinate of the 1st vector
v1z in z coordinate of the 1st vector
v2x in x coordinate of the 2nd vector
v2y in y coordinate of the 2nd vector
v2z in z coordinate of the 2nd vector
rx out x coordinate of the product vector
ry out y coordinate of the product vector
rz out z coordinate of the product vector

Compute the cross-product of two vectors

std::array<real, 3> tk::cross(const std::array<real, 3>& v1, const std::array<real, 3>& v2)

Parameters
v1 in 1st vector
v2 in 2nd vector
Returns Cross-product

Compute the cross-product of two vectors

void tk::crossdiv(real v1x, real v1y, real v1z, real v2x, real v2y, real v2z, real j, real& rx, real& ry, real& rz)

Parameters
v1x in x coordinate of the 1st vector
v1y in y coordinate of the 1st vector
v1z in z coordinate of the 1st vector
v2x in x coordinate of the 2nd vector
v2y in y coordinate of the 2nd vector
v2z in z coordinate of the 2nd vector
in The scalar to divide the product with
rx out x coordinate of the product vector
ry out y coordinate of the product vector
rz out z coordinate of the product vector

Compute the cross-product of two vectors divided by a scalar

std::array<real, 3> tk::crossdiv(const std::array<real, 3>& v1, const std::array<real, 3>& v2, real j)

Parameters
v1 in 1st vector
v2 in 2nd vector
in Scalar to divide each component by
Returns Cross-product divided by scalar

Compute the cross-product of two vectors divided by a scalar

real tk::dot(real v1x, real v1y, real v1z, real v2x, real v2y, real v2z)

Parameters
v1x in x coordinate of 1st vector
v1y in y coordinate of 1st vector
v1z in z coordinate of 1st vector
v2x in x coordinate of 2nd vector
v2y in y coordinate of 2nd vector
v2z in z coordinate of 2ndt vector
Returns Dot-product

Compute the dot-product of two vectors

real tk::dot(const std::array<real, 3>& v1, const std::array<real, 3>& v2)

Parameters
v1 in 1st vector
v2 in 2nd vector
Returns Dot-product

Compute the dot-product of two vectors

real tk::length(real x, real y, real z)

Parameters
in X coordinate of vector
in Y coordinate of vector
in Z coordinate of vector
Returns length

Compute length of a vector

real tk::length(const std::array<real, 3>& v)

Parameters
in vector
Returns length

Compute length of a vector

void tk::unit(std::array<real, 3>& v) noexcept(…)

Parameters
in/out Vector to normalize

Scale vector to unit length

tk::real tk::triple(real v1x, real v1y, real v1z, real v2x, real v2y, real v2z, real v3x, real v3y, real v3z)

Parameters
v1x in x coordinate of the 1st vector
v1y in y coordinate of the 1st vector
v1z in z coordinate of the 1st vector
v2x in x coordinate of the 2nd vector
v2y in y coordinate of the 2nd vector
v2z in z coordinate of the 2nd vector
v3x in x coordinate of the 3rd vector
v3y in y coordinate of the 3rd vector
v3z in z coordinate of the 3rd vector
Returns Scalar value of the triple product

Compute the triple-product of three vectors

real tk::triple(const std::array<real, 3>& v1, const std::array<real, 3>& v2, const std::array<real, 3>& v3)

Parameters
v1 in 1st vector
v2 in 2nd vector
v3 in 3rd vector
Returns Triple-product

Compute the triple-product of three vectors

std::array<real, 3> tk::rotateX(const std::array<real, 3>& v, real angle)

Parameters
in Vector to rotate
angle in Angle to use to rotate with
Returns Rotated vector

Rotate vector about X axis

std::array<real, 3> tk::rotateY(const std::array<real, 3>& v, real angle)

Parameters
in Vector to rotate
angle in Angle to use to rotate with
Returns Rotated vector

Rotate vector about Y axis

std::array<real, 3> tk::rotateZ(const std::array<real, 3>& v, real angle)

Parameters
in Vector to rotate
angle in Angle to use to rotate with
Returns Rotated vector

Rotate vector about Z axis

real tk::normal(real x1, real x2, real x3, real y1, real y2, real y3, real z1, real z2, real z3, real& nx, real& ny, real& nz)

Parameters
x1 in x coordinate of the 1st vertex of the triangle
x2 in x coordinate of the 2nd vertex of the triangle
x3 in x coordinate of the 3rd vertex of the triangle
y1 in y coordinate of the 1st vertex of the triangle
y2 in y coordinate of the 2nd vertex of the triangle
y3 in y coordinate of the 3rd vertex of the triangle
z1 in z coordinate of the 1st vertex of the triangle
z2 in z coordinate of the 2nd vertex of the triangle
z3 in z coordinate of the 3rd vertex of the triangle
nx out x coordinate of the unit normal
ny out y coordinate of the unit normal
nz out z coordinate of the unit normal
Returns Triangle area

Compute the unit normal vector of a triangle

std::pair<int, std::unique_ptr<char[]>> tk::serialize(const std::unordered_map<std::size_t, std::vector<std::size_t>>& d)

Serialize to raw memory stream.

Parameters
in Mesh graph to aggregate
Returns Pair of the length and the raw stream containing the serialized data

CkReductionMsg* tk::mergeGraph(int nmsg, CkReductionMsg** msgs)

Charm++ custom reducer for merging during reduction across PEs.

Parameters
nmsg in Number of messages in msgs
msgs in Charm++ reduction message containing serialized data
Returns Aggregated diagnostics built for further aggregation if needed

std::pair<int, std::unique_ptr<char[]>> tk::serialize(const std::unordered_map<std::size_t, std::size_t>& d)

Serialize to raw memory stream.

Parameters
in Mesh graph to aggregate
Returns Pair of the length and the raw stream containing the serialized data

CkReductionMsg* tk::mergeParts(int nmsg, CkReductionMsg** msgs)

Charm++ custom reducer for merging during reduction across PEs.

Parameters
nmsg in Number of messages in msgs
msgs in Charm++ reduction message containing serialized data
Returns Aggregated diagnostics built for further aggregation if needed

MeshReaderType tk::detectInput(const std::string& filename)

Detect input mesh file type.

Parameters
filename in File to open and detect its type
Returns enum specifying the mesh reader type

MeshWriterType tk::pickOutput(const std::string& filename)

Determine output mesh file type.

Parameters
filename in Filename to pick its type based on extension given
Returns enum specifying the mesh writer type

UnsMesh tk::readUnsMesh(const tk::Print& print, const std::string& filename, std::pair<std::string, tk::real>& timestamp)

Read unstructured mesh from file.

Parameters
print in Pretty printer
filename in Filename to read mesh from
timestamp out A time stamp consisting of a timer label (a string), and a time state (a tk::real in seconds) measuring the mesh read time
Returns Unstructured mesh object

Create unstructured mesh to store mesh

std::vector<std::pair<std::string, tk::real>> tk::writeUnsMesh(const tk::Print& print, const std::string& filename, UnsMesh& mesh, bool reorder)

Write unstructured mesh to file.

Parameters
print in Pretty printer
filename in Filename to write mesh to
mesh in Unstructured mesh object to write from
reorder in Whether to also reorder mesh nodes
Returns Vector of time stamps consisting of a timer label (a string), and a time state (a tk::real in seconds) measuring the renumber and the mesh write time

static std::string tk::workdir()

Returns A stirng containing the current working directory

std::string tk::curtime()

Wrapper for the standard C library's gettimeofday() from.

Returns A stirng containing the current date and time

void tk::echoHeader(HeaderType header)

Echo program header.

Parameters
header in Header type enum indicating which header to print

void tk::echoBuildEnv(const std::string& executable)

Echo build environment.

Parameters
executable in Name of the executable

Echo information read from build_dir/Base/Config.h filled by CMake based on src/Main/Config.h.in.

void tk::echoRunEnv(int argc, char** argv, int quiescence)

Echo runtime environment.

Parameters
argc in Number of command-line arguments to executable
argv in C-style string array to command-line arguments to executable
quiescence in True if quiescence detection is enabled

void tk::finalize(const std::vector<tk::Timer>& timer, std::vector<std::pair<std::string, tk::Timer::Watch>>& timestamp, bool clean)

Finalize function for different executables.

Parameters
timer in Vector of timers, held by the main chare
timestamp in/out Vector of time stamps in h:m:s with labels
clean in True if we should exit with a zero exit code, false to exit with a nonzero exit code

std::size_t tk::npoin_in_graph(const std::vector<std::size_t>& inpoel)

Compute number of points (nodes) in mesh from connectivity.

Parameters
inpoel in Inteconnectivity of points and elements. These are the
Returns Number of mesh points (nodes)

std::pair<std::vector<std::size_t>, std::vector<std::size_t>> tk::genEsup(const std::vector<std::size_t>& inpoel, std::size_t nnpe)

Generate derived data structure, elements surrounding points.

Parameters
inpoel in

Inteconnectivity of points and elements. These are the node ids of each element of an unstructured mesh. Example:

std::vector< std::size_t > inpoel { 12, 14,  9, 11,
                                    10, 14, 13, 12 };

specifies two tetrahedra whose vertices (node ids) are { 12, 14, 9, 11 }, and { 10, 14, 13, 12 }.

nnpe in Number of nodes per element
Returns Linked lists storing elements surrounding points

The data generated here is stored in a linked list, more precisely, two linked arrays (vectors), esup1 and esup2, where esup2 holds the indices at which esup1 holds the element ids surrounding points. Looping over all elements surrounding all points can then be accomplished by the following loop:

for (std::size_t p=0; p<npoin; ++p)
  for (auto i=esup.second[p]+1; i<=esup.second[p+1]; ++i)
     use element id esup.first[i]

To find out the number of points, npoin, the mesh connectivity, inpoel, can be queried:

auto minmax = std::minmax_element( begin(inpoel), end(inpoel) );
Assert( *minmax.first == 0, "node ids should start from zero" );
auto npoin = *minmax.second + 1;

std::pair<std::vector<std::size_t>, std::vector<std::size_t>> tk::genPsup(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup)

Generate derived data structure, points surrounding points.

Parameters
inpoel in

Inteconnectivity of points and elements. These are the node ids of each element of an unstructured mesh. Example:

std::vector< std::size_t > inpoel { 12, 14,  9, 11,
                                    10, 14, 13, 12 };

specifies two tetrahedra whose vertices (node ids) are { 12, 14, 9, 11 }, and { 10, 14, 13, 12 }.

nnpe in Number of nodes per element
esup in Elements surrounding points as linked lists, see tk::genEsup
Returns Linked lists storing points surrounding points

The data generated here is stored in a linked list, more precisely, two linked arrays (vectors), psup1 and psup2, where psup2 holds the indices at which psup1 holds the point ids surrounding points. Looping over all points surrounding all points can then be accomplished by the following loop:

for (std::size_t p=0; p<npoin; ++p)
  for (auto i=psup.second[p]+1; i<=psup.second[p+1]; ++i)
     use point id psup.first[i]

To find out the number of points, npoin, the mesh connectivity, inpoel, can be queried:

auto minmax = std::minmax_element( begin(inpoel), end(inpoel) );
Assert( *minmax.first == 0, "node ids should start from zero" );
auto npoin = *minmax.second + 1;

or the length-1 of the generated index list: auto npoin = psup.second.size()-1;

std::pair<std::vector<std::size_t>, std::vector<std::size_t>> tk::genEdsup(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup)

Generate derived data structure, edges surrounding points.

Parameters
inpoel in

Inteconnectivity of points and elements. These are the node ids of each element of an unstructured mesh. Example:

std::vector< std::size_t > inpoel { 12, 14,  9, 11,
                                    10, 14, 13, 12 };

specifies two tetrahedra whose vertices (node ids) are { 12, 14, 9, 11 }, and { 10, 14, 13, 12 }.

nnpe in Number of nodes per element (3 or 4)
esup in Elements surrounding points as linked lists, see tk::genEsup
Returns Linked lists storing edges (point ids p < q) emanating from points

The data generated here is stored in a linked list, more precisely, two linked arrays (vectors), edsup1 and edsup2, where edsup2 holds the indices at which edsup1 holds the edge-end point ids emanating from points for all points. The generated data structure, linked lists edsup1 and edsup2, are very similar to psup1 and psup2, generated by genPsup(), except here only unique edges are stored, i.e., for edges with point ids p < q, only ids q are stored that are still associated to point p. Looping over all unique edges can then be accomplished by the following loop:

for (std::size_t p=0; p<npoin; ++p)
  for (auto i=edsup.second[p]+1; i<=edsup.second[p+1]; ++i)
    use edge with point ids p < edsup.first[i]

To find out the number of points, npoin, the mesh connectivity, inpoel, can be queried:

auto minmax = std::minmax_element( begin(inpoel), end(inpoel) );
Assert( *minmax.first == 0, "node ids should start from zero" );
auto npoin = *minmax.second + 1;

std::vector<std::size_t> tk::genInpoed(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup)

Generate derived data structure, edge connectivity.

Parameters
inpoel in

Inteconnectivity of points and elements. These are the node ids of each element of an unstructured mesh. Example:

std::vector< std::size_t > inpoel { 12, 14,  9, 11,
                                    10, 14, 13, 12 };

specifies two tetrahedra whose vertices (node ids) are { 12, 14, 9, 11 }, and { 10, 14, 13, 12 }.

nnpe in Number of nodes per element (3 or 4)
esup in Elements surrounding points as linked lists, see tk::genEsup
Returns Linear vector storing edge connectivity (point ids p < q)

The data generated here is stored in a linear vector and is very similar to the linked lists, edsup1 and _edsup2, generated by genEdsup(). The difference is that in the linear vector, inpoed, generated here, both edge point ids are stored as a pair, p < q, as opposed to the linked lists edsup1 and edsup2, in which edsup1 only stores the edge-end point ids (still associated to edge-start point ids when used together with edsup2). The rationale is that while inpoed is larger in memory, it allows direct access to edges (pair of point ids making up an edge), edsup1 and edsup2 are smaller in memory, still allow accessing the same data (edge point id pairs) but only in a linear fashion, not by direct access to particular edges. Accessing all unique edges using the edge connectivity data structure, inpoed, generated here can be accomplished by

for (std::size_t e=0; e<inpoed.size()/2; ++e) {
  use point id p of edge e = inpoed[e*2];
  use point id q of edge e = inpoed[e*2+1];
}

std::pair<std::vector<std::size_t>, std::vector<std::size_t>> tk::genEsupel(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup)

Generate derived data structure, elements surrounding points of elements.

Parameters
inpoel in

Inteconnectivity of points and elements. These are the node ids of each element of an unstructured mesh. Example:

std::vector< std::size_t > inpoel { 12, 14,  9, 11,
                                    10, 14, 13, 12 };

specifies two tetrahedra whose vertices (node ids) are { 12, 14, 9, 11 }, and { 10, 14, 13, 12 }.

nnpe in Number of nodes per element
esup in Elements surrounding points as linked lists, see tk::genEsup
Returns Linked lists storing elements surrounding points of elements

The data generated here is stored in a linked list, more precisely, two linked arrays (vectors), esupel1 and esupel2, where esupel2 holds the indices at which esupel1 holds the element ids surrounding points of elements. Looping over all elements surrounding the points of all elements can then be accomplished by the following loop:

for (std::size_t e=0; e<nelem; ++e)
  for (auto i=esupel.second[e]+1; i<=esupel.second[e+1]; ++i)
     use element id esupel.first[i]

To find out the number of elements, nelem, the size of the mesh connectivity vector, inpoel, can be devided by the number of nodes per elements, nnpe: auto nelem = inpoel.size()/nnpe;

std::pair<std::vector<std::size_t>, std::vector<std::size_t>> tk::genEsuel(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup)

Generate derived data structure, elements surrounding elements.

Parameters
inpoel in

Inteconnectivity of points and elements. These are the node ids of each element of an unstructured mesh. Example:

std::vector< std::size_t > inpoel { 12, 14,  9, 11,
                                    10, 14, 13, 12 };

specifies two tetrahedra whose vertices (node ids) are { 12, 14, 9, 11 }, and { 10, 14, 13, 12 }.

nnpe in Number of nodes per element
esup in Elements surrounding points as linked lists, see tk::genEsup
Returns Linked lists storing elements surrounding elements

The data generated here is stored in a linked list, more precisely, two linked arrays (vectors), esuel1 and esuel2, where esuel2 holds the indices at which esuel1 holds the element ids surrounding elements. Looping over elements surrounding elements can then be accomplished by the following loop:

for (std::size_t e=0; e<nelem; ++e)
  for (auto i=esuel.second[e]+1; i<=esuel.second[e+1]; ++i)
     use element id esuel.first[i]

To find out the number of elements, nelem, the size of the mesh connectivity vector, inpoel, can be devided by the number of nodes per elements, nnpe: auto nelem = inpoel.size()/nnpe;

std::vector<std::size_t> tk::genInedel(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::vector<std::size_t>& inpoed)

Generate derived data structure, edges of elements.

Parameters
inpoel in

Inteconnectivity of points and elements. These are the node ids of each element of an unstructured mesh. Example:

std::vector< std::size_t > inpoel { 12, 14,  9, 11,
                                    10, 14, 13, 12 };

specifies two tetrahedra whose vertices (node ids) are { 12, 14, 9, 11 }, and { 10, 14, 13, 12 }.

nnpe in Number of nodes per element
inpoed in Edge connectivity as linear vector, see tk::genInpoed
Returns Linear vector storing all edge ids * 2 of all elements

The data generated here is stored in a linear vector with all edge ids (as defined by inpoed) of all elements. The edge ids stored in inedel can be directly used to index the vector inpoed. Because the derived data structure generated here, inedel, is intended to be used in conjunction with the linear vector inpoed and not with the linked lists edsup1 and edsup2, this function takes inpoed as an argument. Accessing the edges of element e using the edge of elements data structure, inedel, generated here can be accomplished by

for (std::size_t e=0; e<nelem; ++e) {
  for (std::size_t i=0; i<nepe; ++i) {
    use edge id inedel[e*nepe+i] of element e, or
    use point ids p < q of edge id inedel[e*nepe+i] of element e as
      p = inpoed[ inedel[e*nepe+i]*2 ]
      q = inpoed[ inedel[e*nepe+i]*2+1 ]
  }
}

where nepe denotes the number of edges per elements: 3 for triangles, 6 for tetrahedra. To find out the number of elements, nelem, the size of the mesh connectivity vector, inpoel, can be devided by the number of nodes per elements, nnpe: auto nelem = inpoel.size()/nnpe;

std::unordered_map<UnsMesh::Edge, std::vector<std::size_t>, UnsMesh::Hash<2>, UnsMesh::Eq<2>> tk::genEsued(const std::vector<std::size_t>& inpoel, std::size_t nnpe, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup)

Generate derived data structure, elements surrounding edges.

Parameters
inpoel in

Inteconnectivity of points and elements. These are the node ids of each element of an unstructured mesh. Example:

std::vector< std::size_t > inpoel { 12, 14,  9, 11,
                                    10, 14, 13, 12 };

specifies two tetrahedra whose vertices (node ids) are { 12, 14, 9, 11 }, and { 10, 14, 13, 12 }.

nnpe in Number of nodes per element (3 or 4)
esup in Elements surrounding points as linked lists, see tk::genEsup
Returns Associative container storing elements surrounding edges (value), assigned to edge-end points (key)

Looping over elements surrounding all edges can be accomplished by the following loop:

for (const auto& [edge,surr_elements] : esued) {
  use element edge-end-point ids edge[0] and edge[1]
  for (auto e : surr_elements) {
     use element id e
  }
}

esued.size() equals the number of edges.

std::pair<std::vector<std::size_t>, std::vector<std::size_t>> tk::genEdpas(int mvecl, std::size_t nnpe, std::size_t npoin, const std::vector<std::size_t>& inpoed)

Generate vector-groups for edges.

Parameters
mvecl in Max vector length to target
nnpe in Number of nodes per (super-)edge
npoin in Number mesh points
inpoed in Edge connectivity as linear vector, see tk::genInpoed for nnpe=2
Returns Linked lists storing edge-groups so that any point of a group is accessed only once within a group.

The data generated here is stored in a linked list, more precisely, two linked arrays (vectors), edpas1 and edpas2, where edpas2 holds the indices at which edpas1 holds the edge ids of a vector group. Looping over all groups can then be accomplished by the following loop:

for (std::size_t w=0; w<edpas.second.size()-1; ++w)
  for (auto i=edpas.second[w]+1; i<=edpas.second[w+1]; ++i)
     use edge id edpas.first[i]

std::size_t tk::genNbfacTet(std::size_t tnbfac, const std::vector<std::size_t>& inpoel, const std::vector<std::size_t>& triinpoel_complete, const std::map<int, std::vector<std::size_t>>& bface_complete, const std::unordered_map<std::size_t, std::size_t>& lid, std::vector<std::size_t>& triinpoel, std::map<int, std::vector<std::size_t>>& bface)

Generate number of boundary-faces and triangle boundary-face connectivity.

Parameters
tnbfac in Total number of boundary faces in the entire mesh.
inpoel in Inteconnectivity of points and elements. These are the node ids of each element of an unstructured mesh.
triinpoel_complete in Interconnectivity of points and boundary-face in the entire mesh.
bface_complete in Map of boundary-face lists mapped to corresponding side set ids for the entire mesh.
lid in Mapping between the node indices used in the smaller inpoel connectivity (a subset of the entire triinpoel_complete connectivity), e.g., after mesh partitioning.
triinpoel in/out Interconnectivity of points and boundary-face in this mesh-partition.
bface in/out Map of boundary-face lists mapped to corresponding side set ids for this mesh-partition
Returns Number of boundary-faces on this chare/mesh-partition.

This function takes a mesh by its domain-element (tetrahedron-connectivity) in inpoel and a boundary-face (triangle) connectivity in triinpoel_complete. Based on these two arrays, it searches for those faces of triinpoel_complete that are also in inpoel and as a result it generates (1) the number of boundary faces shared with the mesh in inpoel and (2) the intersection of the triangle element connectivity whose faces are shared with inpoel. An example use case is where triinpoel_complete contains the connectivity for the boundary of the full problem/mesh and inpoel contains the connectivity for only a chunk of an already partitioned mesh. This function then intersects triinpoel_complete with inpoel and returns only those faces that share nodes with inpoel.

std::vector<int> tk::genEsuelTet(const std::vector<std::size_t>& inpoel, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup)

Generate derived data structure, elements surrounding elements as a fixed length data structure as a full vector, including boundary elements as -1.

Parameters
inpoel in

Inteconnectivity of points and elements. These are the node ids of each element of an unstructured mesh. Example:

std::vector< std::size_t > inpoel { 12, 14,  9, 11,
                                    10, 14, 13, 12 };

specifies two tetrahedra whose vertices (node ids) are { 12, 14, 9, 11 }, and { 10, 14, 13, 12 }.

esup in Elements surrounding points as linked lists, see tk::genEsup
Returns Vector storing elements surrounding elements

The data generated here is stored in a single vector, with length nfpe * nelem. Note however, that nelem is not explicitly provided, but calculated from inpoel. For boundary elements, at the boundary face, this esuelTet stores value -1 indicating that this is outside the domain. The convention for numbering the local face (triangle) connectivity is very important, e.g., in generating the inpofa array later. This node ordering convention is stored in tk::lpofa. Thus function is specific to tetrahedra, which is reflected in the fact that nnpe and nfpe are being set here in the function rather than being input arguments. To find out the number of elements, nelem, the size of the mesh connectivity vector, inpoel, can be devided by the number of nodes per elements, nnpe: auto nelem = inpoel.size()/nnpe;

std::size_t tk::genNipfac(std::size_t nfpe, std::size_t nbfac, const std::vector<int>& esuelTet)

Generate number of internal and physical-boundary faces.

Parameters
nfpe in Number of faces per element.
nbfac in Number of boundary faces.
esuelTet in Elements surrounding elements.
Returns Total number of faces in the mesh

The unsigned integer here gives the number of internal and physical-boundary faces in the mesh. The data structure does not include faces that are on partition/chare-boundaries.

std::vector<int> tk::genEsuf(std::size_t nfpe, std::size_t nipfac, std::size_t nbfac, const std::vector<std::size_t>& belem, const std::vector<int>& esuelTet)

Generate derived data structure, elements surrounding faces.

Parameters
nfpe in Number of faces per element.
nipfac in Number of internal and physical-boundary faces.
nbfac in Number of boundary faces.
belem in Boundary element vector.
esuelTet in Elements surrounding elements.
Returns Elements surrounding faces.

The unsigned integer vector gives the IDs of the elements to the

std::vector<std::size_t> tk::genInpofaTet(std::size_t nipfac, std::size_t nbfac, const std::vector<std::size_t>& inpoel, const std::vector<std::size_t>& triinpoel, const std::vector<int>& esuelTet)

Generate derived data structure, points on faces for tetrahedra only.

Parameters
nipfac in Number of internal and physical-boundary faces.
nbfac in Number of boundary faces.
inpoel in Element-node connectivity.
triinpoel in Face-node connectivity.
esuelTet in Elements surrounding elements.
Returns Points surrounding faces. The unsigned integer vector gives the elements to the left and to the right of each face in the mesh.

std::vector<std::size_t> tk::genBelemTet(std::size_t nbfac, const std::vector<std::size_t>& inpofa, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup)

Generate derived data structure, boundary elements.

Parameters
nbfac in Number of boundary faces.
inpofa in Face-node connectivity.
esup in Elements surrounding points as linked lists, see tk::genEsup
Returns Host elements or boundary elements. The unsigned integer vector gives the elements to the left of each boundary face in the mesh.

The data structure generated here contains an array of elements which share one or more of their faces with the physical boundary, i.e., where exodus specifies a side-set for faces. Such elements are sometimes also called host or boundary elements.

bool tk::leakyPartition(const std::vector<int>& esueltet, const std::vector<std::size_t>& inpoel, const std::array<std::vector<real>, 3>& coord)

Perform leak-test on mesh (partition)

Parameters
esueltet in Elements surrounding elements for tetrahedra, see tk::genEsueltet()
inpoel in Element connectivity
coord in Node coordinates
Returns True if partition leaks.

This function computes a surface integral over the boundary of the incoming mesh (partition). A non-zero vector result indicates a leak, e.g., a hole in the mesh (partition), which indicates an error either in the

bool tk::conforming(const std::vector<std::size_t>& inpoel, const std::array<std::vector<real>, 3>& coord, bool cerr, const std::vector<std::size_t>& rid)

Check if mesh (partition) is conforming.

Parameters
inpoel in Element connectivity
coord in Node coordinates
cerr in True if hanging-node edge data should be output to std::cerr (true by default)
rid in AMR Lib node id map std::cerr (true by default)
Returns True if mesh (partition) has no hanging nodes and thus the mesh is conforming, false if non-conforming.

A conforming mesh by definition has no hanging nodes. A node is hanging if an edge of one element coincides with two (or more) edges (of two or more other elements). Thus, testing for conformity relies on checking the coordinates of all vertices: if any vertex coincides with that of a mid-point node of an edge, that is a hanging node. Note that this assumes that hanging nodes can only be at the mid-point of edges. This may happen after a mesh refinement step, due to a problem/bug, within the mesh refinement algorithm given by J. Waltz, Parallel adaptive refinement for unsteady flow calculations on 3D unstructured grids, International Journal for Numerical Methods in Fluids, 46: 37–57, 2004, which always adds/removes vertices at the mid-points of edges of a tetrahedron mesh within a single refinement step. Thus this algorithm is intended for this specific case, i.e., test for conformity after a single refinement step and not after multiple ones or for detecting hanging nodes in an arbitrary mesh.

bool tk::intet(const std::array<std::vector<real>, 3>& coord, const std::vector<std::size_t>& inpoel, const std::vector<real>& p, std::size_t e, std::array<real, 4>& N)

Determine if a point is in a tetrahedron.

Parameters
coord in Mesh node coordinates
inpoel in Mesh element connectivity
in Point coordinates
in Mesh cell index
in/out Shapefunctions evaluated at the point
Returns True if ppoint is in mesh cell

std::array<tk::real, 3> tk::nodegrad(std::size_t node, const std::array<std::vector<tk::real>, 3>& coord, const std::vector<std::size_t>& inpoel, const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& esup, const tk::Fields& U, uint64_t c)

Compute gradient at a mesh node.

Parameters
node in Node id at which to compute gradient
coord in Mesh node coordinates
inpoel in Mesh element connectivity
esup in Linked lists storing elements surrounding points, see tk::genEsup()
in Field vector whose component gradient to compute
in Scalar component to compute gradient of
Returns Gradient of U(c) at mesh node

std::array<tk::real, 3> tk::edgegrad(const std::array<std::vector<tk::real>, 3>& coord, const std::vector<std::size_t>& inpoel, const std::vector<std::size_t>& esued, const tk::Fields& U, uint64_t c)

Compute gradient at a mesh edge.

Parameters
coord in Mesh node coordinates
inpoel in Mesh element connectivity
esued in List of elements surrounding edge, see tk::genEsued()
in Field vector whose component gradient to compute
in Scalar component to compute gradient of
Returns Gradient of U(c) at mesh edge

std::size_t tk::shiftToZero(std::vector<std::size_t>& inpoel)

Shift node IDs to start with zero in element connectivity.

Parameters
inpoel in/out Inteconnectivity of points and elements
Returns Amount shifted

This function implements a simple reordering of the node ids of the element connectivity in inpoel by shifting the node ids so that the smallest is zero.

void tk::remap(std::vector<std::size_t>& ids, const std::vector<std::size_t>& map)

Apply new mapping to vector of indices.

Parameters
ids in/out Vector of integer IDs to remap
map in Array of indices creating a new order

This function applies a mapping (reordering) to the integer IDs passed in using the map passed in. The mapping is expressed between the array index and its value. The function overwrites every value, i, of vector ids with map[i].

void tk::remap(std::vector<tk::real>& r, const std::vector<std::size_t>& map)

Apply new mapping to vector of real numbers.

Parameters
in/out Vector of real numbers to remap
map in Array of indices creating a new order

This function applies a mapping (reordering) to the real values passed in using the map passed in. The mapping is expressed between the array index and its value. The function moves every value r[i] to r[ map[i] ].

std::vector<std::size_t> tk::remap(const std::vector<std::size_t>& ids, const std::vector<std::size_t>& map)

Create remapped vector of indices using a vector.

Parameters
ids in Vector of integer IDs to remap
map in Array of indices creating a new order
Returns Remapped vector of ids

This function applies a mapping (reordering) to the integer IDs passed in using the map passed in. The mapping is expressed between the array index and its value. The function creates and returns a new container with remapped ids of identical size of the origin ids container.

void tk::remap(std::vector<std::size_t>& ids, const std::unordered_map<std::size_t, std::size_t>& map)

In-place remap vector of indices using a map.

Parameters
ids in Vector of integer IDs to remap
map in Hash-map of key->value creating a new order

This function applies a mapping (reordering) to the integer IDs passed in using the map passed in. The mapping is expressed as a hash-map of key->value pairs, where the key is the original and the value is the new ids of the mapping. The function overwrites the ids container with the remapped ids of identical size.

std::vector<std::size_t> tk::remap(const std::vector<std::size_t>& ids, const std::unordered_map<std::size_t, std::size_t>& map)

Create remapped vector of indices using a map.

Parameters
ids in Vector of integer IDs to create new container of ids from
map in Hash-map of key->value creating a new order
Returns Remapped vector of ids

This function applies a mapping (reordering) to the integer IDs passed in using the map passed in. The mapping is expressed as a hash-map of key->value pairs, where the key is the original and the value is the new ids of the mapping. The function creates and returns a new container with the remapped ids of identical size of the original ids container.

std::map<int, std::vector<std::size_t>> tk::remap(const std::map<int, std::vector<std::size_t>>& ids, const std::unordered_map<std::size_t, std::size_t>& map)

Create remapped map of vector of indices using a map.

Parameters
ids in Map of vector of integer IDs to create new container of ids from
map in Hash-map of key->value creating a new order
Returns Remapped vector of ids

This function applies a mapping (reordering) to the map of integer IDs passed in using the map passed in by applying remap(vector,map) on each vector of ids. The keys in the returned map will be the same as in ids.

std::vector<std::size_t> tk::renumber(const std::pair<std::vector<std::size_t>, std::vector<std::size_t>>& psup)

Reorder mesh points with the advancing front technique.

Parameters
psup in Points surrounding points
Returns Mapping created by renumbering (reordering)

std::unordered_map<std::size_t, std::size_t> tk::assignLid(const std::vector<std::size_t>& gid)

Assign local ids to global ids.

Parameters
gid in Global ids
Returns Map associating global ids to local ids

std::tuple<std::vector<std::size_t>, std::vector<std::size_t>, std::unordered_map<std::size_t, std::size_t>> tk::global2local(const std::vector<std::size_t>& ginpoel)

Generate element connectivity of local node IDs from connectivity of global node IDs also returning the mapping between local to global IDs.

Parameters
ginpoel in Element connectivity with global node IDs
Returns Tuple of (1) element connectivity with local node IDs, (2) the vector of unique global node IDs (i.e., the mapping between local to global node IDs), and (3) mapping between global to local node IDs.

bool tk::positiveJacobians(const std::vector<std::size_t>& inpoel, const std::array<std::vector<real>, 3>& coord)

Test positivity of the Jacobian for all cells in a mesh.

Parameters
inpoel in Element connectivity (zero-based, i.e., local if parallel)
coord in Node coordinates
Returns True if Jacobians of all mesh cells are positive

std::map<int, std::vector<std::size_t>> tk::bfacenodes(const std::map<int, std::vector<std::size_t>>& bface, const std::vector<std::size_t>& triinpoel)

Generate nodes of side set faces.

Parameters
bface in Boundary-faces mapped to side set ids
triinpoel in Boundary-face connectivity
Returns Nodes of side set faces for each side set passed in

tk::real tk::count(const std::unordered_map<int, std::unordered_set<std::size_t>>& map, std::size_t node)

Count the number of contributions to a node.

Parameters
map in Node commuinication map to search in
node in Global node id to search for
Returns Count of contributions to node

bool tk::slave(const std::unordered_map<int, std::unordered_set<std::size_t>>& map, std::size_t node, int chare)

Decide if a node is not counted by a chare.

Parameters
map in Node commuinication map to search in
node in Global node id to search for
chare in Caller chare id (but really can be any chare id)
Returns True if the node is a slave (counted by another chare with a lower chare id)

If a node is found in the node communication map and is associated to a lower chare id than the chare id given, it is counted by another chare (and not the caller one), hence a "slave" (for the purpose of this count).

std::pair<int, std::unique_ptr<char[]>> tk::serialize(std::size_t meshid, const std::vector<tk::UniPDF>& u)

Serialize univariate PDF to raw memory stream.

Parameters
meshid in Mesh ID
in Univariate PDFs
Returns Pair of the length and the raw stream containing the serialized PDFs

CkReductionMsg* tk::mergeUniPDFs(int nmsg, CkReductionMsg** msgs)

Charm++ custom reducer for merging a univariate PDF during reduction across PEs.

Parameters
nmsg in Number of messages in msgs
msgs in Charm++ reduction message containing the serialized PDF
Returns Aggregated PDF built for further aggregation if needed

static std::ostream& tk::operator<<(std::ostream& os, const tk::UniPDF& p)

Parameters
os in/out Stream to output to
in PDF to output
Returns Updated stream

Output univariate PDF to output stream

Variable documentation

const uint8_t tk::UnkEqComp

Tags for selecting data layout policies.

const std::array<std::size_t, 2> tk::ExoNnpe

ExodusII mesh cell number of nodes

List of number of nodes per element for different element types supported in the order of tk::ExoElemType

const std::array<std::array<std::size_t, 3>, 4> tk::expofa

ExodusII face-node numbering for tetrahedron side sets

const std::array<UnsMesh::Face, 4> tk::lpofa

Const array defining the node ordering convention for tetrahedron faces

This two-dimensional array stores the naming/ordering convention of the node indices of a tetrahedron (tet) element. The dimensions are 4x3 as a tetrahedron has a total of 4 nodes and each (triangle) face has 3 nodes. Thus the array below associates tet node 0 with nodes {1,2,3}, tet node 1 with {2,0,3}, tet node 2 with {3,0,1}, and tet node 3 with {0,2,1}. Note that not only these mappings are important, but also the order of the nodes within the triplets as this specific order also defines the outwards normal of each face.

const std::array<UnsMesh::Edge, 6> tk::lpoed

Const array defining the node ordering convention for tetrahedron edges.

const std::array<UnsMesh::Edge, 3> tk::lpoet

Const array defining the node ordering convention for triangle edges.

static highwayhash::HH_U64 tk::hh_key constexpr

Highway hash "secret" key

static std::stringstream tk::cerr_quiet

std::tringstream used to quiet std::cerr's stream by redirecting to it

static std::streambuf* tk::cerr_old

std::streambuf used to store state of std::cerr before redirecting it