Sleipnir C++ API
Loading...
Searching...
No Matches
variable_matrix.hpp
1// Copyright (c) Sleipnir contributors
2
3#pragma once
4
5#include <algorithm>
6#include <concepts>
7#include <initializer_list>
8#include <iterator>
9#include <span>
10#include <utility>
11#include <vector>
12
13#include <Eigen/Core>
14#include <Eigen/QR>
15#include <gch/small_vector.hpp>
16
17#include "sleipnir/autodiff/sleipnir_base.hpp"
18#include "sleipnir/autodiff/slice.hpp"
19#include "sleipnir/autodiff/variable.hpp"
20#include "sleipnir/autodiff/variable_block.hpp"
21#include "sleipnir/util/assert.hpp"
22#include "sleipnir/util/concepts.hpp"
23#include "sleipnir/util/empty.hpp"
24#include "sleipnir/util/function_ref.hpp"
25#include "sleipnir/util/symbol_exports.hpp"
26
27namespace slp {
28
32template <typename Scalar_>
34 public:
36 using Scalar = Scalar_;
37
39 VariableMatrix() = default;
40
45 explicit VariableMatrix(int rows) : VariableMatrix{rows, 1} {}
46
51 VariableMatrix(int rows, int cols) : m_rows{rows}, m_cols{cols} {
52 m_storage.reserve(rows * cols);
53 for (int index = 0; index < rows * cols; ++index) {
54 m_storage.emplace_back();
55 }
56 }
57
63 : m_rows{rows}, m_cols{cols} {
64 m_storage.reserve(rows * cols);
65 for (int index = 0; index < rows * cols; ++index) {
66 m_storage.emplace_back(nullptr);
67 }
68 }
69
74 std::initializer_list<std::initializer_list<Variable<Scalar>>> list) {
75 // Get row and column counts for destination matrix
76 m_rows = list.size();
77 m_cols = 0;
78 if (list.size() > 0) {
79 m_cols = list.begin()->size();
80 }
81
82 // Assert all column counts are the same
83 for ([[maybe_unused]]
84 const auto& row : list) {
85 slp_assert(static_cast<int>(row.size()) == m_cols);
86 }
87
88 m_storage.reserve(rows() * cols());
89 for (const auto& row : list) {
90 std::ranges::copy(row, std::back_inserter(m_storage));
91 }
92 }
93
99 // NOLINTNEXTLINE (google-explicit-constructor)
100 VariableMatrix(const std::vector<std::vector<Scalar>>& list) {
101 // Get row and column counts for destination matrix
102 m_rows = list.size();
103 m_cols = 0;
104 if (list.size() > 0) {
105 m_cols = list.begin()->size();
106 }
107
108 // Assert all column counts are the same
109 for ([[maybe_unused]]
110 const auto& row : list) {
111 slp_assert(static_cast<int>(row.size()) == m_cols);
112 }
113
114 m_storage.reserve(rows() * cols());
115 for (const auto& row : list) {
116 std::ranges::copy(row, std::back_inserter(m_storage));
117 }
118 }
119
125 // NOLINTNEXTLINE (google-explicit-constructor)
126 VariableMatrix(const std::vector<std::vector<Variable<Scalar>>>& list) {
127 // Get row and column counts for destination matrix
128 m_rows = list.size();
129 m_cols = 0;
130 if (list.size() > 0) {
131 m_cols = list.begin()->size();
132 }
133
134 // Assert all column counts are the same
135 for ([[maybe_unused]]
136 const auto& row : list) {
137 slp_assert(static_cast<int>(row.size()) == m_cols);
138 }
139
140 m_storage.reserve(rows() * cols());
141 for (const auto& row : list) {
142 std::ranges::copy(row, std::back_inserter(m_storage));
143 }
144 }
145
149 template <typename Derived>
150 // NOLINTNEXTLINE (google-explicit-constructor)
151 VariableMatrix(const Eigen::MatrixBase<Derived>& values)
152 : m_rows{static_cast<int>(values.rows())},
153 m_cols{static_cast<int>(values.cols())} {
154 m_storage.reserve(values.rows() * values.cols());
155 for (int row = 0; row < values.rows(); ++row) {
156 for (int col = 0; col < values.cols(); ++col) {
157 m_storage.emplace_back(values[row, col]);
158 }
159 }
160 }
161
165 template <typename Derived>
166 // NOLINTNEXTLINE (google-explicit-constructor)
167 VariableMatrix(const Eigen::DiagonalBase<Derived>& values)
168 : m_rows{static_cast<int>(values.rows())},
169 m_cols{static_cast<int>(values.cols())} {
170 m_storage.reserve(values.rows() * values.cols());
171 for (int row = 0; row < values.rows(); ++row) {
172 for (int col = 0; col < values.cols(); ++col) {
173 if (row == col) {
174 m_storage.emplace_back(values.diagonal()[row]);
175 } else {
176 m_storage.emplace_back(Scalar(0));
177 }
178 }
179 }
180 }
181
185 // NOLINTNEXTLINE (google-explicit-constructor)
186 VariableMatrix(const Variable<Scalar>& variable) : m_rows{1}, m_cols{1} {
187 m_storage.emplace_back(variable);
188 }
189
193 // NOLINTNEXTLINE (google-explicit-constructor)
194 VariableMatrix(Variable<Scalar>&& variable) : m_rows{1}, m_cols{1} {
195 m_storage.emplace_back(std::move(variable));
196 }
197
201 // NOLINTNEXTLINE (google-explicit-constructor)
203 : m_rows{values.rows()}, m_cols{values.cols()} {
204 m_storage.reserve(rows() * cols());
205 for (int row = 0; row < rows(); ++row) {
206 for (int col = 0; col < cols(); ++col) {
207 m_storage.emplace_back(values[row, col]);
208 }
209 }
210 }
211
215 // NOLINTNEXTLINE (google-explicit-constructor)
217 : m_rows{values.rows()}, m_cols{values.cols()} {
218 m_storage.reserve(rows() * cols());
219 for (int row = 0; row < rows(); ++row) {
220 for (int col = 0; col < cols(); ++col) {
221 m_storage.emplace_back(values[row, col]);
222 }
223 }
224 }
225
229 explicit VariableMatrix(std::span<const Variable<Scalar>> values)
230 : m_rows{static_cast<int>(values.size())}, m_cols{1} {
231 m_storage.reserve(rows() * cols());
232 for (int row = 0; row < rows(); ++row) {
233 for (int col = 0; col < cols(); ++col) {
234 m_storage.emplace_back(values[row * cols() + col]);
235 }
236 }
237 }
238
244 VariableMatrix(std::span<const Variable<Scalar>> values, int rows, int cols)
245 : m_rows{rows}, m_cols{cols} {
246 slp_assert(static_cast<int>(values.size()) == rows * cols);
247 m_storage.reserve(rows * cols);
248 for (int row = 0; row < rows; ++row) {
249 for (int col = 0; col < cols; ++col) {
250 m_storage.emplace_back(values[row * cols + col]);
251 }
252 }
253 }
254
259 template <typename Derived>
260 VariableMatrix& operator=(const Eigen::MatrixBase<Derived>& values) {
261 slp_assert(rows() == values.rows() && cols() == values.cols());
262
263 for (int row = 0; row < values.rows(); ++row) {
264 for (int col = 0; col < values.cols(); ++col) {
265 (*this)[row, col] = values[row, col];
266 }
267 }
268
269 return *this;
270 }
271
279 slp_assert(rows() == 1 && cols() == 1);
280
281 (*this)[0, 0] = value;
282
283 return *this;
284 }
285
289 template <typename Derived>
290 requires std::same_as<typename Derived::Scalar, Scalar>
291 void set_value(const Eigen::MatrixBase<Derived>& values) {
292 slp_assert(rows() == values.rows() && cols() == values.cols());
293
294 for (int row = 0; row < values.rows(); ++row) {
295 for (int col = 0; col < values.cols(); ++col) {
296 (*this)[row, col].set_value(values[row, col]);
297 }
298 }
299 }
300
307 slp_assert(row >= 0 && row < rows());
308 slp_assert(col >= 0 && col < cols());
309 return m_storage[row * cols() + col];
310 }
311
317 const Variable<Scalar>& operator[](int row, int col) const {
318 slp_assert(row >= 0 && row < rows());
319 slp_assert(col >= 0 && col < cols());
320 return m_storage[row * cols() + col];
321 }
322
328 slp_assert(index >= 0 && index < rows() * cols());
329 return m_storage[index];
330 }
331
336 const Variable<Scalar>& operator[](int index) const {
337 slp_assert(index >= 0 && index < rows() * cols());
338 return m_storage[index];
339 }
340
349 int block_rows, int block_cols) {
350 slp_assert(row_offset >= 0 && row_offset <= rows());
351 slp_assert(col_offset >= 0 && col_offset <= cols());
352 slp_assert(block_rows >= 0 && block_rows <= rows() - row_offset);
353 slp_assert(block_cols >= 0 && block_cols <= cols() - col_offset);
355 }
356
365 int col_offset,
366 int block_rows,
367 int block_cols) const {
368 slp_assert(row_offset >= 0 && row_offset <= rows());
369 slp_assert(col_offset >= 0 && col_offset <= cols());
370 slp_assert(block_rows >= 0 && block_rows <= rows() - row_offset);
371 slp_assert(block_cols >= 0 && block_cols <= cols() - col_offset);
373 }
374
386
393 Slice col_slice) const {
394 int row_slice_length = row_slice.adjust(rows());
395 int col_slice_length = col_slice.adjust(cols());
396 return VariableBlock{*this, std::move(row_slice), row_slice_length,
397 std::move(col_slice), col_slice_length};
398 }
399
417
434
441 slp_assert(cols() == 1);
442 slp_assert(offset >= 0 && offset < rows());
443 slp_assert(length >= 0 && length <= rows() - offset);
444 return block(offset, 0, length, 1);
445 }
446
453 int length) const {
454 slp_assert(cols() == 1);
455 slp_assert(offset >= 0 && offset < rows());
456 slp_assert(length >= 0 && length <= rows() - offset);
457 return block(offset, 0, length, 1);
458 }
459
465 slp_assert(row >= 0 && row < rows());
466 return block(row, 0, 1, cols());
467 }
468
474 slp_assert(row >= 0 && row < rows());
475 return block(row, 0, 1, cols());
476 }
477
483 slp_assert(col >= 0 && col < cols());
484 return block(0, col, rows(), 1);
485 }
486
492 slp_assert(col >= 0 && col < cols());
493 return block(0, col, rows(), 1);
494 }
495
500 template <EigenMatrixLike LHS, SleipnirMatrixLike<Scalar> RHS>
502 slp_assert(lhs.cols() == rhs.rows());
503
504 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), rhs.cols());
505
506 for (int i = 0; i < lhs.rows(); ++i) {
507 for (int j = 0; j < rhs.cols(); ++j) {
508 Variable sum{Scalar(0)};
509 for (int k = 0; k < lhs.cols(); ++k) {
510 sum += lhs(i, k) * rhs[k, j];
511 }
512 result[i, j] = sum;
513 }
514 }
515
516 return result;
517 }
518
523 template <SleipnirMatrixLike<Scalar> LHS, EigenMatrixLike RHS>
525 slp_assert(lhs.cols() == rhs.rows());
526
527 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), rhs.cols());
528
529 for (int i = 0; i < lhs.rows(); ++i) {
530 for (int j = 0; j < rhs.cols(); ++j) {
531 Variable sum{Scalar(0)};
532 for (int k = 0; k < lhs.cols(); ++k) {
533 sum += lhs[i, k] * rhs(k, j);
534 }
535 result[i, j] = sum;
536 }
537 }
538
539 return result;
540 }
541
546 template <SleipnirMatrixLike<Scalar> LHS, SleipnirMatrixLike<Scalar> RHS>
548 slp_assert(lhs.cols() == rhs.rows());
549
550 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), rhs.cols());
551
552 for (int i = 0; i < lhs.rows(); ++i) {
553 for (int j = 0; j < rhs.cols(); ++j) {
554 Variable sum{Scalar(0)};
555 for (int k = 0; k < lhs.cols(); ++k) {
556 sum += lhs[i, k] * rhs[k, j];
557 }
558 result[i, j] = sum;
559 }
560 }
561
562 return result;
563 }
564
569 template <EigenMatrixLike LHS>
571 const Variable<Scalar>& rhs) {
572 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), lhs.cols());
573
574 for (int row = 0; row < result.rows(); ++row) {
575 for (int col = 0; col < result.cols(); ++col) {
576 result[row, col] = lhs[row, col] * rhs;
577 }
578 }
579
580 return result;
581 }
582
587 template <SleipnirMatrixLike<Scalar> LHS, ScalarLike RHS>
589 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), lhs.cols());
590
591 for (int row = 0; row < result.rows(); ++row) {
592 for (int col = 0; col < result.cols(); ++col) {
593 result[row, col] = lhs[row, col] * rhs;
594 }
595 }
596
597 return result;
598 }
599
604 template <EigenMatrixLike RHS>
606 const RHS& rhs) {
607 VariableMatrix<Scalar> result(detail::empty, rhs.rows(), rhs.cols());
608
609 for (int row = 0; row < result.rows(); ++row) {
610 for (int col = 0; col < result.cols(); ++col) {
611 result[row, col] = rhs[row, col] * lhs;
612 }
613 }
614
615 return result;
616 }
617
622 template <ScalarLike LHS, SleipnirMatrixLike<Scalar> RHS>
624 VariableMatrix<Scalar> result(detail::empty, rhs.rows(), rhs.cols());
625
626 for (int row = 0; row < result.rows(); ++row) {
627 for (int col = 0; col < result.cols(); ++col) {
628 result[row, col] = rhs[row, col] * lhs;
629 }
630 }
631
632 return result;
633 }
634
640 slp_assert(cols() == rhs.rows() && cols() == rhs.cols());
641
642 for (int i = 0; i < rows(); ++i) {
643 for (int j = 0; j < rhs.cols(); ++j) {
644 Variable sum{Scalar(0)};
645 for (int k = 0; k < cols(); ++k) {
646 sum += (*this)[i, k] * rhs[k, j];
647 }
648 (*this)[i, j] = sum;
649 }
650 }
651
652 return *this;
653 }
654
660 for (int row = 0; row < rows(); ++row) {
661 for (int col = 0; col < rhs.cols(); ++col) {
662 (*this)[row, col] *= rhs;
663 }
664 }
665
666 return *this;
667 }
668
674 template <EigenMatrixLike LHS>
676 const Variable<Scalar>& rhs) {
677 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), lhs.cols());
678
679 for (int row = 0; row < result.rows(); ++row) {
680 for (int col = 0; col < result.cols(); ++col) {
681 result[row, col] = lhs[row, col] / rhs;
682 }
683 }
684
685 return result;
686 }
687
693 template <SleipnirMatrixLike<Scalar> LHS, ScalarLike RHS>
695 friend VariableMatrix<Scalar> operator/(const LHS& lhs, const RHS& rhs) {
696 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), lhs.cols());
697
698 for (int row = 0; row < result.rows(); ++row) {
699 for (int col = 0; col < result.cols(); ++col) {
700 result[row, col] = lhs[row, col] / rhs;
701 }
702 }
703
704 return result;
705 }
706
712 template <SleipnirMatrixLike<Scalar> LHS>
714 const Variable<Scalar>& rhs) {
715 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), lhs.cols());
716
717 for (int row = 0; row < result.rows(); ++row) {
718 for (int col = 0; col < result.cols(); ++col) {
719 result[row, col] = lhs[row, col] / rhs;
720 }
721 }
722
723 return result;
724 }
725
731 for (int row = 0; row < rows(); ++row) {
732 for (int col = 0; col < cols(); ++col) {
733 (*this)[row, col] /= rhs;
734 }
735 }
736
737 return *this;
738 }
739
745 template <EigenMatrixLike LHS, SleipnirMatrixLike<Scalar> RHS>
747 slp_assert(lhs.rows() == rhs.rows() && lhs.cols() == rhs.cols());
748
749 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), lhs.cols());
750
751 for (int row = 0; row < result.rows(); ++row) {
752 for (int col = 0; col < result.cols(); ++col) {
753 result[row, col] = lhs[row, col] + rhs[row, col];
754 }
755 }
756
757 return result;
758 }
759
765 template <SleipnirMatrixLike<Scalar> LHS, EigenMatrixLike RHS>
767 slp_assert(lhs.rows() == rhs.rows() && lhs.cols() == rhs.cols());
768
769 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), lhs.cols());
770
771 for (int row = 0; row < result.rows(); ++row) {
772 for (int col = 0; col < result.cols(); ++col) {
773 result[row, col] = lhs[row, col] + rhs[row, col];
774 }
775 }
776
777 return result;
778 }
779
785 template <SleipnirMatrixLike<Scalar> LHS, SleipnirMatrixLike<Scalar> RHS>
787 slp_assert(lhs.rows() == rhs.rows() && lhs.cols() == rhs.cols());
788
789 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), lhs.cols());
790
791 for (int row = 0; row < result.rows(); ++row) {
792 for (int col = 0; col < result.cols(); ++col) {
793 result[row, col] = lhs[row, col] + rhs[row, col];
794 }
795 }
796
797 return result;
798 }
799
805 slp_assert(rows() == rhs.rows() && cols() == rhs.cols());
806
807 for (int row = 0; row < rows(); ++row) {
808 for (int col = 0; col < cols(); ++col) {
809 (*this)[row, col] += rhs[row, col];
810 }
811 }
812
813 return *this;
814 }
815
821 slp_assert(rows() == 1 && cols() == 1);
822
823 for (int row = 0; row < rows(); ++row) {
824 for (int col = 0; col < cols(); ++col) {
825 (*this)[row, col] += rhs;
826 }
827 }
828
829 return *this;
830 }
831
837 template <EigenMatrixLike LHS, SleipnirMatrixLike<Scalar> RHS>
839 slp_assert(lhs.rows() == rhs.rows() && lhs.cols() == rhs.cols());
840
841 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), lhs.cols());
842
843 for (int row = 0; row < result.rows(); ++row) {
844 for (int col = 0; col < result.cols(); ++col) {
845 result[row, col] = lhs[row, col] - rhs[row, col];
846 }
847 }
848
849 return result;
850 }
851
857 template <SleipnirMatrixLike<Scalar> LHS, EigenMatrixLike RHS>
859 slp_assert(lhs.rows() == rhs.rows() && lhs.cols() == rhs.cols());
860
861 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), lhs.cols());
862
863 for (int row = 0; row < result.rows(); ++row) {
864 for (int col = 0; col < result.cols(); ++col) {
865 result[row, col] = lhs[row, col] - rhs[row, col];
866 }
867 }
868
869 return result;
870 }
871
877 template <SleipnirMatrixLike<Scalar> LHS, SleipnirMatrixLike<Scalar> RHS>
879 slp_assert(lhs.rows() == rhs.rows() && lhs.cols() == rhs.cols());
880
881 VariableMatrix<Scalar> result(detail::empty, lhs.rows(), lhs.cols());
882
883 for (int row = 0; row < result.rows(); ++row) {
884 for (int col = 0; col < result.cols(); ++col) {
885 result[row, col] = lhs[row, col] - rhs[row, col];
886 }
887 }
888
889 return result;
890 }
891
897 slp_assert(rows() == rhs.rows() && cols() == rhs.cols());
898
899 for (int row = 0; row < rows(); ++row) {
900 for (int col = 0; col < cols(); ++col) {
901 (*this)[row, col] -= rhs[row, col];
902 }
903 }
904
905 return *this;
906 }
907
913 slp_assert(rows() == 1 && cols() == 1);
914
915 for (int row = 0; row < rows(); ++row) {
916 for (int col = 0; col < cols(); ++col) {
917 (*this)[row, col] -= rhs;
918 }
919 }
920
921 return *this;
922 }
923
928 const SleipnirMatrixLike<Scalar> auto& lhs) {
929 VariableMatrix<Scalar> result{detail::empty, lhs.rows(), lhs.cols()};
930
931 for (int row = 0; row < result.rows(); ++row) {
932 for (int col = 0; col < result.cols(); ++col) {
933 result[row, col] = -lhs[row, col];
934 }
935 }
936
937 return result;
938 }
939
941 // NOLINTNEXTLINE (google-explicit-constructor)
942 operator Variable<Scalar>() const {
943 slp_assert(rows() == 1 && cols() == 1);
944 return (*this)[0, 0];
945 }
946
951 VariableMatrix<Scalar> result{detail::empty, cols(), rows()};
952
953 for (int row = 0; row < rows(); ++row) {
954 for (int col = 0; col < cols(); ++col) {
955 result[col, row] = (*this)[row, col];
956 }
957 }
958
959 return result;
960 }
961
965 int rows() const { return m_rows; }
966
970 int cols() const { return m_cols; }
971
977 Scalar value(int row, int col) { return (*this)[row, col].value(); }
978
983 Scalar value(int index) { return (*this)[index].value(); }
984
988 Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> value() {
989 Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> result{rows(),
990 cols()};
991
992 for (int row = 0; row < rows(); ++row) {
993 for (int col = 0; col < cols(); ++col) {
994 result[row, col] = value(row, col);
995 }
996 }
997
998 return result;
999 }
1000
1007 const {
1008 VariableMatrix<Scalar> result{detail::empty, rows(), cols()};
1009
1010 for (int row = 0; row < rows(); ++row) {
1011 for (int col = 0; col < cols(); ++col) {
1012 result[row, col] = unary_op((*this)[row, col]);
1013 }
1014 }
1015
1016 return result;
1017 }
1018
1019#ifndef DOXYGEN_SHOULD_SKIP_THIS
1020
1021 class iterator {
1022 public:
1023 using iterator_category = std::bidirectional_iterator_tag;
1024 using value_type = Variable<Scalar>;
1025 using difference_type = std::ptrdiff_t;
1026 using pointer = Variable<Scalar>*;
1027 using reference = Variable<Scalar>&;
1028
1029 constexpr iterator() noexcept = default;
1030
1031 explicit constexpr iterator(
1032 gch::small_vector<Variable<Scalar>>::iterator it) noexcept
1033 : m_it{it} {}
1034
1035 constexpr iterator& operator++() noexcept {
1036 ++m_it;
1037 return *this;
1038 }
1039
1040 constexpr iterator operator++(int) noexcept {
1041 iterator retval = *this;
1042 ++(*this);
1043 return retval;
1044 }
1045
1046 constexpr iterator& operator--() noexcept {
1047 --m_it;
1048 return *this;
1049 }
1050
1051 constexpr iterator operator--(int) noexcept {
1052 iterator retval = *this;
1053 --(*this);
1054 return retval;
1055 }
1056
1057 constexpr bool operator==(const iterator&) const noexcept = default;
1058
1059 constexpr reference operator*() const noexcept { return *m_it; }
1060
1061 private:
1062 gch::small_vector<Variable<Scalar>>::iterator m_it;
1063 };
1064
1065 class const_iterator {
1066 public:
1067 using iterator_category = std::bidirectional_iterator_tag;
1068 using value_type = Variable<Scalar>;
1069 using difference_type = std::ptrdiff_t;
1070 using pointer = Variable<Scalar>*;
1071 using const_reference = const Variable<Scalar>&;
1072
1073 constexpr const_iterator() noexcept = default;
1074
1075 explicit constexpr const_iterator(
1076 gch::small_vector<Variable<Scalar>>::const_iterator it) noexcept
1077 : m_it{it} {}
1078
1079 constexpr const_iterator& operator++() noexcept {
1080 ++m_it;
1081 return *this;
1082 }
1083
1084 constexpr const_iterator operator++(int) noexcept {
1085 const_iterator retval = *this;
1086 ++(*this);
1087 return retval;
1088 }
1089
1090 constexpr const_iterator& operator--() noexcept {
1091 --m_it;
1092 return *this;
1093 }
1094
1095 constexpr const_iterator operator--(int) noexcept {
1096 const_iterator retval = *this;
1097 --(*this);
1098 return retval;
1099 }
1100
1101 constexpr bool operator==(const const_iterator&) const noexcept = default;
1102
1103 constexpr const_reference operator*() const noexcept { return *m_it; }
1104
1105 private:
1106 gch::small_vector<Variable<Scalar>>::const_iterator m_it;
1107 };
1108
1109 using reverse_iterator = std::reverse_iterator<iterator>;
1110 using const_reverse_iterator = std::reverse_iterator<const_iterator>;
1111
1112#endif // DOXYGEN_SHOULD_SKIP_THIS
1113
1117 iterator begin() { return iterator{m_storage.begin()}; }
1118
1122 iterator end() { return iterator{m_storage.end()}; }
1123
1127 const_iterator begin() const { return const_iterator{m_storage.begin()}; }
1128
1132 const_iterator end() const { return const_iterator{m_storage.end()}; }
1133
1137 const_iterator cbegin() const { return const_iterator{m_storage.cbegin()}; }
1138
1142 const_iterator cend() const { return const_iterator{m_storage.cend()}; }
1143
1147 reverse_iterator rbegin() { return reverse_iterator{end()}; }
1148
1152 reverse_iterator rend() { return reverse_iterator{begin()}; }
1153
1157 const_reverse_iterator rbegin() const {
1158 return const_reverse_iterator{end()};
1159 }
1160
1164 const_reverse_iterator rend() const {
1165 return const_reverse_iterator{begin()};
1166 }
1167
1171 const_reverse_iterator crbegin() const {
1172 return const_reverse_iterator{cend()};
1173 }
1174
1178 const_reverse_iterator crend() const {
1179 return const_reverse_iterator{cbegin()};
1180 }
1181
1185 size_t size() const { return m_storage.size(); }
1186
1193 VariableMatrix<Scalar> result{detail::empty, rows, cols};
1194
1195 for (auto& elem : result) {
1196 elem = Scalar(0);
1197 }
1198
1199 return result;
1200 }
1201
1208 VariableMatrix<Scalar> result{detail::empty, rows, cols};
1209
1210 for (auto& elem : result) {
1211 elem = Scalar(1);
1212 }
1213
1214 return result;
1215 }
1216
1217 private:
1218 gch::small_vector<Variable<Scalar>> m_storage;
1219 int m_rows = 0;
1220 int m_cols = 0;
1221};
1222
1223template <typename Derived>
1224VariableMatrix(const Eigen::MatrixBase<Derived>&)
1225 -> VariableMatrix<typename Derived::Scalar>;
1226
1227template <typename Derived>
1228VariableMatrix(const Eigen::DiagonalBase<Derived>&)
1229 -> VariableMatrix<typename Derived::Scalar>;
1230
1237template <typename Scalar>
1238VariableMatrix<Scalar> cwise_reduce(
1239 const VariableMatrix<Scalar>& lhs, const VariableMatrix<Scalar>& rhs,
1240 function_ref<Variable<Scalar>(const Variable<Scalar>& x,
1241 const Variable<Scalar>& y)>
1242 binary_op) {
1243 slp_assert(lhs.rows() == rhs.rows() && lhs.cols() == rhs.cols());
1244
1245 VariableMatrix<Scalar> result{detail::empty, lhs.rows(), lhs.cols()};
1246
1247 for (int row = 0; row < lhs.rows(); ++row) {
1248 for (int col = 0; col < lhs.cols(); ++col) {
1249 result[row, col] = binary_op(lhs[row, col], rhs[row, col]);
1250 }
1251 }
1252
1253 return result;
1254}
1255
1265template <typename Scalar>
1266VariableMatrix<Scalar> block(
1267 std::initializer_list<std::initializer_list<VariableMatrix<Scalar>>> list) {
1268 // Get row and column counts for destination matrix
1269 int rows = 0;
1270 int cols = -1;
1271 for (const auto& row : list) {
1272 if (row.size() > 0) {
1273 rows += row.begin()->rows();
1274 }
1275
1276 // Get number of columns in this row
1277 int latest_cols = 0;
1278 for (const auto& elem : row) {
1279 // Assert the first and latest row have the same height
1280 slp_assert(row.begin()->rows() == elem.rows());
1281
1282 latest_cols += elem.cols();
1283 }
1284
1285 // If this is the first row, record the column count. Otherwise, assert the
1286 // first and latest column counts are the same.
1287 if (cols == -1) {
1288 cols = latest_cols;
1289 } else {
1290 slp_assert(cols == latest_cols);
1291 }
1292 }
1293
1294 VariableMatrix<Scalar> result{detail::empty, rows, cols};
1295
1296 int row_offset = 0;
1297 for (const auto& row : list) {
1298 int col_offset = 0;
1299 for (const auto& elem : row) {
1300 result.block(row_offset, col_offset, elem.rows(), elem.cols()) = elem;
1301 col_offset += elem.cols();
1302 }
1303 if (row.size() > 0) {
1304 row_offset += row.begin()->rows();
1305 }
1306 }
1307
1308 return result;
1309}
1310
1322template <typename Scalar>
1323VariableMatrix<Scalar> block(
1324 const std::vector<std::vector<VariableMatrix<Scalar>>>& list) {
1325 // Get row and column counts for destination matrix
1326 int rows = 0;
1327 int cols = -1;
1328 for (const auto& row : list) {
1329 if (row.size() > 0) {
1330 rows += row.begin()->rows();
1331 }
1332
1333 // Get number of columns in this row
1334 int latest_cols = 0;
1335 for (const auto& elem : row) {
1336 // Assert the first and latest row have the same height
1337 slp_assert(row.begin()->rows() == elem.rows());
1338
1339 latest_cols += elem.cols();
1340 }
1341
1342 // If this is the first row, record the column count. Otherwise, assert the
1343 // first and latest column counts are the same.
1344 if (cols == -1) {
1345 cols = latest_cols;
1346 } else {
1347 slp_assert(cols == latest_cols);
1348 }
1349 }
1350
1351 VariableMatrix<Scalar> result{detail::empty, rows, cols};
1352
1353 int row_offset = 0;
1354 for (const auto& row : list) {
1355 int col_offset = 0;
1356 for (const auto& elem : row) {
1357 result.block(row_offset, col_offset, elem.rows(), elem.cols()) = elem;
1358 col_offset += elem.cols();
1359 }
1360 if (row.size() > 0) {
1361 row_offset += row.begin()->rows();
1362 }
1363 }
1364
1365 return result;
1366}
1367
1374template <typename Scalar>
1375VariableMatrix<Scalar> solve(const VariableMatrix<Scalar>& A,
1376 const VariableMatrix<Scalar>& B) {
1377 // m x n * n x p = m x p
1378 slp_assert(A.rows() == B.rows());
1379
1380 if (A.rows() == 1 && A.cols() == 1) {
1381 // Compute optimal inverse instead of using Eigen's general solver
1382 return B[0, 0] / A[0, 0];
1383 } else if (A.rows() == 2 && A.cols() == 2) {
1384 // Compute optimal inverse instead of using Eigen's general solver
1385 //
1386 // [a b]⁻¹ ___1___ [ d −b]
1387 // [c d] = ad − bc [−c a]
1388
1389 const auto& a = A[0, 0];
1390 const auto& b = A[0, 1];
1391 const auto& c = A[1, 0];
1392 const auto& d = A[1, 1];
1393
1394 VariableMatrix adj_A{{d, -b}, {-c, a}};
1395 auto det_A = a * d - b * c;
1396 return adj_A / det_A * B;
1397 } else if (A.rows() == 3 && A.cols() == 3) {
1398 // Compute optimal inverse instead of using Eigen's general solver
1399 //
1400 // [a b c]⁻¹
1401 // [d e f]
1402 // [g h i]
1403 // 1 [ei − fh ch − bi bf − ce]
1404 // = ------------------------------------ [fg − di ai − cg cd − af]
1405 // a(ei − fh) + b(fg − di) + c(dh − eg) [dh − eg bg − ah ae − bd]
1406 //
1407 // https://www.wolframalpha.com/input?i=inverse+%7B%7Ba%2C+b%2C+c%7D%2C+%7Bd%2C+e%2C+f%7D%2C+%7Bg%2C+h%2C+i%7D%7D
1408
1409 const auto& a = A[0, 0];
1410 const auto& b = A[0, 1];
1411 const auto& c = A[0, 2];
1412 const auto& d = A[1, 0];
1413 const auto& e = A[1, 1];
1414 const auto& f = A[1, 2];
1415 const auto& g = A[2, 0];
1416 const auto& h = A[2, 1];
1417 const auto& i = A[2, 2];
1418
1419 auto ae = a * e;
1420 auto af = a * f;
1421 auto ah = a * h;
1422 auto ai = a * i;
1423 auto bd = b * d;
1424 auto bf = b * f;
1425 auto bg = b * g;
1426 auto bi = b * i;
1427 auto cd = c * d;
1428 auto ce = c * e;
1429 auto cg = c * g;
1430 auto ch = c * h;
1431 auto dh = d * h;
1432 auto di = d * i;
1433 auto eg = e * g;
1434 auto ei = e * i;
1435 auto fg = f * g;
1436 auto fh = f * h;
1437
1438 auto adj_A00 = ei - fh;
1439 auto adj_A10 = fg - di;
1440 auto adj_A20 = dh - eg;
1441
1442 VariableMatrix adj_A{{adj_A00, ch - bi, bf - ce},
1443 {adj_A10, ai - cg, cd - af},
1444 {adj_A20, bg - ah, ae - bd}};
1445 auto det_A = a * adj_A00 + b * adj_A10 + c * adj_A20;
1446 return adj_A / det_A * B;
1447 } else if (A.rows() == 4 && A.cols() == 4) {
1448 // Compute optimal inverse instead of using Eigen's general solver
1449 //
1450 // [a b c d]⁻¹
1451 // [e f g h]
1452 // [i j k l]
1453 // [m n o p]
1454 //
1455 // https://www.wolframalpha.com/input?i=inverse+%7B%7Ba%2C+b%2C+c%2C+d%7D%2C+%7Be%2C+f%2C+g%2C+h%7D%2C+%7Bi%2C+j%2C+k%2C+l%7D%2C+%7Bm%2C+n%2C+o%2C+p%7D%7D
1456
1457 const auto& a = A[0, 0];
1458 const auto& b = A[0, 1];
1459 const auto& c = A[0, 2];
1460 const auto& d = A[0, 3];
1461 const auto& e = A[1, 0];
1462 const auto& f = A[1, 1];
1463 const auto& g = A[1, 2];
1464 const auto& h = A[1, 3];
1465 const auto& i = A[2, 0];
1466 const auto& j = A[2, 1];
1467 const auto& k = A[2, 2];
1468 const auto& l = A[2, 3];
1469 const auto& m = A[3, 0];
1470 const auto& n = A[3, 1];
1471 const auto& o = A[3, 2];
1472 const auto& p = A[3, 3];
1473
1474 auto afk = a * f * k;
1475 auto afl = a * f * l;
1476 auto afo = a * f * o;
1477 auto afp = a * f * p;
1478 auto agj = a * g * j;
1479 auto agl = a * g * l;
1480 auto agn = a * g * n;
1481 auto agp = a * g * p;
1482 auto ahj = a * h * j;
1483 auto ahk = a * h * k;
1484 auto ahn = a * h * n;
1485 auto aho = a * h * o;
1486 auto ajo = a * j * o;
1487 auto ajp = a * j * p;
1488 auto akn = a * k * n;
1489 auto akp = a * k * p;
1490 auto aln = a * l * n;
1491 auto alo = a * l * o;
1492 auto bek = b * e * k;
1493 auto bel = b * e * l;
1494 auto beo = b * e * o;
1495 auto bep = b * e * p;
1496 auto bgi = b * g * i;
1497 auto bgl = b * g * l;
1498 auto bgm = b * g * m;
1499 auto bgp = b * g * p;
1500 auto bhi = b * h * i;
1501 auto bhk = b * h * k;
1502 auto bhm = b * h * m;
1503 auto bho = b * h * o;
1504 auto bio = b * i * o;
1505 auto bip = b * i * p;
1506 auto bjp = b * j * p;
1507 auto bkm = b * k * m;
1508 auto bkp = b * k * p;
1509 auto blm = b * l * m;
1510 auto blo = b * l * o;
1511 auto cej = c * e * j;
1512 auto cel = c * e * l;
1513 auto cen = c * e * n;
1514 auto cep = c * e * p;
1515 auto cfi = c * f * i;
1516 auto cfl = c * f * l;
1517 auto cfm = c * f * m;
1518 auto cfp = c * f * p;
1519 auto chi = c * h * i;
1520 auto chj = c * h * j;
1521 auto chm = c * h * m;
1522 auto chn = c * h * n;
1523 auto cin = c * i * n;
1524 auto cip = c * i * p;
1525 auto cjm = c * j * m;
1526 auto cjp = c * j * p;
1527 auto clm = c * l * m;
1528 auto cln = c * l * n;
1529 auto dej = d * e * j;
1530 auto dek = d * e * k;
1531 auto den = d * e * n;
1532 auto deo = d * e * o;
1533 auto dfi = d * f * i;
1534 auto dfk = d * f * k;
1535 auto dfm = d * f * m;
1536 auto dfo = d * f * o;
1537 auto dgi = d * g * i;
1538 auto dgj = d * g * j;
1539 auto dgm = d * g * m;
1540 auto dgn = d * g * n;
1541 auto din = d * i * n;
1542 auto dio = d * i * o;
1543 auto djm = d * j * m;
1544 auto djo = d * j * o;
1545 auto dkm = d * k * m;
1546 auto dkn = d * k * n;
1547 auto ejo = e * j * o;
1548 auto ejp = e * j * p;
1549 auto ekn = e * k * n;
1550 auto ekp = e * k * p;
1551 auto eln = e * l * n;
1552 auto elo = e * l * o;
1553 auto fio = f * i * o;
1554 auto fip = f * i * p;
1555 auto fkm = f * k * m;
1556 auto fkp = f * k * p;
1557 auto flm = f * l * m;
1558 auto flo = f * l * o;
1559 auto gin = g * i * n;
1560 auto gip = g * i * p;
1561 auto gjm = g * j * m;
1562 auto gjp = g * j * p;
1563 auto glm = g * l * m;
1564 auto gln = g * l * n;
1565 auto hin = h * i * n;
1566 auto hio = h * i * o;
1567 auto hjm = h * j * m;
1568 auto hjo = h * j * o;
1569 auto hkm = h * k * m;
1570 auto hkn = h * k * n;
1571
1572 auto adj_A00 = fkp - flo - gjp + gln + hjo - hkn;
1573 auto adj_A01 = -bkp + blo + cjp - cln - djo + dkn;
1574 auto adj_A02 = bgp - bho - cfp + chn + dfo - dgn;
1575 auto adj_A03 = -bgl + bhk + cfl - chj - dfk + dgj;
1576 auto adj_A10 = -ekp + elo + gip - glm - hio + hkm;
1577 auto adj_A11 = akp - alo - cip + clm + dio - dkm;
1578 auto adj_A12 = -agp + aho + cep - chm - deo + dgm;
1579 auto adj_A13 = agl - ahk - cel + chi + dek - dgi;
1580 auto adj_A20 = ejp - eln - fip + flm + hin - hjm;
1581 auto adj_A21 = -ajp + aln + bip - blm - din + djm;
1582 auto adj_A22 = afp - ahn - bep + bhm + den - dfm;
1583 auto adj_A23 = -afl + ahj + bel - bhi - dej + dfi;
1584 auto adj_A30 = -ejo + ekn + fio - fkm - gin + gjm;
1585 // NOLINTNEXTLINE(build/include_what_you_use)
1586 auto adj_A31 = ajo - akn - bio + bkm + cin - cjm;
1587 auto adj_A32 = -afo + agn + beo - bgm - cen + cfm;
1588 auto adj_A33 = afk - agj - bek + bgi + cej - cfi;
1589
1590 VariableMatrix adj_A{{adj_A00, adj_A01, adj_A02, adj_A03},
1591 {adj_A10, adj_A11, adj_A12, adj_A13},
1592 {adj_A20, adj_A21, adj_A22, adj_A23},
1593 {adj_A30, adj_A31, adj_A32, adj_A33}};
1594 auto det_A = a * adj_A00 + b * adj_A10 + c * adj_A20 + d * adj_A30;
1595 return adj_A / det_A * B;
1596 } else {
1597 using MatrixXv =
1598 Eigen::Matrix<Variable<Scalar>, Eigen::Dynamic, Eigen::Dynamic>;
1599
1600 MatrixXv eigen_A{A.rows(), A.cols()};
1601 for (int row = 0; row < A.rows(); ++row) {
1602 for (int col = 0; col < A.cols(); ++col) {
1603 eigen_A[row, col] = A[row, col];
1604 }
1605 }
1606
1607 MatrixXv eigen_B{B.rows(), B.cols()};
1608 for (int row = 0; row < B.rows(); ++row) {
1609 for (int col = 0; col < B.cols(); ++col) {
1610 eigen_B[row, col] = B[row, col];
1611 }
1612 }
1613
1614 MatrixXv eigen_X = eigen_A.householderQr().solve(eigen_B);
1615
1616 VariableMatrix<Scalar> X{detail::empty, A.cols(), B.cols()};
1617 for (int row = 0; row < X.rows(); ++row) {
1618 for (int col = 0; col < X.cols(); ++col) {
1619 X[row, col] = eigen_X[row, col];
1620 }
1621 }
1622
1623 return X;
1624 }
1625}
1626
1627extern template SLEIPNIR_DLLEXPORT VariableMatrix<double> solve(
1628 const VariableMatrix<double>& A, const VariableMatrix<double>& B);
1629
1630} // namespace slp
Definition intrusive_shared_ptr.hpp:27
Definition sleipnir_base.hpp:9
Represents a sequence of elements in an iterable object.
Definition slice.hpp:25
Definition variable_block.hpp:24
Definition variable_matrix.hpp:33
const_reverse_iterator crbegin() const
Definition variable_matrix.hpp:1171
VariableMatrix(std::initializer_list< std::initializer_list< Variable< Scalar > > > list)
Definition variable_matrix.hpp:73
const_reverse_iterator crend() const
Definition variable_matrix.hpp:1178
Scalar_ Scalar
Scalar type alias.
Definition variable_matrix.hpp:36
iterator end()
Definition variable_matrix.hpp:1122
const Variable< Scalar > & operator[](int row, int col) const
Definition variable_matrix.hpp:317
VariableMatrix(Variable< Scalar > &&variable)
Definition variable_matrix.hpp:194
const VariableBlock< const VariableMatrix > block(int row_offset, int col_offset, int block_rows, int block_cols) const
Definition variable_matrix.hpp:364
size_t size() const
Definition variable_matrix.hpp:1185
static VariableMatrix< Scalar > zero(int rows, int cols)
Definition variable_matrix.hpp:1192
VariableBlock< VariableMatrix > block(int row_offset, int col_offset, int block_rows, int block_cols)
Definition variable_matrix.hpp:348
VariableMatrix(const std::vector< std::vector< Scalar > > &list)
Definition variable_matrix.hpp:100
VariableMatrix & operator=(ScalarLike auto value)
Definition variable_matrix.hpp:278
VariableMatrix & operator-=(const MatrixLike auto &rhs)
Definition variable_matrix.hpp:896
Variable< Scalar > & operator[](int index)
Definition variable_matrix.hpp:327
VariableMatrix & operator=(const Eigen::MatrixBase< Derived > &values)
Definition variable_matrix.hpp:260
friend VariableMatrix< Scalar > operator-(const LHS &lhs, const RHS &rhs)
Definition variable_matrix.hpp:838
VariableMatrix & operator-=(const ScalarLike auto &rhs)
Definition variable_matrix.hpp:912
friend VariableMatrix< Scalar > operator*(const Variable< Scalar > &lhs, const RHS &rhs)
Definition variable_matrix.hpp:605
VariableMatrix & operator*=(const ScalarLike auto &rhs)
Definition variable_matrix.hpp:659
Scalar value(int index)
Definition variable_matrix.hpp:983
VariableMatrix(const VariableBlock< VariableMatrix > &values)
Definition variable_matrix.hpp:202
VariableMatrix(int rows, int cols)
Definition variable_matrix.hpp:51
const_iterator cbegin() const
Definition variable_matrix.hpp:1137
const_iterator begin() const
Definition variable_matrix.hpp:1127
const_iterator cend() const
Definition variable_matrix.hpp:1142
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > value()
Definition variable_matrix.hpp:988
const VariableBlock< const VariableMatrix > operator[](Slice row_slice, Slice col_slice) const
Definition variable_matrix.hpp:392
friend VariableMatrix< Scalar > operator/(const LHS &lhs, const Variable< Scalar > &rhs)
Definition variable_matrix.hpp:675
VariableMatrix(int rows)
Definition variable_matrix.hpp:45
VariableMatrix(std::span< const Variable< Scalar > > values, int rows, int cols)
Definition variable_matrix.hpp:244
VariableBlock< VariableMatrix > col(int col)
Definition variable_matrix.hpp:482
VariableMatrix(const Eigen::MatrixBase< Derived > &values)
Definition variable_matrix.hpp:151
const VariableBlock< const VariableMatrix > operator[](Slice row_slice, int row_slice_length, Slice col_slice, int col_slice_length) const
Definition variable_matrix.hpp:428
const VariableBlock< const VariableMatrix > col(int col) const
Definition variable_matrix.hpp:491
VariableMatrix & operator+=(const ScalarLike auto &rhs)
Definition variable_matrix.hpp:820
friend VariableMatrix< Scalar > operator*(const LHS &lhs, const Variable< Scalar > &rhs)
Definition variable_matrix.hpp:570
VariableBlock< VariableMatrix > segment(int offset, int length)
Definition variable_matrix.hpp:440
Scalar value(int row, int col)
Definition variable_matrix.hpp:977
const VariableBlock< const VariableMatrix > row(int row) const
Definition variable_matrix.hpp:473
VariableMatrix()=default
Constructs an empty VariableMatrix.
static VariableMatrix< Scalar > ones(int rows, int cols)
Definition variable_matrix.hpp:1207
VariableMatrix & operator+=(const MatrixLike auto &rhs)
Definition variable_matrix.hpp:804
VariableMatrix(std::span< const Variable< Scalar > > values)
Definition variable_matrix.hpp:229
friend VariableMatrix< Scalar > operator-(const SleipnirMatrixLike< Scalar > auto &lhs)
Definition variable_matrix.hpp:927
VariableMatrix< Scalar > T() const
Definition variable_matrix.hpp:950
VariableMatrix< Scalar > cwise_transform(function_ref< Variable< Scalar >(const Variable< Scalar > &x)> unary_op) const
Definition variable_matrix.hpp:1005
friend VariableMatrix< Scalar > operator*(const LHS &lhs, const RHS &rhs)
Definition variable_matrix.hpp:501
const_reverse_iterator rbegin() const
Definition variable_matrix.hpp:1157
VariableMatrix(const Variable< Scalar > &variable)
Definition variable_matrix.hpp:186
reverse_iterator rend()
Definition variable_matrix.hpp:1152
VariableMatrix(const VariableBlock< const VariableMatrix > &values)
Definition variable_matrix.hpp:216
int rows() const
Definition variable_matrix.hpp:965
VariableMatrix(const Eigen::DiagonalBase< Derived > &values)
Definition variable_matrix.hpp:167
VariableBlock< VariableMatrix > operator[](Slice row_slice, Slice col_slice)
Definition variable_matrix.hpp:380
const_reverse_iterator rend() const
Definition variable_matrix.hpp:1164
VariableMatrix & operator*=(const MatrixLike auto &rhs)
Definition variable_matrix.hpp:639
VariableMatrix(detail::empty_t, int rows, int cols)
Definition variable_matrix.hpp:62
const_iterator end() const
Definition variable_matrix.hpp:1132
void set_value(const Eigen::MatrixBase< Derived > &values)
Definition variable_matrix.hpp:291
Variable< Scalar > & operator[](int row, int col)
Definition variable_matrix.hpp:306
iterator begin()
Definition variable_matrix.hpp:1117
const VariableBlock< const VariableMatrix > segment(int offset, int length) const
Definition variable_matrix.hpp:452
reverse_iterator rbegin()
Definition variable_matrix.hpp:1147
VariableBlock< VariableMatrix > row(int row)
Definition variable_matrix.hpp:464
friend VariableMatrix< Scalar > operator+(const LHS &lhs, const RHS &rhs)
Definition variable_matrix.hpp:746
VariableMatrix(const std::vector< std::vector< Variable< Scalar > > > &list)
Definition variable_matrix.hpp:126
VariableBlock< VariableMatrix > operator[](Slice row_slice, int row_slice_length, Slice col_slice, int col_slice_length)
Definition variable_matrix.hpp:410
VariableMatrix & operator/=(const ScalarLike auto &rhs)
Definition variable_matrix.hpp:730
const Variable< Scalar > & operator[](int index) const
Definition variable_matrix.hpp:336
int cols() const
Definition variable_matrix.hpp:970
Definition variable.hpp:47
Definition function_ref.hpp:13
Definition concepts.hpp:18
Definition concepts.hpp:24
Definition concepts.hpp:33
Definition concepts.hpp:38
Type tag used to designate an uninitialized VariableMatrix.
Definition empty.hpp:8