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