codac 1.5.6
Loading...
Searching...
No Matches
codac2_arith_div.h
Go to the documentation of this file.
1
9
10#pragma once
11
12#include "codac2_Interval.h"
15#include "codac2_AnalyticType.h"
17#include "codac2_arith_mul.h"
18
19namespace codac2
20{
21 struct DivOp
22 {
23 template<typename X1,typename X2>
24 static std::string str(const X1& x1, const X2& x2)
25 {
26 return x1->str(!x1->is_str_leaf()) + "/" + x2->str(!x2->is_str_leaf());
27 }
28
29 template<typename X1,typename X2>
30 static std::pair<Index,Index> output_shape(const X1& s1, [[maybe_unused]] const X2& s2)
31 {
32 [[maybe_unused]] auto shape2 = s2->output_shape();
33 assert(shape2.first == 1 && shape2.second == 1);
34 return s1->output_shape();
35 }
36
37 static Interval fwd(const Interval& x1, const Interval& x2);
38 static ScalarType fwd_natural(const ScalarType& x1, const ScalarType& x2);
39 static ScalarType fwd_centered(const ScalarType& x1, const ScalarType& x2);
40 static void bwd(const Interval& y, Interval& x1, Interval& x2);
41
42 static IntervalVector fwd(const IntervalVector& x1, const Interval& x2);
43 static VectorType fwd_natural(const VectorType& x1, const ScalarType& x2);
44 static VectorType fwd_centered(const VectorType& x1, const ScalarType& x2);
45 static void bwd(const IntervalVector& y, IntervalVector& x1, Interval& x2);
46
47 static IntervalMatrix fwd(const IntervalMatrix& x1, const Interval& x2);
48 static MatrixType fwd_natural(const MatrixType& x1, const ScalarType& x2);
49 static MatrixType fwd_centered(const MatrixType& x1, const ScalarType& x2);
50 static void bwd(const IntervalMatrix& y, IntervalMatrix& x1, Interval& x2);
51 };
52
53 // operator/
54 // The following functions can be used to build analytic expressions.
55
56 inline ScalarExpr
57 operator/(const ScalarExpr& x1, const ScalarExpr& x2)
58 {
59 return { std::make_shared<AnalyticOperationExpr<DivOp,ScalarType,ScalarType,ScalarType>>(x1,x2) };
60 }
61
62 inline VectorExpr
63 operator/(const VectorExpr& x1, const ScalarExpr& x2)
64 {
65 return { std::make_shared<AnalyticOperationExpr<DivOp,VectorType,VectorType,ScalarType>>(x1,x2) };
66 }
67
68 inline MatrixExpr
69 operator/(const MatrixExpr& x1, const ScalarExpr& x2)
70 {
71 return { std::make_shared<AnalyticOperationExpr<DivOp,MatrixType,MatrixType,ScalarType>>(x1,x2) };
72 }
73
74 // Inline functions
75
76 inline Interval DivOp::fwd(const Interval& x1, const Interval& x2)
77 {
78 return x1 / x2;
79 }
80
81 inline ScalarType DivOp::fwd_natural(const ScalarType& x1, const ScalarType& x2)
82 {
83 return {
84 fwd(x1.a, x2.a),
85 x1.def_domain && x2.def_domain && x2.a != 0. // def domain of the derivative of div
86 };
87 }
88
89 inline ScalarType DivOp::fwd_centered(const ScalarType& x1, const ScalarType& x2)
90 {
91 if(centered_form_not_available_for_args(x1,x2))
92 return fwd_natural(x1,x2);
93
94 assert(x1.da.size() == x2.da.size());
95
96 IntervalMatrix d(1,x1.da.size());
97 for(Index i = 0 ; i < d.size() ; i++)
98 d(0,i) = (x1.da(0,i)*x2.a-x1.a*x2.da(0,i))/sqr(x2.a);
99
100 return {
101 fwd(x1.m, x2.m),
102 fwd(x1.a, x2.a),
103 d,
104 x1.def_domain && x2.def_domain && x2.a != 0. // def domain of the derivative of div
105 };
106 }
107
108 inline void DivOp::bwd(const Interval& y, Interval& x1, Interval& x2)
109 {
110 if((x1 &= y*x2).is_empty())
111 x2.set_empty();
112
113 else
114 {
115 Interval tmp = y;
116 MulOp::bwd(x1, tmp, x2);
117 if(x2.is_empty())
118 x1.set_empty();
119 }
120 }
121
122 inline IntervalVector DivOp::fwd(const IntervalVector& x1, const Interval& x2)
123 {
124 return x1 / x2;
125 }
126
127 inline VectorType DivOp::fwd_natural(const VectorType& x1, const ScalarType& x2)
128 {
129 return {
130 fwd(x1.a, x2.a),
131 x1.def_domain && x2.def_domain && x2.a != 0. // def domain of the derivative of div
132 };
133 }
134
135 inline VectorType DivOp::fwd_centered(const VectorType& x1, const ScalarType& x2)
136 {
137 if(centered_form_not_available_for_args(x1,x2))
138 return fwd_natural(x1,x2);
139
140 assert(x2.da.rows() == 1);
141 assert(x1.da.cols() == x2.da.cols());
142
143 IntervalMatrix d(x1.da.rows(),x1.da.cols());
144 for(Index j = 0 ; j < d.cols() ; j++)
145 for(Index i = 0 ; i < d.rows() ; i++)
146 d(i,j) = x1.da(i,j)/x2.a - x1.a[i]*x2.da(0,j)/sqr(x2.a);
147
148 return {
149 fwd(x1.m, x2.m),
150 fwd(x1.a, x2.a),
151 d,
152 x1.def_domain && x2.def_domain && x2.a != 0. // def domain of the derivative of div
153 };
154 }
155
156 inline void DivOp::bwd(const IntervalVector& y, IntervalVector& x1, Interval& x2)
157 {
158 assert(x1.size() == y.size());
159 for(Index i = 0 ; i < x1.size() ; i++)
160 DivOp::bwd(y[i], x1[i], x2);
161 }
162
163 inline IntervalMatrix DivOp::fwd(const IntervalMatrix& x1, const Interval& x2)
164 {
165 return x1 / x2;
166 }
167
168 inline MatrixType DivOp::fwd_natural(const MatrixType& x1, const ScalarType& x2)
169 {
170 return {
171 fwd(x1.a, x2.a),
172 x1.def_domain && x2.def_domain && x2.a != 0. // def domain of the derivative of div
173 };
174 }
175
176 inline MatrixType DivOp::fwd_centered(const MatrixType& x1, const ScalarType& x2)
177 {
178 if(centered_form_not_available_for_args(x1,x2))
179 return fwd_natural(x1,x2);
180
181 assert(x2.da.rows() == 1);
182 assert(x1.da.cols() == x2.da.cols());
183
184 IntervalMatrix d(x1.da.rows(),x1.da.cols());
185 for(Index j = 0 ; j < d.cols() ; j++)
186 for(Index i = 0 ; i < d.rows() ; i++)
187 d(i,j) = x1.da(i,j)/x2.a -
188 x1.a.reshaped<Eigen::ColMajor>()[i]*x2.da(0,j)/sqr(x2.a);
189
190 return {
191 fwd(x1.m, x2.m),
192 fwd(x1.a, x2.a),
193 d,
194 x1.def_domain && x2.def_domain && x2.a != 0. // def domain of the derivative of div
195 };
196 }
197
198 inline void DivOp::bwd(const IntervalMatrix& y, IntervalMatrix& x1, Interval& x2)
199 {
200 assert(x1.cols() == y.cols() && x1.rows() == y.rows());
201 for(Index j = 0 ; j < x1.cols() ; j++)
202 for(Index i = 0 ; i < x1.rows() ; i++)
203 DivOp::bwd(y(i,j), x1(i,j), x2);
204 }
205}
Interval class, for representing closed and connected subsets of .
Definition codac2_Interval.h:62
Interval operator/(const Interval &x, double y)
Returns with .
Definition codac2_Interval_impl.h:451
Interval sqr(const Interval &x)
Returns .
Definition codac2_Interval_operations_impl.h:21