Medusa  1.1
Coordinate Free Mehless Method implementation
ShapeStorage.hpp
Go to the documentation of this file.
1 #ifndef MEDUSA_BITS_OPERATORS_SHAPESTORAGE_HPP_
2 #define MEDUSA_BITS_OPERATORS_SHAPESTORAGE_HPP_
3 
9 #include "ShapeStorage_fwd.hpp"
10 #include "printShapes.hpp"
11 
12 namespace mm {
13 
14 template <typename Derived, typename vec_t, typename OpFamilies>
16  assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).",
17  node, domain_size_);
18  assert_msg(0 <= j && j < supportSize(node), "Support index %d out of range [0, %d).",
19  j, supportSize(node));
20  assert_msg(access(support_, node)[j] != -1, "Attempted access to shapes and support for node "
21  "%d, which were not computed.", node);
22  return access(support_, node)[j];
23 }
24 
25 template <typename Derived, typename vec_t, typename OpFamilies>
26 Eigen::Map<const Eigen::Matrix<int, Eigen::Dynamic, 1>>
28  assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).",
29  node, domain_size_);
30  assert_msg(*access(support_, node) != -1, "Attempted access to shapes and support for node "
31  "%d, which were not computed.", node);
32  return {access(support_, node), supportSize(node)};
33 }
34 
35 template <typename Derived, typename vec_t, typename OpFamilies>
37  Range<int> sizes(domain_size_, 0);
38  for (int i = 0; i < domain_size_; ++i) {
39  sizes[i] = supportSize(i);
40  }
41  return sizes;
42 }
43 
44 template <typename Derived, typename vec_t, typename OpFamilies>
46  Range<int> sizes = supportSizes();
47  Range<int> res;
48  for (int d = 0; d < dim; ++d) res.append(sizes);
49  for (int& i : res) i *= dim;
50  return res;
51 }
52 
53 template <typename Derived, typename vec_t, typename OpFamilies>
56  assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).",
57  node, domain_size_);
58  assert_msg(0 <= j && j < supportSize(node), "Support index %d out of range [0, %d).",
59  j, supportSize(node));
60  return get<Lap<vec_t::dim>>(node, j);
61 }
62 
63 template <typename Derived, typename vec_t, typename OpFamilies>
65 ShapeStorage<Derived, vec_t, OpFamilies>::d1(int var, int node, int j) const {
66  assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).",
67  node, domain_size_);
68  assert_msg(0 <= j && j < supportSize(node), "Support index %d out of range [0, %d).",
69  j, supportSize(node));
70  assert_msg(0 <= var && var < dim, "Variable index %d out of bounds [%d, %d).", var, 0, dim);
71  return get<Der1s<dim>>(Der1s<dim>::index(Der1<dim>(var)), node, j);
72 }
73 
74 template <typename Derived, typename vec_t, typename OpFamilies>
77  int varmin, int varmax, int node, int j) const {
78  assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).",
79  node, domain_size_);
80  assert_msg(0 <= j && j < supportSize(node), "Support index %d out of range [0, %d).",
81  j, supportSize(node));
82  assert_msg(0 <= varmin && varmin < dim,
83  "Variable varmin %d out of bounds [%d, %d).", varmin, 0, dim);
84  assert_msg(0 <= varmax && varmax < dim,
85  "Variable varmax %d out of bounds [%d, %d).", varmax, 0, dim);
86  assert_msg(varmin <= varmax, "Varmin (%d) must be smaller than varmax (%d).",
87  varmin, varmax);
88  return get<Der2s<dim>>(Der2s<dim>::index(Der2<dim>(varmin, varmax)), node, j);
89 }
90 
91 template <typename Derived, typename vec_t, typename OpFamilies>
93  Eigen::Dynamic, 1>>
95  assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).",
96  node, domain_size_);
97  return getShape<Lap<vec_t::dim>>(node);
98 }
99 
100 template <typename Derived, typename vec_t, typename OpFamilies>
102  Eigen::Dynamic, 1>>
104  assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).",
105  node, domain_size_);
106  assert_msg(0 <= var && var < dim, "Variable index %d out of bounds [%d, %d).", var, 0, dim);
107  return getShape<Der1s<dim>>(Der1s<dim>::index(Der1<dim>(var)), node);
108 }
109 
110 template <typename Derived, typename vec_t, typename OpFamilies>
112  Eigen::Dynamic, 1>>
113 ShapeStorage<Derived, vec_t, OpFamilies>::d2(int varmin, int varmax, int node) const {
114  assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).",
115  node, domain_size_);
116  assert_msg(0 <= varmin && varmin < dim,
117  "Variable varmin %d out of bounds [%d, %d).", varmin, 0, dim);
118  assert_msg(0 <= varmax && varmax < dim,
119  "Variable varmax %d out of bounds [%d, %d).", varmax, 0, dim);
120  assert_msg(varmin <= varmax, "Varmin (%d) must be smaller than varmax (%d).",
121  varmin, varmax);
122  return getShape<Der2s<dim>>(Der2s<dim>::index(Der2<dim>(varmin, varmax)), node);
123 }
124 
125 template <typename Derived, typename vec_t, typename OpFamilies>
127  int node, const std::vector<int>& support) {
128  assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).",
129  node, domain_size_);
130  assert_msg(static_cast<int>(support.size()) == supportSize(node), "Support of unexpected "
131  "length %d, expected %d.", support.size(), supportSize(node));
132  std::memcpy(access(support_, node), support.data(),
133  support.size()*sizeof(int));
134 }
135 
136 template <typename Derived, typename vec_t, typename OpFamilies>
138  int node, const Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>& shape) {
139  assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).",
140  node, domain_size_);
141  assert_msg(shape.size() == supportSize(node), "Laplace shape of unexpected length %d, "
142  "expected %d.", shape.size(), supportSize(node));
143  setShape<Lap<vec_t::dim>>(node, shape);
144 }
145 
146 template <typename Derived, typename vec_t, typename OpFamilies>
148  int var, int node, const Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>& shape) {
149  assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).",
150  node, domain_size_);
151  assert_msg(0 <= var && var < dim, "Variable index %d out of bounds [%d, %d).", var, 0, dim);
152  assert_msg(shape.size() == supportSize(node), "D1 shape of unexpected length %d, expected %d.",
153  shape.size(), supportSize(node));
154  setShape<Der1s<dim>>(var, node, shape);
155 }
156 
157 template <typename Derived, typename vec_t, typename OpFamilies>
159  int varmin, int varmax, int node, const Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>& shape) {
160  assert_msg(0 <= node && node <= domain_size_, "Node index %d out or range [0, %d).",
161  node, domain_size_);
162  assert_msg(0 <= varmin && varmin < dim,
163  "Variable varmin %d out of bounds [%d, %d).", varmin, 0, dim);
164  assert_msg(0 <= varmax && varmax < dim,
165  "Variable varmax %d out of bounds [%d, %d).", varmax, 0, dim);
166  assert_msg(varmin <= varmax, "Varmin (%d) must be smaller than varmax (%d).",
167  varmin, varmax);
168  assert_msg(shape.size() == supportSize(node), "D2 shape of unexpected length %d, expected %d.",
169  shape.size(), supportSize(node));
170  int idx = varmax*(varmax+1)/2 + varmin;
171  setShape<Der2s<dim>>(idx, node, shape);
172 }
173 
175 template <typename D, typename V, typename O>
176 std::ostream& operator<<(std::ostream& os, const ShapeStorage<D, V, O>& shapes) {
177  os << "Shape storage:\n";
178  return shapes_internal::printShapes(shapes, os);
179 }
180 
181 template <typename Derived, typename vec_t, typename OpFamilies>
183  size_t size = sizeof(*this);
184  for (int i = 0; i < num_operators; ++i) {
185  size += mem_used(shapes_[i]);
186  }
187  return size;
188 }
189 
191 namespace resize_internal {
192 
194 // A helper for resizing storage -- recursive case
195 template <typename OpFamilies, size_t i>
196 struct resize_storage_ {
197  template <typename Shapes>
198  static void resize(Shapes& shapes, int base_size) {
199  resize_storage_<OpFamilies, i-1>::template resize<Shapes>(shapes, base_size);
200  std::get<i-1>(shapes).resize(
201  base_size*std::tuple_element<i-1, OpFamilies>::type::size(), NaN);
202  }
203 };
204 
205 // A helper for resizing storage -- base case
206 template <typename OpFamilies>
207 struct resize_storage_<OpFamilies, 0> {
208  template <typename Shapes> static void resize(Shapes&, int) {}
209 };
211 
212 } // namespace resize_internal
213 
214 } // namespace mm
215 
216 #endif // MEDUSA_BITS_OPERATORS_SHAPESTORAGE_HPP_
mm::ShapeStorage::setD1
SINL void setD1(int var, int node, const Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > &shape)
Sets shape for derivative wrt. variable var for node to shape.
Definition: ShapeStorage.hpp:147
mm
Root namespace for the whole library.
Definition: Gaussian.hpp:14
scalar_t
Scalar scalar_t
Type of the elements, alias of Scalar.
Definition: MatrixBaseAddons.hpp:16
ShapeStorage_fwd.hpp
printShapes.hpp
dim
@ dim
Number of elements of this matrix.
Definition: MatrixBaseAddons.hpp:14
mm::NaN
static const double NaN
Not-a-number floating point value.
Definition: Config.hpp:46
mm::ShapeStorage::d1
SINL scalar_t d1(int var, int node, int j) const
Return j-th shape coefficient for derivative wrt. variable var in node.
Definition: ShapeStorage.hpp:65
mm::operator<<
std::ostream & operator<<(std::ostream &os, const Gaussian< S > &b)
Output basic information about given Gaussian RBF.
Definition: Gaussian.hpp:37
mm::ShapeStorage
Shape storage base class.
Definition: ShapeStorage_fwd.hpp:41
assert_msg
#define assert_msg(cond,...)
Assert with better error reporting.
Definition: assert.hpp:75
mm::Der1s::index
static int index(Der1< dim > d)
Get index of operator d in the family.
Definition: Operators_fwd.hpp:109
mm::ShapeStorage::setLaplace
SINL void setLaplace(int node, const Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > &shape)
Sets the laplace shape for node to shape.
Definition: ShapeStorage.hpp:137
mm::ShapeStorage::scalar_t
vec_t::scalar_t scalar_t
Scalar type used.
Definition: ShapeStorage_fwd.hpp:44
mm::ShapeStorage::laplace
SINL scalar_t laplace(int node, int j) const
Return j-th laplace shape coefficient for node-th node.
Definition: ShapeStorage.hpp:55
mm::ShapeStorage::support
SINL int support(int node, int j) const
Returns index of j-th neighbour of node-th node, i.e. support[node][j], but possibly faster.
Definition: ShapeStorage.hpp:15
mm::Der2
Represents a second derivative wrt. var1 and var2.
Definition: Monomials_fwd.hpp:24
mm::ShapeStorage::supportSizes
Range< int > supportSizes() const
Returns a vector of support sizes for all nodes, useful for matrix space preallocation.
Definition: ShapeStorage.hpp:36
mm::ShapeStorage::setSupport
void setSupport(int node_idx, const std::vector< int > &support)
Sets support of node-th node to support.
Definition: ShapeStorage.hpp:126
mm::ShapeStorage::setD2
SINL void setD2(int varmin, int varmax, int node, const Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > &shape)
Sets shape for mixed derivative wrt. variables varmin and varmax for node to shape.
Definition: ShapeStorage.hpp:158
mm::ShapeStorage::supportSizesVec
Range< int > supportSizesVec() const
Returns a dim*N vector of dim*support_size, useful for matrix space preallocation for vector equation...
Definition: ShapeStorage.hpp:45
mm::Range::append
Range & append(const Range &rng)
Append all elements of rng to self.
Definition: Range.hpp:66
mm::ShapeStorage::memoryUsed
size_t memoryUsed() const
Returns the approximate memory used (in bytes).
Definition: ShapeStorage.hpp:182
mm::ShapeStorage::d2
SINL scalar_t d2(int varmin, int varmax, int node, int j) const
Return j-th shape coefficient for mixed derivative wrt.
Definition: ShapeStorage.hpp:76
mm::Der1
Represents a first derivative wrt. var.
Definition: Monomials_fwd.hpp:22
mm::Range< int >
mm::mem_used
std::size_t mem_used(const container_t &v)
Returns number of bytes the container uses in memory.
Definition: memutils.hpp:83
mm::Der2s::index
static int index(Der2< dim > d)
Get index of operator d in the family.
Definition: Operators_fwd.hpp:132