codac 1.5.6
Loading...
Searching...
No Matches
codac2_matrices.h
Go to the documentation of this file.
1
19
20#pragma once
21
22#include <type_traits>
23#include "codac2_Interval.h"
25#include "codac2_assert.h"
26
27namespace Eigen
28{
29 template<int R,int C>
30 concept IsVectorOrRow = (C == 1 || R == 1);
31
32 template<typename Scalar>
33 concept IsIntervalDomain = std::is_same_v<Scalar,codac2::Interval>;
34}
35
36#define EIGEN_MATRIXBASE_PLUGIN "codac2_MatrixBase_addons_include.h"
37#define EIGEN_MATRIX_PLUGIN "codac2_Matrix_addons_include.h"
38
39#ifndef EIGEN_NO_DEBUG
40/* Disables Eigen's assertions if defined.
41 * Not defined by default, unless the NDEBUG macro is defined
42 * (this is a standard C++ macro which disables all asserts).
43 * https://eigen.tuxfamily.org/dox/TopicPreprocessorDirectives.html
44 */
45//#define EIGEN_NO_DEBUG // uncomment to disable Eigen's assertions
46#endif
47
48#include <Eigen/Core>
49#include <Eigen/Dense>
50
51namespace Eigen
52{
53 template<> struct NumTraits<codac2::Interval>
54 : NumTraits<double> // permits to get the epsilon, dummy_precision, lowest, highest functions
55 {
56 typedef codac2::Interval Real;
57 typedef codac2::Interval NonInteger;
58 typedef codac2::Interval Nested;
59 typedef codac2::Interval Scalar;
60 typedef double RealScalar;
61
62 enum {
63 IsComplex = 0,
64 IsInteger = 0,
65 IsSigned = 1,
66 RequireInitialization = 1,
67 ReadCost = 1,
68 AddCost = 3,
69 MulCost = 3
70 };
71 };
72
73 template<typename BinOp>
74 struct ScalarBinaryOpTraits<codac2::Interval,double,BinOp>
75 { typedef codac2::Interval ReturnType; };
76
77 template<typename BinOp>
78 struct ScalarBinaryOpTraits<double,codac2::Interval,BinOp>
79 { typedef codac2::Interval ReturnType; };
80}
81
82namespace codac2
83{
84 using Eigen::Dynamic;
85
86 inline const Interval& conj(const Interval& x) { return x; }
87 inline const Interval& real(const Interval& x) { return x; }
88 inline Interval imag(const Interval&) { return 0.; }
89 //inline Interval abs(const Interval& x) { return codac2::abs(x); }
90 inline Interval abs2(const Interval& x) { return codac2::sqr(x); }
91
92 template<typename Scalar,int R,int C>
93 using Mat = Eigen::Matrix<Scalar,R,C>;
94
95 template<int R,int C>
96 struct is_interval_based<Eigen::Matrix<double,R,C>> : std::false_type {};
97
98 template<int R,int C>
99 struct is_interval_based<Eigen::Matrix<Interval,R,C>> : std::true_type {};
100
101 template<typename T,int R,int C>
102 struct is_ctc<Eigen::Matrix<T,R,C>> : std::false_type {};
103
104 template<typename T,int R,int C>
105 struct is_sep<Eigen::Matrix<T,R,C>> : std::false_type {};
106}
107
108namespace codac2
109{
110 template<typename OtherDerived>
111 inline auto abs(const Eigen::MatrixBase<OtherDerived>& x)
112 {
113 using M = Eigen::MatrixBase<OtherDerived>;
114 Eigen::Matrix<typename M::Scalar,M::RowsAtCompileTime,M::ColsAtCompileTime> a(x.rows(),x.cols());
115
116 for(Index i = 0 ; i < x.rows() ; i++)
117 for(Index j = 0 ; j < x.cols() ; j++)
118 {
119 if constexpr(std::is_same_v<typename M::Scalar,double>)
120 a(i,j) = fabs(x(i,j));
121 else
122 a(i,j) = abs(x(i,j));
123 }
124
125 return a;
126 }
127
128 template<typename OtherDerived>
129 requires (!Eigen::IsIntervalDomain<typename OtherDerived::Scalar>)
130 inline auto floor(const Eigen::MatrixBase<OtherDerived>& x)
131 {
132 return x.array().floor().matrix();
133 }
134
135 template<typename OtherDerived>
136 requires (!Eigen::IsIntervalDomain<typename OtherDerived::Scalar>)
137 inline auto ceil(const Eigen::MatrixBase<OtherDerived>& x)
138 {
139 return x.array().ceil().matrix();
140 }
141
142 template<typename OtherDerived>
143 requires (!Eigen::IsIntervalDomain<typename OtherDerived::Scalar>)
144 inline auto round(const Eigen::MatrixBase<OtherDerived>& x)
145 {
146 return x.array().round().matrix();
147 }
148
149 inline Eigen::IOFormat codac_row_fmt()
150 {
151 return Eigen::IOFormat(Eigen::StreamPrecision, Eigen::DontAlignCols, " ", "", "", "", "[ ", " ]");
152 }
153
154 inline Eigen::IOFormat codac_vector_fmt()
155 {
156 return Eigen::IOFormat(Eigen::StreamPrecision, Eigen::DontAlignCols, "", " ; ", "", "", "[ ", " ]");
157 }
158
159 inline Eigen::IOFormat codac_matrix_fmt()
160 {
161 return Eigen::IOFormat(Eigen::StreamPrecision, 0, " , ", "\n", "[ ", " ]", "[", "]");
162 }
163}
Interval class, for representing closed and connected subsets of .
Definition codac2_Interval.h:62
Interval ceil(const Interval &x)
Returns ceil of .
Definition codac2_Interval_operations_impl.h:294
Interval floor(const Interval &x)
Returns floor of .
Definition codac2_Interval_operations_impl.h:289
Interval sqr(const Interval &x)
Returns .
Definition codac2_Interval_operations_impl.h:21
Interval abs(const Interval &x)
Returns .
Definition codac2_Interval_operations_impl.h:264