#ifndef MEDUSA_BITS_OPERATORS_SHAPESTORAGE_HPP_ #define MEDUSA_BITS_OPERATORS_SHAPESTORAGE_HPP_ /** * @file * Implementations of the ShapeStorage interface. */ #include "ShapeStorage_fwd.hpp" #include "printShapes.hpp" namespace mm { template int ShapeStorage::support(int node, int j) const { assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).", node, domain_size_); assert_msg(0 <= j && j < supportSize(node), "Support index %d out of range [0, %d).", j, supportSize(node)); assert_msg(access(support_, node)[j] != -1, "Attempted access to shapes and support for node " "%d, which were not computed.", node); return access(support_, node)[j]; } template Eigen::Map> ShapeStorage::support(int node) const { assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).", node, domain_size_); assert_msg(*access(support_, node) != -1, "Attempted access to shapes and support for node " "%d, which were not computed.", node); return {access(support_, node), supportSize(node)}; } template Range ShapeStorage::supportSizes() const { Range sizes(domain_size_, 0); for (int i = 0; i < domain_size_; ++i) { sizes[i] = supportSize(i); } return sizes; } template Range ShapeStorage::supportSizesVec() const { Range sizes = supportSizes(); Range res; for (int d = 0; d < dim; ++d) res.append(sizes); for (int& i : res) i *= dim; return res; } template typename ShapeStorage::scalar_t ShapeStorage::laplace(int node, int j) const { assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).", node, domain_size_); assert_msg(0 <= j && j < supportSize(node), "Support index %d out of range [0, %d).", j, supportSize(node)); return get>(node, j); } template typename ShapeStorage::scalar_t ShapeStorage::d1(int var, int node, int j) const { assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).", node, domain_size_); assert_msg(0 <= j && j < supportSize(node), "Support index %d out of range [0, %d).", j, supportSize(node)); assert_msg(0 <= var && var < dim, "Variable index %d out of bounds [%d, %d).", var, 0, dim); return get>(Der1s::index(Der1(var)), node, j); } template typename ShapeStorage::scalar_t ShapeStorage::d2( int varmin, int varmax, int node, int j) const { assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).", node, domain_size_); assert_msg(0 <= j && j < supportSize(node), "Support index %d out of range [0, %d).", j, supportSize(node)); assert_msg(0 <= varmin && varmin < dim, "Variable varmin %d out of bounds [%d, %d).", varmin, 0, dim); assert_msg(0 <= varmax && varmax < dim, "Variable varmax %d out of bounds [%d, %d).", varmax, 0, dim); assert_msg(varmin <= varmax, "Varmin (%d) must be smaller than varmax (%d).", varmin, varmax); return get>(Der2s::index(Der2(varmin, varmax)), node, j); } template Eigen::Map::scalar_t, Eigen::Dynamic, 1>> ShapeStorage::laplace(int node) const { assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).", node, domain_size_); return getShape>(node); } template Eigen::Map::scalar_t, Eigen::Dynamic, 1>> ShapeStorage::d1(int var, int node) const { assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).", node, domain_size_); assert_msg(0 <= var && var < dim, "Variable index %d out of bounds [%d, %d).", var, 0, dim); return getShape>(Der1s::index(Der1(var)), node); } template Eigen::Map::scalar_t, Eigen::Dynamic, 1>> ShapeStorage::d2(int varmin, int varmax, int node) const { assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).", node, domain_size_); assert_msg(0 <= varmin && varmin < dim, "Variable varmin %d out of bounds [%d, %d).", varmin, 0, dim); assert_msg(0 <= varmax && varmax < dim, "Variable varmax %d out of bounds [%d, %d).", varmax, 0, dim); assert_msg(varmin <= varmax, "Varmin (%d) must be smaller than varmax (%d).", varmin, varmax); return getShape>(Der2s::index(Der2(varmin, varmax)), node); } template void ShapeStorage::setSupport( int node, const std::vector& support) { assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).", node, domain_size_); assert_msg(static_cast(support.size()) == supportSize(node), "Support of unexpected " "length %d, expected %d.", support.size(), supportSize(node)); std::memcpy(access(support_, node), support.data(), support.size()*sizeof(int)); } template void ShapeStorage::setLaplace( int node, const Eigen::Matrix& shape) { assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).", node, domain_size_); assert_msg(shape.size() == supportSize(node), "Laplace shape of unexpected length %d, " "expected %d.", shape.size(), supportSize(node)); setShape>(node, shape); } template void ShapeStorage::setD1( int var, int node, const Eigen::Matrix& shape) { assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).", node, domain_size_); assert_msg(0 <= var && var < dim, "Variable index %d out of bounds [%d, %d).", var, 0, dim); assert_msg(shape.size() == supportSize(node), "D1 shape of unexpected length %d, expected %d.", shape.size(), supportSize(node)); setShape>(var, node, shape); } template void ShapeStorage::setD2( int varmin, int varmax, int node, const Eigen::Matrix& shape) { assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).", node, domain_size_); assert_msg(0 <= varmin && varmin < dim, "Variable varmin %d out of bounds [%d, %d).", varmin, 0, dim); assert_msg(0 <= varmax && varmax < dim, "Variable varmax %d out of bounds [%d, %d).", varmax, 0, dim); assert_msg(varmin <= varmax, "Varmin (%d) must be smaller than varmax (%d).", varmin, varmax); assert_msg(shape.size() == supportSize(node), "D2 shape of unexpected length %d, expected %d.", shape.size(), supportSize(node)); int idx = varmax*(varmax+1)/2 + varmin; setShape>(idx, node, shape); } /// Output basic info about this shape storage. template std::ostream& operator<<(std::ostream& os, const ShapeStorage& shapes) { os << "Shape storage:\n"; return shapes_internal::printShapes(shapes, os); } template size_t ShapeStorage::memoryUsed() const { size_t size = sizeof(*this); for (int i = 0; i < num_operators; ++i) { size += mem_used(shapes_[i]); } return size; } /// Contains a helper for resizing storage. namespace resize_internal { /// @cond // A helper for resizing storage -- recursive case template struct resize_storage_ { template static void resize(Shapes& shapes, int base_size) { resize_storage_::template resize(shapes, base_size); std::get(shapes).resize( base_size*std::tuple_element::type::size(), NaN); } }; // A helper for resizing storage -- base case template struct resize_storage_ { template static void resize(Shapes&, int) {} }; /// @endcond } // namespace resize_internal } // namespace mm #endif // MEDUSA_BITS_OPERATORS_SHAPESTORAGE_HPP_