378 #include "polymodel.hpp"
380 #include "mclapack.hpp"
382 #undef MC__TMODEL_DEBUG
383 #undef MC__TMODEL_DEBUG_SCALE
384 #define MC__TMODEL_CHECK
385 #undef MC__TMODEL_CHECK_PMODEL
386 #undef MC__TVAR_DISPLAY_EXT
387 #undef MC__TVAR_DEBUG_EXP
388 #undef MC__TVAR_DEBUG_BERSTEIN
389 #undef MC__TVAR_HYBRID_EIGEN
390 #define MC__TMODEL_TIGHT_REMAINDER
395 template <
typename T>
class TVar;
405 template <
typename T>
409 friend class TVar<T>;
410 template <
typename U>
friend class TModel;
412 template <
typename U>
friend TVar<U> inv
414 template <
typename U>
friend TVar<U> sqrt
416 template <
typename U>
friend TVar<U> log
418 template <
typename U>
friend TVar<U> exp
420 template <
typename U>
friend TVar<U> pow
430 (
const unsigned int nvar,
const unsigned int nord )
444 {
return _prodmon; };
448 {
return _refpoint; }
484 return "mc::TModel\t Division by zero scalar";
486 return "mc::TModel\t Inverse operation with zero in range";
488 return "mc::TModel\t Log operation with non-positive numbers in range";
490 return "mc::TModel\t Square-root operation with negative numbers in range";
492 return "mc::TModel\t Sine/Cosine inverse operation with range outside [-1,1]";
494 return "mc::TModel\t Range bounder with eigenvalue decomposition failed";
496 return "mc::TModel\t Bernstein remainder evalution failed";
498 return "mc::TModel\t Taylor variable initialization failed";
500 return "mc::TModel\t Inconsistent bounds with template parameter arithmetic";
502 return "mc::TModel\t Operation between Taylor variables in different Taylor model environment not allowed";
504 return "mc::TModel\t Feature not yet implemented in mc::TModel class";
506 return "mc::TModel\t Undocumented error";
540 template <
typename U>
Options&
operator =
592 unsigned int **_prodmon;
624 (
const unsigned int ivar,
const T&X,
const double Xref,
635 (
const unsigned int ivar,
double*coef )
const;
638 typedef double (puniv)
639 (
const double x,
const double*rusr,
const int*iusr );
642 typedef double (punivopt)
643 (
const double x,
const double*rusr,
const int*iusr, puniv df,
644 const T&I,
const std::pair<unsigned int,const double*>&bern );
648 (
const T&x,
const double*rusr,
const int*iusr );
652 (
const TVar<T>&TV, puniv f, puniv df, punivext If,
653 punivext d2If,
const double*rusr,
const int*iusr );
656 static double _gap_bernstein
657 (
const double x,
const double*rusr,
const int*iusr, puniv df,
658 const T&I,
const std::pair<unsigned int,const double*>&bern );
661 static double _dgap_bernstein
662 (
const double x,
const double*rusr,
const int*iusr, puniv df,
663 const T&I,
const std::pair<unsigned int,const double*>&bern );
667 (
const double xL,
const double xU, punivopt fopt,
const double*rusr,
668 const int*iusr, puniv df,
const T&Ix,
669 const std::pair<unsigned int,const double*>&bern );
672 double _goldsect_iter
673 (
const bool init,
const double a,
const double fa,
const double b,
674 const double fb,
const double c,
const double fc, punivopt fopt,
675 const double*rusr,
const int*iusr, puniv df,
const T&Ix,
676 const std::pair<unsigned int,const double*>&bern );
680 (
const TVar<T>&TV,
const int n );
716 = {
"NAIVE",
"LSB",
"EIGEN",
"BERNSTEIN",
"HYBRID" };
724 template <
typename T>
725 class TVar:
public PolyVar<T>
728 template <
typename U>
friend class TVar;
729 template <
typename U>
friend class TModel;
731 template <
typename U>
friend TVar<U>
operator+
733 template <
typename U,
typename V>
friend TVar<U>
operator+
734 (
const TVar<U>&,
const TVar<V>& );
735 template <
typename U>
friend TVar<U>
operator+
736 (
const TVar<U>&,
const U& );
737 template <
typename U>
friend TVar<U>
operator+
738 (
const U&,
const TVar<U>& );
739 template <
typename U>
friend TVar<U>
operator+
740 (
const double,
const TVar<U>& );
741 template <
typename U>
friend TVar<U>
operator+
742 (
const TVar<U>&,
const double );
743 template <
typename U>
friend TVar<U>
operator-
745 template <
typename U,
typename V>
friend TVar<U>
operator-
746 (
const TVar<U>&,
const TVar<V>& );
747 template <
typename U>
friend TVar<U>
operator-
748 (
const TVar<U>&,
const U& );
749 template <
typename U>
friend TVar<U>
operator-
750 (
const U&,
const TVar<U>& );
751 template <
typename U>
friend TVar<U>
operator-
752 (
const double,
const TVar<U>& );
753 template <
typename U>
friend TVar<U>
operator-
754 (
const TVar<U>&,
const double );
755 template <
typename U>
friend TVar<U>
operator*
756 (
const TVar<U>&,
const TVar<U>& );
757 template <
typename U>
friend TVar<U>
operator*
758 (
const double,
const TVar<U>& );
759 template <
typename U>
friend TVar<U>
operator*
760 (
const TVar<U>&,
const double );
761 template <
typename U>
friend TVar<U>
operator*
762 (
const U&,
const TVar<U>& );
763 template <
typename U>
friend TVar<U>
operator*
764 (
const TVar<U>&,
const U& );
765 template <
typename U>
friend TVar<U>
operator/
766 (
const TVar<U>&,
const TVar<U>& );
767 template <
typename U>
friend TVar<U>
operator/
768 (
const double,
const TVar<U>& );
769 template <
typename U>
friend TVar<U>
operator/
770 (
const TVar<U>&,
const double );
771 template <
typename U>
friend std::ostream&
operator<<
772 ( std::ostream&,
const TVar<U>& );
774 template <
typename U>
friend TVar<U> inv
776 template <
typename U>
friend TVar<U> sqr
778 template <
typename U>
friend TVar<U> sqrt
780 template <
typename U>
friend TVar<U> exp
782 template <
typename U>
friend TVar<U> log
784 template <
typename U>
friend TVar<U> xlog
786 template <
typename U>
friend TVar<U> pow
787 (
const TVar<U>&,
const int );
788 template <
typename U>
friend TVar<U> pow
789 (
const TVar<U>&,
const double );
790 template <
typename U>
friend TVar<U> pow
791 (
const double,
const TVar<U>& );
792 template <
typename U>
friend TVar<U> pow
793 (
const TVar<U>&,
const TVar<U>& );
794 template <
typename U>
friend TVar<U> monomial
795 (
const unsigned int,
const TVar<U>*,
const int* );
796 template <
typename U>
friend TVar<U> cos
798 template <
typename U>
friend TVar<U> sin
800 template <
typename U>
friend TVar<U> tan
802 template <
typename U>
friend TVar<U> acos
804 template <
typename U>
friend TVar<U> asin
806 template <
typename U>
friend TVar<U> atan
808 template <
typename U>
friend TVar<U> hull
809 (
const TVar<U>&,
const TVar<U>& );
810 template <
typename U>
friend bool inter
811 ( TVar<U>&,
const TVar<U>&,
const TVar<U>& );
814 using PolyVar<T>::_coefmon;
815 using PolyVar<T>::_bndord;
816 using PolyVar<T>::_bndord_uptd;
817 using PolyVar<T>::_bndrem;
818 using PolyVar<T>::_bndpol;
819 using PolyVar<T>::_set_bndpol;
820 using PolyVar<T>::_unset_bndpol;
821 using PolyVar<T>::nvar;
822 using PolyVar<T>::nord;
823 using PolyVar<T>::nmon;
824 using PolyVar<T>::_posord;
825 using PolyVar<T>::_expmon;
826 using PolyVar<T>::_loc_expmon;
827 using PolyVar<T>::_get_binom;
828 using PolyVar<T>::_resize;
829 using PolyVar<T>::_set;
830 using PolyVar<T>::set;
831 using PolyVar<T>::get;
832 using PolyVar<T>::bound;
833 using PolyVar<T>::bndord;
834 using PolyVar<T>::_center;
849 {
return _TM->_TV; };
852 unsigned int _prodmon
853 (
const unsigned int imon,
const unsigned int jmon )
const
854 {
return _TM->prodmon()[imon][jmon]; };
858 (
const unsigned int imon )
const
859 {
return _TM->bndmon()[imon]; };
863 (
const unsigned int ivar )
const
864 {
return _TM->_refpoint[ivar]; };
868 (
const unsigned int ivar )
const
869 {
return _TM->_scaling[ivar]; };
875 TModel<T>*
env()
const
881 (
const double d=0. );
889 ( TModel<T>*TM,
const unsigned int ix,
const T&X,
const double Xref );
893 ( TModel<T>*TM,
const unsigned int ix,
const T&X );
896 template <
typename U> TVar
897 ( TModel<T>*&TM,
const TVar<U>&TV );
900 template <
typename U> TVar
901 ( TModel<T>*&TM,
const TVar<U>&TV,
const T& (U::*method)()
const );
904 template <
typename U> TVar
905 ( TModel<T>*&TM,
const TVar<U>&TV, T (*method)(
const U& ) );
912 #ifdef MC__TMODEL_CHECK_PMODEL
933 (
const unsigned int ivar,
const T&X,
const double Xref );
941 (
TModel<T>*TM,
const unsigned int ix,
const T&X,
const double Xref )
942 {
set( TM ); _TM = TM; _set( ix, X, Xref );
return *
this; }
956 (
const double*x )
const;
960 (
const double*x )
const
967 {
TVar<T> var = *
this; *(var._bndrem) = 0.;
return var; }
985 (
const bool reset=
false );
992 (
const unsigned int ivar,
const bool reset=
false );
1001 template <
typename U>
TVar<T>&
operator +=
1003 template <
typename U>
TVar<T>&
operator +=
1007 template <
typename U>
TVar<T>&
operator -=
1009 template <
typename U>
TVar<T>&
operator -=
1026 void _update_bndord()
const;
1030 (
const int type )
const;
1035 {
return _polybound( _TM? _TM->options.BOUNDER_TYPE: 0 ); }
1038 T _polybound_naive()
const;
1041 T _polybound_LSB()
const;
1044 T _polybound_eigen()
const;
1047 T _polybound_bernstein()
const;
1050 double _coef_bernstein
1051 (
const double*
coefmon,
const unsigned int*jexp,
1052 const unsigned int jmon,
const unsigned int maxord )
const;
1057 template <
typename T>
inline void
1062 _bndpow =
new T*[_nvar];
1063 for(
unsigned int i=0; i<_nvar; i++ ) _bndpow[i] = 0;
1064 _bndmon =
new T[_nmon];
1066 _refpoint =
new double[_nvar];
1067 _scaling =
new double[_nvar];
1068 _cbern =
new double[_nord+1];
1070 _TV =
new TVar<T>( this );
1073 template <
typename T>
inline void
1076 for(
unsigned int i=0; i<_nvar; i++ ){
1077 delete[] _bndpow[i];
1082 template <
typename T>
inline void
1083 TModel<T>::_cleanup()
1085 for(
unsigned int i=0; i<_nmon; i++ )
delete[] _prodmon[i];
1087 for(
unsigned int i=0; i<_nvar; i++ )
delete[] _bndpow[i];
1096 template <
typename T>
inline void
1097 TModel<T>::_set_prodmon()
1099 _prodmon =
new unsigned int*[_nmon];
1100 _prodmon[0] =
new unsigned int[_nmon+1];
1101 _prodmon[0][0] = _nmon;
1102 for(
unsigned int i=1; i<=_nmon; i++ ) _prodmon[0][i] = i-1;
1103 #ifdef MC__TMODEL_DEBUG
1104 std::ostringstream ohead0;
1105 ohead0 <<
"_prodmon[" << 0 <<
"]";
1106 mc::display( 1, _nmon+1, _prodmon[0], 1, ohead0.str(), std::cout );
1108 if( !_nord )
return;
1110 unsigned int *iexp =
new unsigned int[_nvar];
1111 for(
unsigned int i=1; i<_nord; i++ ){
1112 for(
unsigned int j=_posord[i]; j<_posord[i+1]; j++ ){
1113 _prodmon[j] =
new unsigned int [_posord[_nord+1-i]+1];
1114 _prodmon[j][0] = _posord[_nord+1-i];
1115 for(
unsigned int k=0; k<_posord[_nord+1-i]; k++ ){
1116 for(
unsigned int in=0; in<_nvar; in++ )
1117 iexp[in] = _expmon[j*_nvar+in] + _expmon[k*_nvar+in] ;
1118 _prodmon[j][k+1] = _loc_expmon( iexp );
1120 #ifdef MC__TMODEL_DEBUG
1121 std::ostringstream oheadj;
1122 oheadj <<
"_prodmon[" << j <<
"]";
1123 mc::display( 1, _posord[_nord+1-i]+1, _prodmon[j], 1, oheadj.str(), std::cout );
1129 for(
unsigned int i=_posord[_nord]; i<_nmon; i++ ){
1130 _prodmon[i] =
new unsigned int[2];
1133 #ifdef MC__TMODEL_DEBUG
1134 std::ostringstream oheadi;
1135 oheadi <<
"_prodmon[" << i <<
"]";
1136 mc::display( 1, 2, _prodmon[i], 1, oheadi.str(), std::cout );
1141 template <
typename T>
inline void
1142 TModel<T>::_set_bndpow
1143 (
const unsigned int ivar,
const T&X,
const double Xref,
1144 const double scaling )
1146 if( ivar>=_nvar )
throw Exceptions( Exceptions::INIT );
1148 delete[] _bndpow[ivar];
1149 _bndpow[ivar] =
new T [_nord+1];
1150 _refpoint[ivar] = Xref/scaling;
1151 _scaling[ivar] = scaling;
1152 T Xr = X/scaling - _refpoint[ivar];
1153 _bndpow[ivar][0] = 1.;
1154 for(
unsigned int i=1; i<=_nord; i++ ){
1155 _bndpow[ivar][i] = Op<T>::pow(Xr,(
int)i);
1160 template <
typename T>
inline void
1161 TModel<T>::_set_bndmon()
1163 if( !_modvar )
return;
1166 for(
unsigned int i=1; i<_nmon; i++ ){
1168 for(
unsigned int j=0; j<_nvar; j++)
1169 if( _bndpow[j] ) _bndmon[i] *= _bndpow[j][_expmon[i*_nvar+j]];
1173 #ifdef MC__TMODEL_DEBUG
1174 mc::display( 1, _nmon, _bndmon, 1,
"_bndmon", std::cout );
1178 template <
typename T>
inline void
1180 (
const unsigned int ivar,
double*coef )
const
1182 if( ivar>=_nvar || !_nord )
return;
1183 double dscale = Op<T>::diam(_bndpow[ivar][1]);
1184 double fscale = Op<T>::l(_bndpow[ivar][1]) / dscale;
1185 unsigned int iexp[_nvar];
1186 for(
unsigned int imon=0; imon<_nmon; imon++ ){
1187 unsigned int iord = 0;
1188 for(
unsigned int jvar=0; jvar<_nvar; jvar++ ){
1189 iexp[jvar] = _expmon[imon*_nvar+jvar];
1192 double coefmod = coef[imon] * std::pow(dscale,iexp[ivar]);
1194 for(
unsigned int k=1; k<=_nord-iord; k++, iexp[ivar]++ ){
1195 coefmod += coef[_loc_expmon(iexp)]
1196 * _get_binom(iexp[ivar],k) * std::pow(fscale,k)
1197 * std::pow(dscale,iexp[ivar]);
1199 #ifdef MC__TMODEL_DEBUG_SCALE
1200 std::cout <<
" a" << std::left << std::setw(3) << imon <<
" = "
1201 << std::right << std::setw(12) << coef[imon]
1202 << std::right << std::setw(14) << coefmod << std::endl;
1204 coef[imon] = coefmod;
1209 template <
typename T>
inline TVar<T>
1210 TModel<T>::_univ_bernstein
1211 (
const TVar<T>&TV, puniv f, puniv df, punivext If,
1212 punivext d2If,
const double*rusr,
const int*iusr )
1214 assert( TV._TM ==
this );
1216 const T B( TV.B() );
1217 const TVar<T> TVs( (TV-Op<T>::l(B))/Op<T>::diam(B) );
1218 TVar<T> TV2(
this, 0. ), MON( 1. );
1219 for(
unsigned int j=0; j<=_nord; j++ ){
1220 double sign = ( j%2? -1.: 1. );
1221 _cbern[j] = sign * _get_binom(_nord,j) * f( Op<T>::l(B), rusr, iusr );
1222 for(
unsigned int i=1; i<=j; i++ ){
1224 _cbern[j] += sign * _get_binom(_nord,i) * _get_binom(_nord-i,j-i)
1225 * f( Op<T>::l(B)+(
double)i/(
double)_nord*Op<T>::diam(B), rusr, iusr );
1227 TV2 += _cbern[j] * MON;
1231 T R( - d2If( B, rusr, iusr ) * sqr(Op<T>::diam(B)) / 2. / _nord );
1232 if( Op<T>::l( R ) * Op<T>::u( R ) < 0 )
1233 return TV2 + Op<T>::zeroone() * R;
1235 if( !options.BERNSTEIN_OPT ){
1236 T R0( If(B,rusr,iusr) );
1237 R = Op<T>::zeroone() * ( Op<T>::u(R) <= 0.?
1238 -std::min( Op<T>::diam(R0), -Op<T>::l(R) ):
1239 std::min( Op<T>::diam(R0), Op<T>::u(R) ) );
1243 std::pair<unsigned int,const double*> bern = std::make_pair(_nord+1,_cbern);
1244 double xopt = _goldsect( 0., 1., _dgap_bernstein, rusr, iusr, df, B, bern );
1245 R = Op<T>::zeroone() * _gap_bernstein( xopt, rusr, iusr, f, B, bern );
1249 template <
typename T>
inline double
1250 TModel<T>::_dgap_bernstein
1251 (
const double x,
const double*rusr,
const int*iusr, puniv df,
1252 const T&B,
const std::pair<unsigned int,const double*>&bern )
1254 double phi = df(Op<T>::l(B)+x*Op<T>::diam(B),rusr,iusr)*Op<T>::diam(B);
1255 for(
unsigned int j=1; j<bern.first; j++ )
1256 phi -= j * bern.second[j] * std::pow(x,j-1);
1260 template <
typename T>
inline double
1261 TModel<T>::_gap_bernstein
1262 (
const double x,
const double*rusr,
const int*iusr, puniv f,
1263 const T&B,
const std::pair<unsigned int,const double*>&bern )
1265 double phi = f(Op<T>::l(B)+x*Op<T>::diam(B),rusr,iusr);
1266 for(
unsigned int j=0; j<bern.first; j++ )
1267 phi -= bern.second[j] * std::pow(x,j);
1271 template <
typename T>
inline double
1272 TModel<T>::_goldsect
1273 (
const double xL,
const double xU, punivopt fopt,
const double*rusr,
1274 const int*iusr, puniv df,
const T&Ix,
1275 const std::pair<unsigned int,const double*>&bern )
1277 const double phi = 2.-(1.+std::sqrt(5.))/2.;
1278 const double fL = fopt( xL, rusr, iusr, df, Ix, bern ),
1279 fU = fopt( xU, rusr, iusr, df, Ix, bern );
1280 if( fL*fU > 0 )
throw Exceptions( Exceptions::BERNSTEIN );
1281 const double xm = xU-phi*(xU-xL),
1282 fm = fopt( xm, rusr, iusr, df, Ix, bern );
1283 return _goldsect_iter(
true, xL, fL, xm, fm, xU, fU, fopt, rusr, iusr,
1287 template <
typename T>
inline double
1288 TModel<T>::_goldsect_iter
1289 (
const bool init,
const double a,
const double fa,
const double b,
1290 const double fb,
const double c,
const double fc, punivopt fopt,
1291 const double*rusr,
const int*iusr, puniv df,
const T&Ix,
1292 const std::pair<unsigned int,const double*>&bern )
1296 static unsigned int iter;
1297 iter = ( init? 1: iter+1 );
1298 const double phi = 2.-(1.+std::sqrt(5.))/2.;
1299 bool b_then_x = ( c-b > b-a );
1300 double x = ( b_then_x? b+phi*(c-b): b-phi*(b-a) );
1301 if( std::fabs(c-a) < options.BERNSTEIN_TOL*(std::fabs(b)+std::fabs(x))
1302 || iter > options.BERNSTEIN_MAXIT )
return (c+a)/2.;
1303 double fx = fopt( x, rusr, iusr, df, Ix, bern );
1306 _goldsect_iter(
false, a, fa, b, fb, x, fx, fopt, rusr, iusr, df, Ix, bern ):
1307 _goldsect_iter(
false, b, fb, x, fx, c, fc, fopt, rusr, iusr, df, Ix, bern ) );
1309 _goldsect_iter(
false, a, fa, x, fx, b, fb, fopt, rusr, iusr, df, Ix, bern ):
1310 _goldsect_iter(
false, x, fx, b, fb, c, fc, fopt, rusr, iusr, df, Ix, bern ) );
1315 template <
typename T>
inline TVar<T>&
1317 (
const TVar<T>&TV )
1321 #ifdef MC__TMODEL_CHECK_PMODEL
1322 if( _TM !=
dynamic_cast< TModel<T>*
>( PolyVar<T>::_PM ) ) assert(
false );
1327 template <
typename T>
inline
1337 template <
typename T>
inline TVar<T>&
1341 if( _TM ){ _TM = 0; _resize( _TM ); }
1348 template <
typename T>
inline
1358 template <
typename T>
inline TVar<T>&
1362 if( _TM ){ _TM = 0; _resize( _TM ); }
1369 template <
typename T>
inline
1371 ( TModel<T>*TM,
const double d )
1372 : PolyVar<T>( TM ), _TM( TM )
1381 template <
typename T>
inline
1383 ( TModel<T>*TM,
const T&B )
1384 : PolyVar<T>( TM ), _TM( TM )
1387 for(
unsigned int i=0; i<
nmon(); i++ )
_coefmon[i] = 0.;
1388 for(
unsigned int i=0; i<=
nord(); i++)
_bndord[i] = 0.;
1396 template <
typename T>
template <
typename U>
inline
1402 _coefmon[0] = TVtrunc._coefmon[0];
1403 TVtrunc._coefmon[0] = 0. ;
1404 for(
unsigned int i=1; TM && i<nmon(); i++ ){
1405 if( TVtrunc.TM && i < TVtrunc.nmon() ){
1406 _coefmon[i] = TVtrunc._coefmon[i];
1407 TVtrunc._coefmon[i] = 0.;
1412 TVtrunc._bndord_uptd =
false;
1413 *_bndrem = T( TVtrunc.B() );
1414 _bndord_uptd =
false;
1418 template <
typename T>
template <
typename U>
inline
1425 _coefmon[0] = TVtrunc._coefmon[0];
1426 TVtrunc._coefmon[0] = 0. ;
1427 for(
unsigned int i=1; TM && i<nmon(); i++ ){
1428 if( TVtrunc.TM && i < TVtrunc.nmon() ){
1429 _coefmon[i] = TVtrunc._coefmon[i];
1430 TVtrunc._coefmon[i] = 0.;
1435 TVtrunc._bndord_uptd =
false;
1436 *_bndrem = (*method)( TVtrunc.B() );
1437 _bndord_uptd =
false;
1441 template <
typename T>
template <
typename U>
inline
1448 _coefmon[0] = TVtrunc._coefmon[0];
1449 TVtrunc._coefmon[0] = 0. ;
1450 for(
unsigned int i=1; TM && i<nmon(); i++ ){
1451 if( TVtrunc.TM && i < TVtrunc.nmon() ){
1452 _coefmon[i] = TVtrunc._coefmon[i];
1453 TVtrunc._coefmon[i] = 0.;
1458 TVtrunc._bndord_uptd =
false;
1459 *_bndrem = (TVtrunc.B().*method)();
1460 _bndord_uptd =
false;
1464 template <
typename T>
inline TVar<T>&
1466 (
const unsigned int ivar,
const T&X,
const double Xref )
1471 double scaling = ( _TM->options.SCALE_VARIABLES?
Op<T>::diam(X)/2.: 1. );
1472 if( isequal( scaling, 0. ) ) scaling = 1.;
1473 _TM->_set_bndpow( ivar, X, Xref, scaling );
1478 for(
unsigned int i=1; i<nmon(); i++ ) _coefmon[i] = 0.;
1479 if( nord() > 0 ) _coefmon[nvar()-ivar] = scaling;
1482 _bndord[0] = _coefmon[0];
1484 _bndord[1] = X-_coefmon[0];
1485 for(
unsigned int i=2; i<nord()+2; i++) _bndord[i] = 0.;
1489 *_bndrem = X-_coefmon[0];
1490 _set_bndpol( _coefmon[0] );
1492 _bndord_uptd =
true;
1497 template <
typename T>
inline
1505 template <
typename T>
inline
1507 (
TModel<T>*TM,
const unsigned int ivar,
const T&X,
const double Xref )
1510 _set( ivar, X, Xref );
1513 template <
typename T>
inline void
1516 if( !_TM || _bndord_uptd )
return;
1518 _bndord[0] = _coefmon[0];
1519 for(
unsigned int i=1; i<=nord(); i++ ){
1521 for(
unsigned int j=_posord(i); j<_posord(i+1); j++ )
1522 _bndord[i] += _coefmon[j] * _bndmon(j);
1524 _bndord_uptd =
true;
1527 template <
typename T>
inline T
1528 TVar<T>::_polybound_eigen()
const
1530 static const double TOL = 1e-8;
1532 T bndpol = _coefmon[0];
1533 if( nord() == 1 ) bndpol += _bndord[1];
1535 else if( nord() > 1 ){
1536 double*U =
new double[nvar()*nvar()];
1537 for(
unsigned int i=0; i<nvar(); i++ ){
1538 for(
unsigned int j=0; j<i; j++ ){
1539 U[nvar()*(nvar()-i-1)+nvar()-j-1] = 0.;
1540 U[nvar()*(nvar()-j-1)+nvar()-i-1] = _coefmon[_prodmon(i+1,j+2)]/2.;
1542 U[(nvar()+1)*(nvar()-i-1)] = _coefmon[_prodmon(i+1,i+2)];
1544 #ifdef MC__TVAR_DEBUG_EIGEN
1545 display( nvar(), nvar(), U, nvar(),
"Matrix U", std::cout );
1547 double*D = mc::dsyev_wrapper( nvar(), U,
true );
1553 #ifdef MC__TVAR_HYBRID_EIGEN
1557 for(
unsigned int i=0; i<nvar(); i++ ){
1560 for(
unsigned int k=0; k<nvar(); k++ ){
1561 linaux += U[i*nvar()+k] * _coefmon[nvar()-k];
1562 bndaux += U[i*nvar()+k] * _bndmon(nvar()-k);
1564 #ifdef MC__TVAR_DEBUG_EIGEN
1565 std::cout << i <<
": LINAUX = " << linaux
1566 <<
" BNDAUX = " << bndaux << std::endl;
1568 #ifdef MC__TVAR_HYBRID_EIGEN
1569 bndtype1 += _coefmon[i+1] * _bndmon(i+1) + D[i] * Op<T>::sqr( bndaux );
1571 #ifdef MC__TVAR_DEBUG_EIGEN
1572 std::cout << std::endl <<
"BNDTYPE1: " << bndtype1 << std::endl;
1574 if( std::fabs(D[i]) > TOL )
1575 bndtype2 += D[i] * Op<T>::sqr( linaux/D[i]/2. + bndaux )
1576 - linaux*linaux/D[i]/4.;
1579 bndtype2 += linaux * bndaux + D[i] * Op<T>::sqr( bndaux );
1580 #ifdef MC__TVAR_DEBUG_EIGEN
1581 std::cout <<
"BNDTYPE2: " << bndtype2 << std::endl;
1587 #ifdef MC__TVAR_HYBRID_EIGEN
1588 if( !Op<T>::inter( bndtype1, bndtype1, bndtype2 ) ){
1590 #ifdef MC__TVAR_DEBUG_EIGEN
1591 std::cout <<
"BNDTYPE3: " << bndtype1 << std::endl;
1599 #ifdef MC__TVAR_DEBUG_EIGEN
1600 int tmp; std::cin >> tmp;
1603 for(
unsigned int i=3; i<=nord(); i++ ) bndpol += _bndord[i];
1607 template <
typename T>
inline T
1608 TVar<T>::_polybound_LSB()
const
1610 static const double TOL = 1e-8;
1612 T bndpol = _coefmon[0];
1613 if( nord() == 1 ) bndpol += _bndord[1];
1614 else if( nord() > 1 ){
1615 for(
unsigned int i=1; i<=nvar(); i++ ){
1617 unsigned int ii = _prodmon(i,i+1);
1618 if( std::fabs(_coefmon[ii]) > TOL )
1619 bndpol += _coefmon[ii] * Op<T>::sqr( _coefmon[i]/_coefmon[ii]/2.
1620 + _bndmon(i) ) - _coefmon[i]*_coefmon[i]/_coefmon[ii]/4.;
1622 bndpol += _coefmon[i] * _bndmon(i) + _coefmon[ii] * _bndmon(ii);
1624 for(
unsigned int k=i+1; k<=nvar(); k++ ){
1625 unsigned int ik = _prodmon(i,k+1) ;
1626 bndpol += _coefmon[ik] * _bndmon(ik);
1631 for(
unsigned int i=3; i<=nord(); i++ ) bndpol += _bndord[i];
1635 template <
typename T>
inline T
1636 TVar<T>::_polybound_bernstein()
const
1639 double *coeftrans =
new double[nmon()];
1640 for(
unsigned int imon=0; imon<nmon(); imon++ )
1641 coeftrans[imon] = _coefmon[imon];
1642 for(
unsigned int ivar=0; ivar<nvar(); ivar++ )
1643 _TM->_scale( ivar, coeftrans );
1646 const unsigned int maxord = (_TM->options.BOUNDER_ORDER>nord()?
1647 _TM->options.BOUNDER_ORDER: nord() );
1648 const unsigned int maxmon = std::pow(maxord+1,nvar());
1649 _TM->_ext_expmon( maxord,
true );
1652 T bndpol = coeftrans[0];
1653 #ifdef MC__TVAR_DEBUG_BERSTEIN
1654 std::cout <<
"\n0: " << bndpol << std::endl;
1656 for(
unsigned int jmon=1; jmon<maxmon; jmon++ ){
1657 const unsigned int*jexp = _expmon(jmon);
1658 const double coefbern = _coef_bernstein( coeftrans, jexp, jmon, maxord );
1659 bndpol = Op<T>::hull( bndpol, coefbern );
1660 #ifdef MC__TVAR_DEBUG_BERNSTEIN
1661 std::cout << jmon <<
" [";
1662 for(
unsigned int ivar=0; ivar<nvar(); ivar++ )
1663 std::cout << std::setw(3) << jexp[ivar];
1664 std::cout <<
"] : " << coefbern <<
" : " << bndpol << std::endl;
1671 template <
typename T>
inline double
1672 TVar<T>::_coef_bernstein
1673 (
const double*coefmon,
const unsigned int*jexp,
1674 const unsigned int jmon,
const unsigned int maxord )
const
1677 double coefbern = coefmon[0];
1678 for(
unsigned int imon=1; imon<=std::min(jmon,nmon()-1); imon++ ){
1679 const unsigned int*iexp = _expmon(imon);
1681 bool inrange =
true;
1682 for(
unsigned int ivar=0; ivar<nvar() && inrange; ivar++ )
1683 if( iexp[ivar] > jexp[ivar] ) inrange =
false;
1684 if( !inrange )
continue;
1685 double termbern = coefmon[imon];
1686 #ifdef MC__TVAR_DEBUG_BERSTEIN
1688 for(
unsigned int ivar=0; ivar<nvar(); ivar++ )
1689 std::cout << std::setw(3) << iexp[ivar];
1692 for(
unsigned int ivar=0; ivar<nvar(); ivar++ )
1693 termbern *= (
double)_get_binom(jexp[ivar],iexp[ivar])
1694 / (double)_get_binom(maxord,iexp[ivar]);
1695 coefbern += termbern;
1697 #ifdef MC__TVAR_DEBUG_BERSTEIN
1698 std::cout << std::endl;
1703 template <
typename T>
inline T
1704 TVar<T>::_polybound_naive()
const
1706 T bndpol = _coefmon[0];
1707 for(
unsigned int i=1; i<=nord(); i++ ) bndpol += _bndord[i];
1711 template <
typename T>
inline T
1713 (
const int type )
const
1715 if( !_TM )
return _coefmon[0];
1720 return _polybound_naive();
1722 return _polybound_bernstein();
1724 return _polybound_eigen();
1726 T bndpol, bndpol_LSB = _polybound_LSB();
1727 if( !Op<T>::inter( bndpol, bndpol_LSB, _polybound_eigen() ) )
1733 return _polybound_LSB();
1737 template <
typename T>
inline double
1739 (
const double*x )
const
1741 if( !_TM )
return _coefmon[0];
1742 double Pval = _coefmon[0];
1743 for(
unsigned int i=1; i<nmon(); i++ ){
1745 for(
unsigned int k=0; k<nvar(); k++ )
1746 valmon *= std::pow( x[k]/_scaling(k)-_refpoint(k), _expmon(i)[k] );
1747 Pval += _coefmon[i] * valmon;
1752 template <
typename T>
inline double
1754 (
const bool reset )
1756 const double coefconst = _coefmon[0];
1757 if( reset ) *
this -= coefconst;
1761 template <
typename T>
inline double*
1764 if( !nvar() || !nord() )
return 0;
1766 double*plin =
new double[nvar()];
1767 for(
unsigned int i=0; i<nvar(); i++ )
1768 plin[i] = _coefmon[nvar()-i] / _scaling(i);
1772 template <
typename T>
inline double
1774 (
const unsigned int ivar,
const bool reset )
1776 if( ivar>=nvar() || !nord() )
return 0.;
1777 const double coeflin = _coefmon[nvar()-ivar] / _scaling(ivar);
1778 if( reset ){ _coefmon[nvar()-ivar] = 0.; _bndord_uptd =
false; _unset_bndpol(); }
1782 template <
typename T>
inline std::ostream&
1784 ( std::ostream&out,
const TVar<T>&TV )
1787 << std::scientific << std::setprecision(5)
1790 #ifdef MC__TVAR_DISPLAY_EXT
1791 out <<
" TM: " << TV._TM <<
" PM: " << TV._PM << std::endl;
1796 out <<
" a0 = " << std::right << std::setw(12) << TV._coefmon[0]
1798 <<
" R = " << *(TV._bndrem) << std::endl;
1803 out << std::setprecision(TV._TM->options.DISPLAY_DIGITS);
1804 for(
unsigned int i=0; i<TV.nmon(); i++ ){
1806 for(
unsigned int k=0; k<TV.nvar(); k++ )
1807 scal *= std::pow( TV._scaling(k), TV._expmon(i)[k] );
1808 out <<
" a" << std::left << std::setw(4) << i <<
" = "
1809 << std::right << std::setw(TV._TM->options.DISPLAY_DIGITS+7)
1810 << TV._coefmon[i]*scal <<
" ";
1811 for(
unsigned int k=0; k<TV.nvar(); k++ )
1812 out << std::setw(3) << TV._expmon(i)[k];
1815 #ifdef MC__TVAR_DISPLAY_EXT
1817 for(
unsigned int j=0; j<=TV.nord(); j++ )
1818 out << std::right <<
" B" << std::left << std::setw(2) << j <<
" = "
1819 << std::right << TV.bndord(j) << std::endl;
1822 out << std::right <<
" R = " << *(TV._bndrem)
1827 out << std::right <<
" B = " << TV.B()
1833 template <
typename T>
inline TVar<T>
1835 (
const TVar<T>&TV )
1840 template <
typename T>
template <
typename U>
inline TVar<T>&
1841 TVar<T>::operator +=
1842 (
const TVar<U>&TV )
1845 _coefmon[0] += TV._coefmon[0];
1846 if( _TM && _bndord_uptd ) _bndord[0] += TV._coefmon[0];
1847 *_bndrem += *(TV._bndrem);
1848 if( _bndpol ) *_bndpol += TV._coefmon[0];
1852 *
this = TV; *
this += TV2;
1857 for(
unsigned int i=0; i<nmon(); i++ )
1858 _coefmon[i] += TV._coefmon[i];
1859 if( TV._bndord_uptd )
1860 for(
unsigned int i=0; _bndord_uptd && i<=nord(); i++ )
1861 _bndord[i] += TV._bndord[i];
1862 else _bndord_uptd =
false;
1863 *_bndrem += *(TV._bndrem);
1864 if( _bndpol && TV._bndpol ) *_bndpol += *(TV._bndpol);
1865 else _unset_bndpol();
1867 if( _TM && _TM->options.CENTER_REMAINDER ) _center();
1871 template <
typename T,
typename U>
inline TVar<T>
1873 (
const TVar<T>&TV1,
const TVar<U>&TV2 )
1880 template <
typename T>
inline TVar<T>&
1881 TVar<T>::operator +=
1885 if( _TM && _bndord_uptd ) _bndord[0] += c;
1886 if( _bndpol ) *_bndpol += c;
1887 if( _TM && _TM->options.CENTER_REMAINDER ) _center();
1891 template <
typename T>
inline TVar<T>
1893 (
const TVar<T>&TV1,
const double c )
1900 template <
typename T>
inline TVar<T>
1902 (
const double c,
const TVar<T>&TV2 )
1909 template <
typename T>
template <
typename U>
inline TVar<T>&
1910 TVar<T>::operator +=
1914 if( _TM && _TM->options.CENTER_REMAINDER ) _center();
1918 template <
typename T>
inline TVar<T>
1920 (
const TVar<T>&TV1,
const T&I )
1927 template <
typename T>
inline TVar<T>
1929 (
const T&I,
const TVar<T>&TV2 )
1936 template <
typename T>
inline TVar<T>
1938 (
const TVar<T>&TV )
1942 TV2._coefmon[0] = -TV._coefmon[0];
1943 *TV2._bndrem = - *TV._bndrem;
1944 if( TV._bndpol ) TV2._set_bndpol( - *TV._bndpol );
1947 TVar<T>& TV2 = *TV._TV();
1948 TV2._unset_bndpol(); TV2._bndord_uptd =
false;
1949 for(
unsigned int i=0; i<TV.nmon(); i++ ) TV2._coefmon[i] = -TV._coefmon[i];
1950 *TV2._bndrem = - *TV._bndrem;
1951 if( TV._bndord_uptd ){
1952 for(
unsigned int i=0; i<=TV.nord(); i++ ) TV2._bndord[i] = -TV._bndord[i];
1953 TV2._bndord_uptd =
true;
1955 if( TV._bndpol ) TV2._set_bndpol( - *TV._bndpol );
1956 if( TV._TM->options.CENTER_REMAINDER ) TV2._center();
1960 template <
typename T>
template <
typename U>
inline TVar<T>&
1961 TVar<T>::operator -=
1962 (
const TVar<U>&TV )
1965 _coefmon[0] -= TV._coefmon[0];
1966 if( _TM && _bndord_uptd ) _bndord[0] -= TV._coefmon[0];
1967 *_bndrem -= *(TV._bndrem);
1968 if( _bndpol ) *_bndpol -= TV._coefmon[0];
1972 *
this = -TV; *
this += TV2;
1977 for(
unsigned int i=0; i<nmon(); i++ )
1978 _coefmon[i] -= TV._coefmon[i];
1979 if( TV._bndord_uptd )
1980 for(
unsigned int i=0; _bndord_uptd && i<=nord(); i++ )
1981 _bndord[i] -= TV._bndord[i];
1982 else _bndord_uptd =
false;
1983 *_bndrem -= *(TV._bndrem);
1984 if( _bndpol && TV._bndpol ) *_bndpol -= *(TV._bndpol);
1985 else _unset_bndpol();
1987 if( _TM && _TM->options.CENTER_REMAINDER ) _center();
1991 template <
typename T,
typename U>
inline TVar<T>
1993 (
const TVar<T>&TV1,
const TVar<U>&TV2 )
2000 template <
typename T>
inline TVar<T>&
2001 TVar<T>::operator -=
2005 if( _TM && _bndord_uptd ) _bndord[0] -= c;
2006 if( _bndpol ) *_bndpol -= c;
2007 if( _TM && _TM->options.CENTER_REMAINDER ) _center();
2011 template <
typename T>
inline TVar<T>
2013 (
const TVar<T>&TV1,
const double c )
2020 template <
typename T>
inline TVar<T>
2022 (
const double c,
const TVar<T>&TV2 )
2024 TVar<T> TV3( -TV2 );
2029 template <
typename T>
template <
typename U>
inline TVar<T>&
2030 TVar<T>::operator -=
2034 if( _TM && _TM->options.CENTER_REMAINDER ) _center();
2038 template <
typename T>
inline TVar<T>
2040 (
const TVar<T>&TV1,
const T&I )
2047 template <
typename T>
inline TVar<T>
2049 (
const T&I,
const TVar<T>&TV2 )
2051 TVar<T> TV3( -TV2 );
2056 template <
typename T>
inline TVar<T>&
2057 TVar<T>::operator *=
2058 (
const TVar<T>&TV )
2060 TVar<T> TV2( *
this );
2065 template <
typename T>
inline TVar<T>
2067 (
const TVar<T>&TV1,
const TVar<T>&TV2 )
2069 if( !TV2._TM )
return( TV1 * TV2._coefmon[0] + TV1 * *(TV2._bndrem) );
2070 else if( !TV1._TM )
return( TV2 * TV1._coefmon[0] + TV2 * *(TV1._bndrem) );
2071 else if( &TV1 == &TV2 )
return sqr(TV1);
2073 if( TV1._TM != TV2._TM )
2075 TVar<T>& TV3 = *TV1._TV();
2076 for(
unsigned int i=0; i<TV3.nmon(); i++ ) TV3._coefmon[i] = 0.;
2079 for(
unsigned int i=0; i<TV3._posord(TV3.nord()/2+1); i++){
2080 TV3._coefmon[TV3._prodmon(i,i+1)] += TV1._coefmon[i] * TV2._coefmon[i];
2081 for(
unsigned int j=i+1; j<TV3._prodmon(i,0); j++ )
2082 TV3._coefmon[TV3._prodmon(i,j+1)] += TV1._coefmon[i] * TV2._coefmon[j]
2083 + TV1._coefmon[j] * TV2._coefmon[i];
2086 TV1._update_bndord(); TV2._update_bndord();
2088 for(
unsigned int i=0; i<=TV3.nord()+1; i++ ){
2090 for(
unsigned int j=TV3.nord()+1-i; j<=TV3.nord()+1; j++ ){
2092 r2 += TV2._bndord[j];
2094 s1 += TV2._bndord[i] * r1 ;
2095 s2 += TV1._bndord[i] * r2 ;
2097 if( !Op<T>::inter( *(TV3._bndrem), s1, s2) ){
2098 *(TV3._bndrem) = s1;
2102 TV3._unset_bndpol();
2103 TV3._bndord_uptd =
false;
2104 if( TV3._TM->options.CENTER_REMAINDER ) TV3._center();
2108 template <
typename T>
inline TVar<T>
2110 (
const TVar<T>&TV )
2114 TV2._coefmon[0] *= TV2._coefmon[0];
2115 *(TV2._bndrem) *= 2. + *(TV2._bndrem);
2116 TV2._unset_bndpol();
2121 TVar<T> TV2( TV._TM, 0. );
2122 for(
unsigned int i=0; i<TV2._posord(TV2.nord()/2+1); i++){
2123 TV2._coefmon[TV2._prodmon(i,i+1)] += TV._coefmon[i] * TV._coefmon[i];
2124 for(
unsigned int j=i+1; j<TV2._prodmon(i,0); j++ )
2125 TV2._coefmon[TV2._prodmon(i,j+1)] += TV._coefmon[i] * TV._coefmon[j] * 2.;
2128 TV._update_bndord();
2130 for(
unsigned int i=0; i<=TV2.nord()+1; i++ ){
2131 unsigned int k = std::max(TV2.nord()+1-i, i+1);
2133 for(
unsigned int j=k; j<=TV2.nord()+1; j++ )
2135 s += TV._bndord[i] * r;
2139 for(
unsigned int i=TV2.nord()/2+1; i<=TV2.nord()+1; i++ )
2140 r += Op<T>::sqr(TV._bndord[i]) ;
2141 *(TV2._bndrem) = 2. * s + r;
2144 TV2._unset_bndpol();
2145 TV2._bndord_uptd =
false;
2146 if( TV2._TM->options.CENTER_REMAINDER ) TV2._center();
2150 template <
typename T>
inline TVar<T>&
2151 TVar<T>::operator *=
2159 for(
unsigned int i=0; i<nmon(); i++ ) _coefmon[i] *= c;
2160 for(
unsigned int i=0; _bndord_uptd && i<=nord(); i++ ) _bndord[i] *= c;
2163 if( _bndpol ) *_bndpol *= c;
2167 template <
typename T>
inline TVar<T>
2169 (
const TVar<T>&TV1,
const double c )
2176 template <
typename T>
inline TVar<T>
2178 (
const double c,
const TVar<T>&TV2 )
2185 template <
typename T>
inline TVar<T>&
2186 TVar<T>::operator *=
2190 *(_bndrem) += _coefmon[0];
2195 const double Imid = Op<T>::mid(I);
2197 for(
unsigned int i=0; i<nmon(); i++ ) _coefmon[i] *= Imid;
2198 for(
unsigned int i=0; _bndord_uptd && i<=nord(); i++ ) _bndord[i] *= Imid;
2200 *_bndrem += (I-Imid)*Icur;
2204 if( _TM && _TM->options.CENTER_REMAINDER ) _center();
2208 template <
typename T>
inline TVar<T>
2210 (
const TVar<T>&TV1,
const T&I )
2217 template <
typename T>
inline TVar<T>
2219 (
const T&I,
const TVar<T>&TV2 )
2226 template <
typename T>
inline TVar<T>&
2227 TVar<T>::operator /=
2228 (
const TVar<T>&TV )
2234 template <
typename T>
inline TVar<T>
2236 (
const TVar<T>&TV1,
const TVar<T>&TV2 )
2238 return TV1 * inv(TV2);
2241 template <
typename T>
inline TVar<T>&
2242 TVar<T>::operator /=
2245 if ( isequal( c, 0. ))
2251 template <
typename T>
inline TVar<T>
2253 (
const TVar<T>&TV,
const double c )
2255 if ( isequal( c, 0. ))
2260 template <
typename T>
inline TVar<T>
2262 (
const double c,
const TVar<T>&TV )
2267 template <
typename T>
inline TVar<T>
2269 (
const TVar<T>&TV )
2272 return TVar<T>( Op<T>::inv(TV._coefmon[0] + *(TV._bndrem)) );
2274 TVar<T> TV2 = TV._TM->_inv_taylor( TV );
2275 if( TV2._TM->options.BERNSTEIN_USE ){
2276 TVar<T> TV2b = TV._TM->_inv_bernstein( TV );
2277 if( Op<T>::diam(*(TV2b._bndrem)) < Op<T>::diam(*(TV2._bndrem)) )
2281 TV2._bndord_uptd =
false;
2282 TV2._unset_bndpol();
2283 if( TV2._TM->options.CENTER_REMAINDER ) TV2._center();
2287 template <
typename T>
inline TVar<T>
2288 TModel<T>::_inv_taylor
2289 (
const TVar<T>&TV )
2291 const T I( TV.B() );
2292 if ( Op<T>::l(TV.B()) <= 0. && Op<T>::u(TV.B()) >= 0. )
2294 double x0 = ( TV._TM->options.REF_MIDPOINT? Op<T>::mid(I):
2296 const TVar<T> TVmx0( TV - x0 );
2297 const T Imx0( I - x0 );
2299 TVar<T> TV2( TV._TM, 1. ), MON( 1. );
2300 #ifdef MC__TMODEL_TIGHT_REMAINDER
2301 double R[3] = { 0., 1./Op<T>::l(I)-1./x0, 1./Op<T>::u(I)-1./x0 };
2302 double aL = 1./x0, aU = 1./x0;
2304 for(
unsigned int i=1; i<=TV.nord(); i++ ){
2305 MON *= TVmx0 / (-x0);
2307 #ifdef MC__TMODEL_TIGHT_REMAINDER
2308 aL *= (Op<T>::l(I)-x0) / (-x0);
2310 aU *= (Op<T>::u(I)-x0) / (-x0);
2315 #ifdef MC__TMODEL_TIGHT_REMAINDER
2316 TV2 += max(3,R) + Op<T>::zeroone() * ( min(3,R) - max(3,R) );
2318 TV2 += Op<T>::pow( -Imx0, (
int)TV2.nord()+1 )
2319 / Op<T>::pow( Op<T>::zeroone()*Imx0+x0, (int)TV2.nord()+2 );
2324 template <
typename T>
inline TVar<T>
2325 TModel<T>::_inv_bernstein
2326 (
const TVar<T>&TV )
2328 const T I( TV.B() );
2329 if ( Op<T>::l(TV.B()) <= 0. && Op<T>::u(TV.B()) >= 0. )
2333 (
const double x,
const double*rusr,
const int*iusr )
2336 (
const double x,
const double*rusr,
const int*iusr )
2337 {
return -1./mc::sqr(x); }
2339 (
const T&x,
const double*rusr,
const int*iusr )
2340 {
return Op<T>::inv(x); }
2342 (
const T&x,
const double*rusr,
const int*iusr )
2343 {
return 2.*Op<T>::pow( Op<T>::inv(x), 3 ); }
2345 return TV._TM->_univ_bernstein( TV, loc::inv, loc::dinv, loc::Iinv,
2346 loc::d2Iinv, 0, 0 );
2349 template <
typename T>
inline TVar<T>
2351 (
const TVar<T>&TV )
2354 return TVar<T>( Op<T>::sqrt(TV._coefmon[0] + *(TV._bndrem)) );
2356 TVar<T> TV2 = TV._TM->_sqrt_taylor( TV );
2357 if( TV2._TM->options.BERNSTEIN_USE ){
2358 TVar<T> TV2b = TV._TM->_sqrt_bernstein( TV );
2359 if( Op<T>::diam(*(TV2b._bndrem)) < Op<T>::diam(*(TV2._bndrem)) )
2363 TV2._bndord_uptd =
false;
2364 TV2._unset_bndpol();
2365 if( TV2._TM->options.CENTER_REMAINDER ) TV2._center();
2369 template <
typename T>
inline TVar<T>
2370 TModel<T>::_sqrt_taylor
2371 (
const TVar<T>&TV )
2373 const T I( TV.B() );
2374 if ( Op<T>::l(I) < 0. )
2376 double x0 = ( TV._TM->options.REF_MIDPOINT? Op<T>::mid(I):
2378 const TVar<T> TVmx0( TV - x0 );
2379 const T Imx0( I - x0 );
2382 TVar<T> TV2( TV._TM, 1. ), MON( 1. );
2383 #ifdef MC__TMODEL_TIGHT_REMAINDER
2384 double R[3] = { 0., std::sqrt(Op<T>::l(I))-std::sqrt(x0), std::sqrt(Op<T>::u(I))-std::sqrt(x0) };
2385 double aL = std::sqrt(x0), aU = std::sqrt(x0);
2387 for(
unsigned int i=1; i<=TV.nord(); i++ ){
2390 #ifdef MC__TMODEL_TIGHT_REMAINDER
2391 aL *= (Op<T>::l(I)-x0) / x0;
2393 aU *= (Op<T>::u(I)-x0) / x0;
2396 s *= -(2.*i-1.)/(2.*i+2.);
2398 TV2 *= std::sqrt(x0);
2399 #ifdef MC__TMODEL_TIGHT_REMAINDER
2400 TV2 += max(3,R) + Op<T>::zeroone() * ( min(3,R) - max(3,R) );
2402 TV2 += s * Op<T>::pow( Imx0, (
int)TV2.nord()+1 )
2403 / Op<T>::pow( Op<T>::zeroone()*Imx0+x0, (int)TV2.nord()+1/2 );
2408 template <
typename T>
inline TVar<T>
2409 TModel<T>::_sqrt_bernstein
2410 (
const TVar<T>&TV )
2412 const T I( TV.B() );
2413 if ( Op<T>::l(I) < 0. )
2417 (
const double x,
const double*rusr,
const int*iusr )
2418 {
return std::sqrt(x); }
2420 (
const double x,
const double*rusr,
const int*iusr )
2421 {
return 1./(2.*std::sqrt(x)); }
2423 (
const T&x,
const double*rusr,
const int*iusr )
2424 {
return Op<T>::sqrt(x); }
2426 (
const T&x,
const double*rusr,
const int*iusr )
2427 {
return -Op<T>::inv( 4.*Op<T>::sqrt(Op<T>::pow(x,3)) ); }
2429 return TV._TM->_univ_bernstein( TV, loc::sqrt, loc::dsqrt, loc::Isqrt,
2430 loc::d2Isqrt, 0, 0 );
2433 template <
typename T>
inline TVar<T>
2435 (
const TVar<T>&TV )
2438 return TVar<T>( Op<T>::exp(TV._coefmon[0] + *(TV._bndrem)) );
2440 TVar<T> TV2 = TV._TM->_exp_taylor( TV );
2441 if( TV2._TM->options.BERNSTEIN_USE ){
2442 TVar<T> TV2b = TV._TM->_exp_bernstein( TV );
2444 if( Op<T>::diam(*(TV2b._bndrem)) < Op<T>::diam(*(TV2._bndrem)) )
2448 TV2._bndord_uptd =
false;
2449 TV2._unset_bndpol();
2450 if( TV2._TM->options.CENTER_REMAINDER ) TV2._center();
2454 template <
typename T>
inline TVar<T>
2455 TModel<T>::_exp_taylor
2456 (
const TVar<T>&TV )
2458 const T I( TV.B() );
2459 double x0 = ( TV._TM->options.REF_MIDPOINT? Op<T>::mid(I):
2461 const TVar<T> TVmx0( TV - x0 );
2462 const T Imx0( I - x0 );
2465 TVar<T> TV2( TV._TM, 1. ), MON( 1. );
2466 #ifdef MC__TMODEL_TIGHT_REMAINDER
2467 double R[3] = { 0., std::exp(Op<T>::l(I))-std::exp(x0), std::exp(Op<T>::u(I))-std::exp(x0) };
2468 double aL = std::exp(x0), aU = std::exp(x0);
2470 for(
unsigned int i=1; i<=TV.nord(); i++ ){
2473 #ifdef MC__TMODEL_TIGHT_REMAINDER
2474 aL *= (Op<T>::l(I)-x0);
2476 aU *= (Op<T>::u(I)-x0);
2481 #ifdef MC__TMODEL_TIGHT_REMAINDER
2482 TV2 *= std::exp(x0);
2483 TV2 += max(3,R) + Op<T>::zeroone() * ( min(3,R) - max(3,R) );
2485 TV2 += s * Op<T>::pow( Imx0, (
int)TV2.nord()+1 )
2486 * Op<T>::exp( Op<T>::zeroone()*Imx0 );
2487 TV2 *= std::exp(x0);
2492 template <
typename T>
inline TVar<T>
2493 TModel<T>::_exp_bernstein
2494 (
const TVar<T>&TV )
2498 (
const double x,
const double*rusr,
const int*iusr )
2499 {
return std::exp(x); }
2501 (
const double x,
const double*rusr,
const int*iusr )
2502 {
return std::exp(x); }
2504 (
const T&x,
const double*rusr,
const int*iusr )
2505 {
return Op<T>::exp(x); }
2507 (
const T&x,
const double*rusr,
const int*iusr )
2508 {
return Op<T>::exp(x); }
2510 return TV._TM->_univ_bernstein( TV, loc::exp, loc::dexp, loc::Iexp,
2511 loc::d2Iexp, 0, 0 );
2514 template <
typename T>
inline TVar<T>
2516 (
const TVar<T>&TV )
2519 return TVar<T>( Op<T>::log(TV._coefmon[0] + *(TV._bndrem)) );
2521 TVar<T> TV2 = TV._TM->_log_taylor( TV );
2522 if( TV2._TM->options.BERNSTEIN_USE ){
2523 TVar<T> TV2b = TV._TM->_log_bernstein( TV );
2524 if( Op<T>::diam(*(TV2b._bndrem)) < Op<T>::diam(*(TV2._bndrem)) )
2528 TV2._bndord_uptd =
false;
2529 TV2._unset_bndpol();
2530 if( TV2._TM->options.CENTER_REMAINDER ) TV2._center();
2534 template <
typename T>
inline TVar<T>
2535 TModel<T>::_log_taylor
2536 (
const TVar<T>&TV )
2538 const T I( TV.B() );
2539 if ( Op<T>::l(I) <= 0. )
2541 double x0 = ( TV._TM->options.REF_MIDPOINT? Op<T>::mid(I):
2543 const TVar<T> TVmx0( TV - x0 );
2544 const T Imx0( I - x0 );
2546 TVar<T> TV2( TV._TM, 0. ), MON( -1. );
2547 #ifdef MC__TMODEL_TIGHT_REMAINDER
2548 double R[3] = { 0., std::log(Op<T>::l(I))-std::log(x0), std::log(Op<T>::u(I))-std::log(x0) };
2549 double aL = -1., aU = -1.;
2551 for(
unsigned int i=1; i<=TV.nord(); i++ ){
2552 MON *= TVmx0 / (-x0);
2553 TV2 += MON / (double)i;
2554 #ifdef MC__TMODEL_TIGHT_REMAINDER
2555 aL *= (Op<T>::l(I)-x0) / (-x0);
2556 R[1] -= aL / (double)i;
2557 aU *= (Op<T>::u(I)-x0) / (-x0);
2558 R[2] -= aU / (double)i;
2561 TV2._coefmon[0] += std::log(x0);
2562 #ifdef MC__TMODEL_TIGHT_REMAINDER
2563 TV2 += max(3,R) + Op<T>::zeroone() * ( min(3,R) - max(3,R) );
2565 TV2 -= Op<T>::pow( - Imx0 / ( Op<T>::zeroone()*Imx0+x0 ),
2566 (
int)TV2.nord()+1 ) / ( TV2.nord()+1. );
2571 template <
typename T>
inline TVar<T>
2572 TModel<T>::_log_bernstein
2573 (
const TVar<T>&TV )
2575 const T I( TV.B() );
2576 if ( Op<T>::l(I) <= 0. )
2580 (
const double x,
const double*rusr,
const int*iusr )
2581 {
return std::log(x); }
2583 (
const double x,
const double*rusr,
const int*iusr )
2586 (
const T&x,
const double*rusr,
const int*iusr )
2587 {
return Op<T>::log(x); }
2589 (
const T&x,
const double*rusr,
const int*iusr )
2590 {
return -Op<T>::sqr( Op<T>::inv(x) ); }
2592 return TV._TM->_univ_bernstein( TV, loc::log, loc::dlog, loc::Ilog,
2593 loc::d2Ilog, 0, 0 );
2596 template <
typename T>
inline TVar<T>
2598 (
const TVar<T>&TV )
2600 return TV * log( TV );
2603 template <
typename T>
inline TVar<T>
2605 (
const TVar<T>&TV,
const int n )
2608 return TVar<T>( Op<T>::pow(TV._coefmon[0] + *(TV._bndrem), n) );
2610 if( n < 0 )
return pow( inv( TV ), -n );
2611 TVar<T> TV2( TV._TM->_intpow( TV, n ) );
2612 if( TV2._TM && TV2._TM->options.CENTER_REMAINDER ) TV2._center();
2616 template <
typename T>
inline TVar<T>
2618 (
const TVar<T>&TV,
const int n )
2620 if( n == 0 )
return 1.;
2621 else if( n == 1 )
return TV;
2622 return n%2 ? sqr( _intpow( TV, n/2 ) ) * TV : sqr( _intpow( TV, n/2 ) );
2625 template <
typename T>
inline TVar<T>
2627 (
const TVar<T> &TV,
const double a )
2629 return exp( a * log( TV ) );
2632 template <
typename T>
inline TVar<T>
2634 (
const TVar<T> &TV1,
const TVar<T> &TV2 )
2636 return exp( TV2 * log( TV1 ) );
2639 template <
typename T>
inline TVar<T>
2641 (
const double a,
const TVar<T> &TV )
2643 return exp( TV * std::log( a ) );
2646 template <
typename T>
inline TVar<T>
2648 (
const unsigned int n,
const TVar<T>*TV,
const int*k)
2654 return pow( TV[0], k[0] );
2656 return pow( TV[0], k[0] ) * monomial( n-1, TV+1, k+1 );
2659 template <
typename T>
inline TVar<T>
2661 (
const TVar<T> &TV )
2664 return TVar<T>( Op<T>::cos(TV._coefmon[0] + *(TV._bndrem)) );
2666 const T I( TV.B() );
2667 double x0 = ( TV._TM->options.REF_MIDPOINT? Op<T>::mid(I):
2669 const TVar<T> TVmx0( TV - x0 );
2670 const T Imx0( I - x0 );
2673 TVar<T> TV2( TV._TM, 0. ), MON( 1. );
2674 for(
unsigned int i=1; i<=TV.nord(); i++ ){
2676 case 0: c = std::cos(x0);
break;
2677 case 1: c = -std::sin(x0);
break;
2678 case 2: c = -std::cos(x0);
break;
2680 default: c = std::sin(x0);
break;
2686 TV2._coefmon[0] += std::cos(x0);
2688 switch( (TV2.nord()+1)%4 ){
2689 case 0: TV2 += s * Op<T>::pow( Imx0, (
int)TV2.nord()+1 )
2690 * Op<T>::cos( Op<T>::zeroone()*Imx0+x0 );
break;
2691 case 1: TV2 -= s * Op<T>::pow( Imx0, (
int)TV2.nord()+1 )
2692 * Op<T>::sin( Op<T>::zeroone()*Imx0+x0 );
break;
2693 case 2: TV2 -= s * Op<T>::pow( Imx0, (
int)TV2.nord()+1 )
2694 * Op<T>::cos( Op<T>::zeroone()*Imx0+x0 );
break;
2695 case 3: TV2 += s * Op<T>::pow( Imx0, (
int)TV2.nord()+1 )
2696 * Op<T>::sin( Op<T>::zeroone()*Imx0+x0 );
break;
2699 TV2._bndord_uptd =
false;
2700 TV2._unset_bndpol();
2701 if( TV2._TM->options.CENTER_REMAINDER ) TV2._center();
2705 template <
typename T>
inline TVar<T>
2707 (
const TVar<T> &TV )
2709 return cos( TV - PI/2. );
2712 template <
typename T>
inline TVar<T>
2714 (
const TVar<T> &TV )
2717 return TVar<T>( Op<T>::asin(TV._coefmon[0] + *(TV._bndrem)) );
2718 if ( Op<T>::l(TV.B()) < -1. && Op<T>::u(TV.B()) > 1. )
2721 double x0 = ( TV._TM->options.REF_MIDPOINT? Op<T>::mid(TV.B()):
2723 double s = 1., t = 1.;
2724 TVar<T> G = TV * std::sqrt(1-x0*x0) - x0 * sqrt(1.-sqr(TV));
2725 TVar<T> TV2( G ), MON( G ), sqrG( sqr(G) );
2726 T IG( G.B() ), IG0( IG*Op<T>::zeroone() ), sqrIG0( Op<T>::sqr(IG0) ) ;
2727 T ASIN1, ASIN2( 1./Op<T>::sqrt(1-sqrIG0) ),
2728 ASIN3( IG0*Op<T>::pow(Op<T>::sqrt(1-sqrIG0),-3) );
2729 for(
unsigned int i=1; i<=TV.nord(); i++ ){
2731 s *= (double)(2*i-1)*(double)(2*i-1)/(double)(2*i)/(double)(2*i+1);
2734 ASIN1 = ASIN2; ASIN2 = ASIN3;
2735 ASIN3 = ((2*i-1)*IG0*ASIN2+(i-1)*(i-1)*ASIN1)/(1-sqrIG0);
2737 TV2._coefmon[0] += std::asin(x0);
2738 TV2 += Op<T>::pow( IG, (
int)TV2.nord()+1 ) / t * ASIN2;
2740 TV2._bndord_uptd =
false;
2741 TV2._unset_bndpol();
2742 if( TV2._TM->options.CENTER_REMAINDER ) TV2._center();
2746 template <
typename T>
inline TVar<T>
2748 (
const TVar<T> &TV )
2750 return PI/2. - asin( TV );
2753 template <
typename T>
inline TVar<T>
2755 (
const TVar<T> &TV )
2757 return sin(TV) / cos(TV);
2760 template <
typename T>
inline TVar<T>
2762 (
const TVar<T> &TV )
2765 return TVar<T>( Op<T>::atan(TV._coefmon[0] + *(TV._bndrem)) );
2767 double x0 = ( TV._TM->options.REF_MIDPOINT? Op<T>::mid(TV.B()):
2769 TVar<T> G = ( TV - x0 ) / ( 1. + x0 * TV );
2770 TVar<T> TV2( G ), MON( G ), msqrG( -sqr(G) );
2771 T IG( G.B() ), IG0( IG*Op<T>::zeroone() );
2772 for(
unsigned int i=1; i<=TV.nord(); i++ ){
2774 TV2 += MON / (2*i+1);
2776 TV2._coefmon[0] += std::atan(x0);
2777 TV2 += Op<T>::pow( IG * Op<T>::cos( Op<T>::atan(IG0) ),
2778 (
int)TV2.nord()+1 ) / (
double)(TV2.nord()+1)
2779 * sin( (TV2.nord()+1) * (atan(IG0)+PI/2.) );
2781 TV2._bndord_uptd =
false;
2782 TV2._unset_bndpol();
2783 if( TV2._TM->options.CENTER_REMAINDER ) TV2._center();
2787 template <
typename T>
inline TVar<T>
2789 (
const TVar<T>&TV1,
const TVar<T>&TV2 )
2792 if( !TV1._TM && !TV2._TM ){
2793 T R1 = TV1._coefmon[0] + *(TV1._bndrem);
2794 T R2 = TV2._coefmon[0] + *(TV2._bndrem);
2795 return Op<T>::hull(R1, R2);
2800 return hull( TV2, TV1 );
2803 else if( !TV2._TM ){
2804 TVar<T> TVR = TV1.P();
2805 return TVR + Op<T>::hull( TV1.R(), TV2._coefmon[0]+*(TV2._bndrem)-TVR.B() );
2809 else if( TV1._TM != TV2._TM )
2813 TVar<T> TV1C( TV1 ), TV2C( TV2 );
2814 const double eta = TV1._TM->options.REF_POLY;
2815 T R1C = TV1C.
C().R(), R2C = TV2C.C().R();
2818 T BTVD = (TV1C-TV2C).B();
2819 TVar<T> TVR = (1.-eta)*TV1C + eta*TV2C + Op<T>::hull( R1C+eta*BTVD, R2C+(eta-1.)*BTVD );
2823 if( TVR._TM->options.CENTER_REMAINDER ) TVR._center();
2827 template <
typename T>
inline bool
2829 ( TVar<T>&TVR,
const TVar<T>&TV1,
const TVar<T>&TV2 )
2832 if( !TV1._TM && !TV2._TM ){
2833 T R1 = TV1._coefmon[0] + *TV1._bndrem;
2834 T R2 = TV2._coefmon[0] + *TV2._bndrem;
2836 bool flag = Op<T>::inter(RR, R1, R2);
2843 return inter( TVR, TV2, TV1 );
2846 else if( !TV2._TM ){
2847 T R1 = TV1.R(), B2 = TV2.B();
2849 if( !Op<T>::inter(*(TVR._bndrem), R1, B2-TVR.B()) )
2852 if( TVR._TM->options.CENTER_REMAINDER ) TVR._center();
2857 else if( TV1._TM != TV2._TM )
2861 TVar<T> TV1C( TV1 ), TV2C( TV2 );
2862 const double eta = TV1._TM->options.REF_POLY;
2863 T R1C = TV1C.
C().R(), R2C = TV2C.C().R();
2866 TVR = (1.-eta)*TV1C + eta*TV2C;
2869 if( !Op<T>::inter( *(TVR._bndrem), R1C+eta*BTVD, R2C+(eta-1.)*BTVD ) )
2872 TVR._bndord_uptd =
false;
2873 TVR._unset_bndpol();
2874 if( TVR._TM->options.CENTER_REMAINDER ) TVR._center();
2886 template <>
template<
typename T>
struct Op< mc::
TVar<T> >
2889 static TV point(
const double c ) {
return TV(c); }
2891 static void I(TV& x,
const TV&y) { x = y; }
2892 static double l(
const TV& x) {
return mc::Op<T>::l(x.
B()); }
2893 static double u(
const TV& x) {
return mc::Op<T>::u(x.
B()); }
2897 static TV inv (
const TV& x) {
return mc::inv(x); }
2898 static TV sqr (
const TV& x) {
return mc::sqr(x); }
2899 static TV sqrt(
const TV& x) {
return mc::sqrt(x); }
2900 static TV log (
const TV& x) {
return mc::log(x); }
2901 static TV xlog(
const TV& x) {
return x*mc::log(x); }
2903 static TV exp (
const TV& x) {
return mc::exp(x); }
2904 static TV sin (
const TV& x) {
return mc::sin(x); }
2905 static TV cos (
const TV& x) {
return mc::cos(x); }
2906 static TV tan (
const TV& x) {
return mc::tan(x); }
2907 static TV asin(
const TV& x) {
return mc::asin(x); }
2908 static TV acos(
const TV& x) {
return mc::acos(x); }
2909 static TV atan(
const TV& x) {
return mc::atan(x); }
2914 static TV hull(
const TV& x,
const TV& y) {
return mc::hull(x,y); }
2915 static TV min (
const TV& x,
const TV& y) {
return mc::Op<T>::min(x.
B(),y.
B()); }
2916 static TV max (
const TV& x,
const TV& y) {
return mc::Op<T>::max(x.
B(),y.
B()); }
2917 static TV arh (
const TV& x,
const double k) {
return mc::exp(-k/x); }
2918 template <
typename X,
typename Y>
static TV pow(
const X& x,
const Y& y) {
return mc::pow(x,y); }
2919 static TV monomial (
const unsigned int n,
const T* x,
const int* k) {
return mc::monomial(n,x,k); }
2920 static bool inter(TV& xIy,
const TV& x,
const TV& y) {
return mc::inter(xIy,x,y); }
2921 static bool eq(
const TV& x,
const TV& y) {
return mc::Op<T>::eq(x.
B(),y.
B()); }
2922 static bool ne(
const TV& x,
const TV& y) {
return mc::Op<T>::ne(x.
B(),y.
B()); }
2923 static bool lt(
const TV& x,
const TV& y) {
return mc::Op<T>::lt(x.
B(),y.
B()); }
2924 static bool le(
const TV& x,
const TV& y) {
return mc::Op<T>::le(x.
B(),y.
B()); }
2925 static bool gt(
const TV& x,
const TV& y) {
return mc::Op<T>::gt(x.
B(),y.
B()); }
2926 static bool ge(
const TV& x,
const TV& y) {
return mc::Op<T>::ge(x.
B(),y.
B()); }
~TVar()
Destructor of Taylor variable.
Definition: tmodel.hpp:918
C++ class for Taylor model computation of factorable function - Taylor model propagation.
Definition: tmodel.hpp:395
bool CENTER_REMAINDER
Whether to center the remainder term during Taylor model propagation.
Definition: tmodel.hpp:572
const T * bndmon() const
Const pointer to array of size nmon() with bounds on each monomial term in Taylor model...
Definition: tmodel.hpp:439
Eigenvalue decomposition-based bounder.
Definition: tmodel.hpp:559
Exceptions of mc::TModel.
Definition: tmodel.hpp:459
int ierr()
Error flag.
Definition: tmodel.hpp:479
unsigned int nmon() const
Total number of monomial terms in polynomial model.
Definition: polymodel.hpp:463
Lin & Stadtherr range bounder.
Definition: tmodel.hpp:558
C++ base class for the computation of polynomial models for factorable functions: Environment...
Definition: polymodel.hpp:36
unsigned int BOUNDER_ORDER
Order of Bernstein polynomial for Taylor model range bounding (no less than Taylor model order!)...
Definition: tmodel.hpp:566
std::pair< unsigned int, const double * > coefmon() const
Get pair of size of, and const pointer to, array of (possibly scaled) monomial coefficients in multiv...
Definition: polymodel.hpp:870
bool REF_MIDPOINT
Whether to take the midpoint of the inner range as the reference in the outer composition with a univ...
Definition: tmodel.hpp:574
TVar< T > & C()
Shortcut to mc::TVar::center.
Definition: tmodel.hpp:980
TModel< T > * env() const
Get pointer to linked Taylor model environment.
Definition: tmodel.hpp:876
std::string what()
Error description.
Definition: tmodel.hpp:481
double constant(const bool reset=false)
Get coefficient of constant term in Taylor variable. The value of this coefficient is reset to 0 if r...
Definition: tmodel.hpp:1754
virtual void _center()
Center remainder error term _bndrem
Definition: polymodel.hpp:850
TVar< T > P() const
Shortcut to mc::TVar::polynomial.
Definition: tmodel.hpp:971
Options of mc::TModel.
Definition: tmodel.hpp:515
Square-root operation with negative numbers in range.
Definition: tmodel.hpp:467
const unsigned int *const * prodmon() const
Const pointer to double array of size (_nmon+1,<=_nmon) with indices of monomial terms from product o...
Definition: tmodel.hpp:443
bool SCALE_VARIABLES
Whether to scale the variable ranges to [-1,1] internally.
Definition: tmodel.hpp:570
~TModel()
Destructor of Taylor model environment.
Definition: tmodel.hpp:435
T * _bndrem
Pointer to remainder bound of variable (possibly NULL if not computed)
Definition: polymodel.hpp:443
Options()
Constructor of mc::TModel::Options.
Definition: tmodel.hpp:518
bool BERNSTEIN_USE
Whether to compute a Berstein model [Stancu, 1963] of the outer function in a univariate composition ...
Definition: tmodel.hpp:578
Operation between Taylor variables linked to different Taylor models.
Definition: tmodel.hpp:473
void reset()
Reset the bounds on powers of (possibly scaled) variable ranges.
Definition: tmodel.hpp:455
C++ class for the computation of Taylor models of factorable function - Taylor model environment...
Definition: tmodel.hpp:406
BOUNDER
Taylor model range bounder option.
Definition: tmodel.hpp:556
Hybrid LSB + EIGEN range bounder.
Definition: tmodel.hpp:561
int BOUNDER_TYPE
Taylor model range bounder - See How are the options set for the computation of a Taylor model...
Definition: tmodel.hpp:564
Exceptions(TYPE ierr)
Constructor for error ierr
Definition: tmodel.hpp:477
unsigned int BERNSTEIN_MAXIT
Maximum number of iterations for determination of the exact remainder bounds in a Berstein model of c...
Definition: tmodel.hpp:582
const double * reference() const
Get const pointer to array of size nvar with references for all variables in Taylor model...
Definition: tmodel.hpp:447
Sine/Cosine inverse operation with range outside [-1,1].
Definition: tmodel.hpp:468
TVar< T > & center()
Center remainder term of Taylor variable.
Definition: tmodel.hpp:976
const double * scaling() const
Get const pointer to array of size nvar with scaling coefficients in Taylor model.
Definition: tmodel.hpp:451
Log operation with non-positive numbers in range.
Definition: tmodel.hpp:466
virtual T B() const
Shortcut to mc::PolyVar::bound.
Definition: polymodel.hpp:651
bool BERNSTEIN_OPT
Whether to compute exact remainder bounds for Berstein models of convex/concave univariates exp...
Definition: tmodel.hpp:580
Division by zero scalar.
Definition: tmodel.hpp:464
virtual PolyVar< T > & set(PolyModel *env)
Set polynomial model environment in variable as env
Definition: polymodel.hpp:573
static const std::string BOUNDER_NAME[5]
Array of Taylor model range bounder names (for display)
Definition: tmodel.hpp:568
C++ base class for the computation of polynomial models for factorable functions: Variable...
Definition: polymodel.hpp:425
virtual T B(const int type) const
Shortcut to mc::PolyVar::bound.
Definition: polymodel.hpp:645
unsigned int nord() const
Order of polynomial model.
Definition: polymodel.hpp:451
virtual void _set_bndpol(const T &bndpol)
Set polynomial bound in variable as bndpol
Definition: polymodel.hpp:841
Failed to construct Taylor variable.
Definition: tmodel.hpp:471
bool _bndord_uptd
Whether the bounds in bndord are up-to-date.
Definition: polymodel.hpp:440
unsigned int nord() const
Order of polynomial model environment.
Definition: polymodel.hpp:62
unsigned int DISPLAY_DIGITS
Number of digits in output stream for Taylor model coefficients.
Definition: tmodel.hpp:586
double BERNSTEIN_TOL
Termination tolerance for determination of the exact remainder bounds in a Berstein model of convex/c...
Definition: tmodel.hpp:584
TVar< T > & set(TModel< T > *TM, const unsigned int ix, const T &X, const double Xref)
Set Taylor variable with index ix (starting from 0), bounded by X and with reference point Xref...
Definition: tmodel.hpp:941
Inverse operation with zero in range.
Definition: tmodel.hpp:465
Feature not yet implemented in mc::TModel.
Definition: tmodel.hpp:474
PolyModel(const unsigned int nvar, const unsigned int nord)
Constructor of polynomial model environment for nvar variables and order nord
Definition: polymodel.hpp:46
TYPE
Enumeration type for TModel exception handling.
Definition: tmodel.hpp:463
double * linear() const
Get pointer to array of size nvar with coefficients of linear term in Taylor variable.
Definition: tmodel.hpp:1762
double * _coefmon
Array of size _nmon with monomial coefficients of variable.
Definition: polymodel.hpp:434
Taylor model bound does not intersect with bound in template parameter arithmetic.
Definition: tmodel.hpp:472
Naive polynomial range bounder.
Definition: tmodel.hpp:557
T * _bndord
Array of size _nord+2 with bounds for all terms of degrees iord=0,...,_nord as well as the remainder ...
Definition: polymodel.hpp:437
TVar< T > polynomial() const
Return new Taylor variable with same multivariate polynomial part but zero remainder.
Definition: tmodel.hpp:965
unsigned int nvar() const
Number of variables in polynomial model environment.
Definition: polymodel.hpp:56
double REF_POLY
Scalar in related to the choice of the polynomial part in the overloaded functions mc::inter and mc:...
Definition: tmodel.hpp:576
Bernstein range bounder.
Definition: tmodel.hpp:560
Failed to compute the maximum gap between a univariate term and its Bernstein model.
Definition: tmodel.hpp:470
C++ structure to allow usage of various template parameters in the types mc::McCormick, mc::TModel, mc::TVar, and mc::SpecBnd of MC++ _ Specialization of this structure is required for the template parameters can be found in the header files defining the types mc::McCormick, mc::TModel, mc::TVar, and mc::SpecBnd.
Definition: mcop.hpp:13
Failed to compute eigenvalue decomposition in range bounder TModel::Options::EIGEN.
Definition: tmodel.hpp:469