293 #ifndef MC__SPECBND_H
294 #define MC__SPECBND_H
299 #include "mclapack.hpp"
304 #undef MC__SPECBND_DEBUG_SPECTRUM
305 #undef MC__SPECBND_DEBUG_HESSBND
316 template <
typename T>
320 template <
typename U>
friend class Specbnd;
355 template <
class U>
friend std::ostream& operator<<(std::ostream&, const Specbnd<U>&);
366 static T _LambdaS(
const fadbad::F<T> &a,
const unsigned int n );
368 static T _LambdaT(
const fadbad::F<T> &a,
const fadbad::F<T> &b,
const unsigned int n );
371 static std::pair<double,double> _gershgorin_bound
372 (
const fadbad::B< fadbad::F< T > >* D2X );
374 static std::pair<double,double> _hertzrohn_bound
375 (
const fadbad::B< fadbad::F< T > >* D2X );
377 static std::pair<double,double> _gershgorin_bound
378 (
const fadbad::F< fadbad::F< T > >& D2F );
380 static std::pair<double,double> _hertzrohn_bound
381 (
const fadbad::F< fadbad::F< T > >& D2F );
446 return "mc::Specbnd\t Computation of Hessian spectrum failed";
448 return "mc::Specbnd\t Computation of Hessian spectral bound failed";
450 return "mc::Specbnd\t Operation between variables with different numbers of dependents";
452 return "mc::Specbnd\t Feature not yet implemented in mc::Specbnd class";
454 return "mc::Specbnd\t Undocumented error";
464 _n(0), _FI(0), _spec(0.)
469 _n(0), _FI(c), _spec(0.)
474 _n(0), _FI(B), _spec(0.)
478 Specbnd(
const T &B,
const unsigned int i,
const unsigned int n ):
479 _n(n), _FI(B), _spec(0.)
486 _n(x._n), _FI(x._FI), _spec(x._spec)
528 const fadbad::F<T>&
FI()
const
540 static std::pair<double,double>
spectrum
541 (
const fadbad::B< fadbad::F< double > >* D2X );
544 static std::pair<double,double>
spectrum
545 (
const fadbad::F< fadbad::F< double > >& D2X );
549 (
const fadbad::B< fadbad::F< T > >* D2X );
553 (
const fadbad::F< fadbad::F< T > >& D2X );
559 template <
typename T>
typename Specbnd<T>::Options Specbnd<T>::options;
561 template <
class T>
inline T
563 (
const fadbad::F<T> &a,
unsigned int n )
566 if( n == 1 )
return a.size()? Op<T>::sqr( a[0] ): 0.;
568 for (
unsigned int i=0; i<n; i++)
569 if( a.size() ) upbnd += Op<T>::u( Op<T>::sqr( a[i] ) );
570 return Op<T>::zeroone() * upbnd;
573 template <
class T>
inline T
575 (
const fadbad::F<T> &a,
const fadbad::F<T> &b,
unsigned int n )
578 if( n == 1 )
return (a.size() && b.size())? 2.*a[0]*b[0]: 0.;
579 double upbnda=0., upbndb=0.;
580 for (
unsigned int i=0; i<n; i++){
581 if( a.size() ) upbnda += Op<T>::u( Op<T>::sqr( a[i] ) );
582 if( b.size() ) upbndb += Op<T>::u( Op<T>::sqr( b[i] ) );
584 T lamb = 2.*(Op<T>::zeroone()-0.5) * std::sqrt(upbnda*upbndb);
585 for (
unsigned int i=0; i<n; i++)
586 if( a.size() && b.size() ) lamb += a[i]*b[i];
590 template <
typename T>
inline std::pair<double,double>
592 (
const fadbad::B< fadbad::F< double > >* D2X )
595 const unsigned int N = D2X->val().size();
597 double*H =
new double[N*N];
599 for(
unsigned int j=0; j<N; j++ )
600 for(
unsigned int i=j; i<N; i++ )
601 H[j*N+i] = H[i*N+j] = ( i==j? D2X[i].deriv(0).deriv(i):
602 0.5*(D2X[i].deriv(0).deriv(j)+D2X[j].deriv(0).deriv(i)) );
603 #ifdef MC__SPECBND_DEBUG_SPECTRUM
604 mc::display( N, N, H, N,
"\nMatrix H", std::cout );
607 double*D = mc::dsyev_wrapper( N, H );
610 std::pair<double,double> spbnd = std::make_pair( min(N,D), max(N,D) );
615 template <
typename T>
inline std::pair<double,double>
617 (
const fadbad::F< fadbad::F< double > >& D2F )
619 const unsigned int N = D2F.size();
621 double*H =
new double[N*N];
623 for(
unsigned int j=0; j<N; j++ )
624 for(
unsigned int i=j; i<N; i++ )
625 H[j*N+i] = H[i*N+j] = ( i==j? D2F.deriv(i).deriv(i):
626 0.5*(D2F.deriv(i).deriv(j)+D2F.deriv(j).deriv(i)) );
627 #ifdef MC__SPECBND_DEBUG_SPECTRUM
628 mc::display( N, N, H, N,
"\nMatrix H", std::cout );
631 double*D = mc::dsyev_wrapper( N, H );
634 std::pair<double,double> spbnd = std::make_pair( min(N,D), max(N,D) );
639 template <
typename T>
inline std::pair<double,double>
641 (
const fadbad::B< fadbad::F< T > >* D2X )
644 switch( options.HESSBND ){
645 case Options::GERSHGORIN:
646 return _gershgorin_bound( D2X );
647 case Options::HERTZROHN:
default:
648 return _hertzrohn_bound( D2X );
652 template <
typename T>
inline std::pair<double,double>
654 (
const fadbad::B< fadbad::F< T > >* D2X )
656 const unsigned int N = D2X->val().size();
659 std::pair<double,double> spbnd;
660 for(
unsigned int i=0; i<N; i++ ){
662 for(
unsigned int j=0; j<N; j++ ){
663 if( j == i )
continue;
665 if( !inter( D2Fij, D2X[i].deriv(0).deriv(j), D2X[i].deriv(0).deriv(j) ) )
666 D2Fij = D2X[i].deriv(0).deriv(j);
671 spbnd.first = Op<T>::l(D2X[i].deriv(0).deriv(i)) - ri;
672 spbnd.second = Op<T>::u(D2X[i].deriv(0).deriv(i)) + ri;
675 spbnd.first = std::min( spbnd.first, Op<T>::l(D2X[i].deriv(0).deriv(i)) - ri );
676 spbnd.second = std::max( spbnd.second, Op<T>::u(D2X[i].deriv(0).deriv(i)) + ri );
682 template <
typename T>
inline std::pair<double,double>
683 Specbnd<T>::_hertzrohn_bound
684 (
const fadbad::B< fadbad::F< T > >* D2X )
686 const unsigned int N = D2X->val().size();
690 unsigned int col_S = 2;
691 short *S =
new short[col_S];
693 #ifdef MC__SPECBND_DEBUG_HESSBND
694 mc::display( 1, col_S, S, 1,
"\nMatrix S1", std::cout );
696 for(
unsigned int k=1; k<N; k++, col_S*=2 ){
697 short *Sprev =
new short[k*col_S];
698 for(
unsigned int i=0; i<k; i++ )
699 for(
unsigned int j=0; j<col_S; j++ )
700 Sprev[j*k+i] = S[j*k+i];
701 delete[] S; S =
new short[(k+1)*2*col_S];
702 for(
unsigned int i=0; i<k; i++ ){
703 for(
unsigned int j=0; j<col_S; j++ )
704 S[j*(k+1)+i] = S[(col_S+j)*(k+1)+i] = Sprev[j*k+i];
706 for(
unsigned int j=0; j<col_S; j++ ){
707 S[j*(k+1)+k] = 1; S[(col_S+j)*(k+1)+k] = -1;
710 #ifdef MC__SPECBND_DEBUG_HESSBND
711 mc::display( k+1, 2*col_S, S, k+1,
"\nMatrix Sk", std::cout );
716 std::pair<double,double> spbnd;
717 double *Lk =
new double[N*N], *Uk =
new double[N*N];
718 for(
unsigned int k=0; k<col_S; k++ ){
720 for(
unsigned int j=0; j<N; j++ )
721 for(
unsigned int i=j; i<N; i++ ){
723 Lk[i*N+i] = Op<T>::l( D2X[i].deriv(0).deriv(i) );
724 Uk[i*N+i] = Op<T>::u( D2X[i].deriv(0).deriv(i) );
728 if( !inter( D2Fij, D2X[i].deriv(0).deriv(j), D2X[j].deriv(0).deriv(i) ) )
729 D2Fij = D2X[i].deriv(0).deriv(j);
730 Lk[i*N+j] = Lk[j*N+i] = S[k*N+i]*S[k*N+j]==1? Op<T>::l( D2Fij ): Op<T>::u( D2Fij );
731 Uk[i*N+j] = Uk[j*N+i] = S[k*N+i]*S[k*N+j]==1? Op<T>::u( D2Fij ): Op<T>::l( D2Fij );
733 #ifdef MC__SPECBND_DEBUG_HESSBND
734 mc::display( N, N, Lk, N,
"\nMatrix Lk", std::cout );
735 mc::display( N, N, Uk, N,
"\nMatrix Uk", std::cout );
738 double*DLk = mc::dsyev_wrapper( N, Lk );
740 delete[] Lk;
delete[] Uk;
delete[] S;
743 #ifdef MC__SPECBND_DEBUG_HESSBND
744 mc::display( 1, N, DLk, 1,
"\nMatrix DLk", std::cout );
746 spbnd.first = k? std::min( spbnd.first, min(N,DLk) ): min(N,DLk);
749 double*DUk = mc::dsyev_wrapper( N, Uk );
751 delete[] Lk;
delete[] Uk;
delete[] S;
754 #ifdef MC__SPECBND_DEBUG_HESSBND
755 mc::display( 1, N, DUk, 1,
"\nMatrix DUk", std::cout );
757 spbnd.second = k? std::max( spbnd.second, max(N,DUk) ): max(N,DUk);
761 delete[] Lk;
delete[] Uk;
delete[] S;
765 template <
typename T>
inline std::pair<double,double>
767 (
const fadbad::F< fadbad::F< T > >& D2F )
769 switch( options.HESSBND ){
770 case Options::GERSHGORIN:
771 return _gershgorin_bound( D2F );
772 case Options::HERTZROHN:
default:
773 return _hertzrohn_bound( D2F );
777 template <
typename T>
inline std::pair<double,double>
779 (
const fadbad::F< fadbad::F< T > >& D2F )
781 const unsigned int N = D2F.size();
784 std::pair<double,double> spbnd;
785 for(
unsigned int i=0; i<N; i++ ){
787 for(
unsigned int j=0; j<N; j++ ){
788 if( j == i )
continue;
790 if( !inter( D2Fij, D2F.deriv(i).deriv(j), D2F.deriv(i).deriv(j) ) )
791 D2Fij = D2F.deriv(i).deriv(j);
796 spbnd.first = Op<T>::l(D2F.deriv(i).deriv(i)) - ri;
797 spbnd.second = Op<T>::u(D2F.deriv(i).deriv(i)) + ri;
800 spbnd.first = std::min( spbnd.first, Op<T>::l(D2F.deriv(i).deriv(i)) - ri );
801 spbnd.second = std::max( spbnd.second, Op<T>::u(D2F.deriv(i).deriv(i)) + ri );
807 template <
typename T>
inline std::pair<double,double>
808 Specbnd<T>::_hertzrohn_bound
809 (
const fadbad::F< fadbad::F< T > >& D2F )
811 const unsigned int N = D2F.size();
815 unsigned int col_S = 2;
816 short *S =
new short[col_S];
818 #ifdef MC__SPECBND_DEBUG_HESSBND
819 mc::display( 1, col_S, S, 1,
"\nMatrix S1", std::cout );
821 for(
unsigned int k=1; k<N; k++, col_S*=2 ){
822 short *Sprev =
new short[k*col_S];
823 for(
unsigned int i=0; i<k; i++ )
824 for(
unsigned int j=0; j<col_S; j++ )
825 Sprev[j*k+i] = S[j*k+i];
826 delete[] S; S =
new short[(k+1)*2*col_S];
827 for(
unsigned int i=0; i<k; i++ ){
828 for(
unsigned int j=0; j<col_S; j++ )
829 S[j*(k+1)+i] = S[(col_S+j)*(k+1)+i] = Sprev[j*k+i];
831 for(
unsigned int j=0; j<col_S; j++ ){
832 S[j*(k+1)+k] = 1; S[(col_S+j)*(k+1)+k] = -1;
835 #ifdef MC__SPECBND_DEBUG_HESSBND
836 mc::display( k+1, 2*col_S, S, k+1,
"\nMatrix Sk", std::cout );
840 #ifdef MC__SPECBND_DEBUG_HESSBND
842 for(
unsigned int j=0; j<N; j++ )
843 for(
unsigned int i=j; i<N; i++ ){
845 H[i*N+i] = D2F.deriv(i).deriv(i);
849 if( !inter( D2Fij, D2F.deriv(i).deriv(j), D2F.deriv(j).deriv(i) ) )
850 D2Fij = D2F.deriv(i).deriv(j);
851 H[i*N+j] = H[j*N+i] = D2Fij;
853 mc::display( N, N, H, N,
"\nMatrix H", std::cout );
858 std::pair<double,double> spbnd;
859 double *Lk =
new double[N*N], *Uk =
new double[N*N];
860 for(
unsigned int k=0; k<col_S; k++ ){
862 for(
unsigned int j=0; j<N; j++ )
863 for(
unsigned int i=j; i<N; i++ ){
865 Lk[i*N+i] = Op<T>::l( D2F.deriv(i).deriv(i) );
866 Uk[i*N+i] = Op<T>::u( D2F.deriv(i).deriv(i) );
870 if( !inter( D2Fij, D2F.deriv(i).deriv(j), D2F.deriv(j).deriv(i) ) )
871 D2Fij = D2F.deriv(i).deriv(j);
872 Lk[i*N+j] = Lk[j*N+i] = S[k*N+i]*S[k*N+j]==1? Op<T>::l( D2Fij ): Op<T>::u( D2Fij );
873 Uk[i*N+j] = Uk[j*N+i] = S[k*N+i]*S[k*N+j]==1? Op<T>::u( D2Fij ): Op<T>::l( D2Fij );
875 #ifdef MC__SPECBND_DEBUG_HESSBND
876 mc::display( N, N, Lk, N,
"\nMatrix Lk", std::cout );
877 mc::display( N, N, Uk, N,
"\nMatrix Uk", std::cout );
880 double*DLk = mc::dsyev_wrapper( N, Lk );
882 delete[] Lk;
delete[] Uk;
delete[] S;
885 #ifdef MC__SPECBND_DEBUG_HESSBND
886 mc::display( 1, N, DLk, 1,
"\nMatrix DLk", std::cout );
888 spbnd.first = k? std::min( spbnd.first, min(N,DLk) ): min(N,DLk);
891 double*DUk = mc::dsyev_wrapper( N, Uk );
893 delete[] Lk;
delete[] Uk;
delete[] S;
896 #ifdef MC__SPECBND_DEBUG_HESSBND
897 mc::display( 1, N, DUk, 1,
"\nMatrix DUk", std::cout );
899 spbnd.second = k? std::max( spbnd.second, max(N,DUk) ): max(N,DUk);
903 delete[] Lk;
delete[] Uk;
delete[] S;
907 template <
class T>
inline std::ostream&
909 ( std::ostream &out,
const Specbnd<T> &y )
911 out <<
" " << y._spec << std::endl
912 <<
" " << y._FI.val() << std::endl;
913 for(
unsigned int i=0; i<y._n; i++ )
914 out <<
" ( " << y._FI.deriv(i) <<
" )" << std::endl;
918 template <
class T>
inline Specbnd<T>&
919 Specbnd<T>::operator=
928 template <
class T>
inline Specbnd<T>&
929 Specbnd<T>::operator=
938 template <
class T>
inline Specbnd<T>&
940 (
const fadbad::F<T>&FB,
const T&SB )
958 template <
class T>
inline Specbnd<T>&
959 Specbnd<T>::operator+=
966 template <
class T>
inline Specbnd<T>&
967 Specbnd<T>::operator+=
968 (
const Specbnd<T> &y )
970 if( _n && y._n && _n != y._n )
978 template <
class T>
inline Specbnd<T>
980 (
const Specbnd<T> &y )
985 template <
class T> Specbnd<T>
987 (
const Specbnd<T> &x,
const Specbnd<T> &y )
989 if( x._n && y._n && x._n != y._n )
992 z._FI = x._FI + y._FI;
994 z._spec = x._spec + y._spec;
998 template <
class T> Specbnd<T>
1000 (
const Specbnd<T> &y,
const double c )
1009 template <
class T> Specbnd<T>
1011 (
const double c,
const Specbnd<T> &y )
1016 template <
class T>
inline Specbnd<T>&
1017 Specbnd<T>::operator-=
1024 template <
class T>
inline Specbnd<T>&
1025 Specbnd<T>::operator-=
1026 (
const Specbnd<T> &y )
1028 if( _n && y._n && _n != y._n )
1036 template <
class T>
inline Specbnd<T>
1038 (
const Specbnd<T> &y )
1047 template <
class T>
inline Specbnd<T>
1049 (
const Specbnd<T> &x,
const Specbnd<T> &y )
1051 if( x._n && y._n && x._n != y._n )
1054 z._FI = x._FI - y._FI;
1055 z._n = z._FI.size();
1056 z._spec = x._spec - y._spec;
1060 template <
class T>
inline Specbnd<T>
1062 (
const Specbnd<T> &y,
const double c )
1067 template <
class T>
inline Specbnd<T>
1069 (
const double c,
const Specbnd<T> &y )
1074 template <
typename T>
inline Specbnd<T>&
1075 Specbnd<T>::operator*=
1078 Specbnd<T> z = c * (*this);
1083 template <
typename T>
inline Specbnd<T>&
1084 Specbnd<T>::operator*=
1085 (
const Specbnd<T> &x )
1087 Specbnd<T> z = x * (*this);
1092 template <
class T>
inline Specbnd<T>
1094 (
const Specbnd<T> &x,
const Specbnd<T> &y )
1096 if( x._n && y._n && x._n != y._n )
1099 z._FI = x._FI * y._FI;
1100 z._n = z._FI.size();
1101 z._spec = x._FI.val() * y._spec + x._spec * y._FI.val()
1102 + Specbnd<T>::_LambdaT( x._FI, y._FI, z._n );
1106 template <
class T>
inline Specbnd<T>
1108 (
const double c,
const Specbnd<T> &y )
1113 z._spec = c * y._spec;
1117 template <
class T>
inline Specbnd<T>
1119 (
const Specbnd<T> &y,
const double c )
1124 template <
class T>
inline Specbnd<T>
1126 (
const Specbnd<T> &y )
1130 z._FI = fadbad::sqr( y._FI );
1131 z._spec = 2.*( y._FI.val() * y._spec + Specbnd<T>::_LambdaS( y._FI, z._n ) );
1135 template <
class T>
inline Specbnd<T>
1137 (
const Specbnd<T> &y,
const int m )
1140 if( m == 1 )
return y;
1141 if( m == 2 )
return sqr(y);
1142 if( m == -1 )
return inv( y );
1143 if( m < -1 )
return inv( pow( y, -m ) );
1146 z._FI = fadbad::pow( y._FI, m );
1147 z._spec = (double)m * Op<T>::pow( y._FI.val(), m-2 )
1148 * ( y._FI.val() * y._spec + (m-1) * Specbnd<T>::_LambdaS( y._FI, z._n ) );
1152 template <
class T>
inline Specbnd<T>
1154 (
const Specbnd<T> &y )
1159 z._spec = Op<T>::sqr( z._FI.val() )
1160 * ( 2. * z._FI.val() * Specbnd<T>::_LambdaS( y._FI, z._n ) - y._spec );
1164 template <
typename T>
inline Specbnd<T>&
1165 Specbnd<T>::operator/=
1171 template <
typename T>
inline Specbnd<T>&
1172 Specbnd<T>::operator/=
1173 (
const Specbnd<T> &x )
1178 template <
class T>
inline Specbnd<T>
1180 (
const Specbnd<T> &x,
const Specbnd<T> &y )
1185 template <
class T>
inline Specbnd<T>
1187 (
const Specbnd<T> &y,
const double c )
1192 template <
class T>
inline Specbnd<T>
1194 (
const double c,
const Specbnd<T> &y )
1199 template <
class T>
inline Specbnd<T>
1201 (
const Specbnd<T> &y )
1205 z._FI = fadbad::sqrt( y._FI );
1206 z._spec = 1./(2.*Op<T>::sqr( z._FI.val() ))
1207 * ( y._spec - Specbnd<T>::_LambdaS( y._FI, z._n ) / (2. * y._FI.val()) );
1211 template <
class T>
inline Specbnd<T>
1213 (
const Specbnd<T> &y )
1217 z._FI = fadbad::exp( y._FI );
1218 z._spec = z._FI.val() * ( y._spec + Specbnd<T>::_LambdaS( y._FI, z._n ) );
1222 template <
class T>
inline Specbnd<T>
1224 (
const Specbnd<T> &y )
1228 z._FI = fadbad::log( y._FI );
1229 z._spec = 1./y._FI.val() * ( y._spec - Specbnd<T>::_LambdaS( y._FI, z._n ) / y._FI.val() );
1233 template <
class T>
inline Specbnd<T>
1235 (
const Specbnd<T> &y )
1239 z._FI = y._FI*fadbad::log( y._FI );
1240 z._spec = (Op<T>::log(y._FI.val())+1.)*y._spec
1241 + Specbnd<T>::_LambdaS( y._FI, z._n ) / y._FI.val();
1245 template <
class T>
inline Specbnd<T>
1247 (
const Specbnd<T> &x,
const Specbnd<T> &y )
1249 return exp( y * log( x ) );
1252 template <
class T>
inline Specbnd<T>
1254 (
const double c,
const Specbnd<T> &y )
1256 return exp( y * std::log( c ) );
1259 template <
class T>
inline Specbnd<T>
1261 (
const Specbnd<T> &x,
const double c )
1263 return exp( c * log( x ) );
1266 template <
class T>
inline Specbnd<T>
1268 (
const unsigned int n,
const Specbnd<T>*x,
const int*k )
1274 return pow( x[0], k[0] );
1276 return pow( x[0], k[0] ) * monomial( n-1, x+1, k+1 );
1279 template <
class T>
inline Specbnd<T>
1281 (
const Specbnd<T> &y )
1286 template <
class T>
inline Specbnd<T>
1288 (
const Specbnd<T> &y )
1292 z._FI = fadbad::cos( y._FI );
1293 z._spec = - z._FI.val() * ( y._spec + Specbnd<T>::_LambdaS( y._FI, z._n ) );
1297 template <
class T>
inline Specbnd<T>
1299 (
const Specbnd<T> &y )
1303 z._FI = fadbad::sin( y._FI );
1304 z._spec = z._FI.val() * ( y._spec - Specbnd<T>::_LambdaS( y._FI, z._n ) );
1308 template <
class T>
inline Specbnd<T>
1310 (
const Specbnd<T> &y )
1314 z._FI = fadbad::tan( y._FI );
1315 z._spec = ( Op<T>::sqr(z._FI.val()) + 1. ) * ( y._spec
1316 + 2. * y._FI.val() * Specbnd<T>::_LambdaS( y._FI, z._n ) );
1320 template <
class T>
inline Specbnd<T>
1322 (
const Specbnd<T> &y )
1326 z._FI = fadbad::acos( y._FI );
1327 z._spec = - 1./Op<T>::sqrt(1.-Op<T>::sqr(y._FI.val())) * ( y._spec
1328 + y._FI.val()/(1.-Op<T>::sqr(y._FI.val()))
1329 *Specbnd<T>::_LambdaS( y._FI, z._n ) );
1333 template <
class T>
inline Specbnd<T>
1335 (
const Specbnd<T> &y )
1339 z._FI = fadbad::asin( y._FI );
1340 z._spec = 1./Op<T>::sqrt(1.-Op<T>::sqr(y._FI.val())) * ( y._spec
1341 + y._FI.val()/(1.-Op<T>::sqr(y._FI.val()))
1342 *Specbnd<T>::_LambdaS( y._FI, z._n ) );
1346 template <
class T>
inline Specbnd<T>
1348 (
const Specbnd<T> &y )
1352 z._FI = fadbad::atan( y._FI );
1353 z._spec = 1./(Op<T>::sqr(y._FI.val())+1.) * ( y._spec
1354 - 2.*y._FI.val()/(Op<T>::sqr(y._FI.val())+1.)
1355 *Specbnd<T>::_LambdaS( y._FI, z._n ) );
1359 template <
class T>
inline Specbnd<T>
1361 (
const Specbnd<T> &y )
1372 template <
class T>
inline Specbnd<T>
1374 (
const Specbnd<T> &y )
1394 template <>
template<
typename T>
struct Op< mc::
Specbnd<T> >
1397 static SB point(
const double c ) {
return SB(c); }
1399 static void I(SB& x,
const SB&y) { x = y; }
1400 static double l(
const SB& x) {
return Op<T>::l(x.
I()); }
1401 static double u(
const SB& x) {
return Op<T>::u(x.
I()); }
1405 static SB inv (
const SB& x) {
return mc::inv(x); }
1406 static SB sqr (
const SB& x) {
return mc::sqr(x); }
1407 static SB sqrt(
const SB& x) {
return mc::sqrt(x); }
1408 static SB log (
const SB& x) {
return mc::log(x); }
1409 static SB xlog(
const SB& x) {
return x*mc::log(x); }
1410 static SB fabs(
const SB& x) {
return mc::fabs(x); }
1411 static SB exp (
const SB& x) {
return mc::exp(x); }
1412 static SB sin (
const SB& x) {
return mc::sin(x); }
1413 static SB cos (
const SB& x) {
return mc::cos(x); }
1414 static SB tan (
const SB& x) {
return mc::tan(x); }
1415 static SB asin(
const SB& x) {
return mc::asin(x); }
1416 static SB acos(
const SB& x) {
return mc::acos(x); }
1417 static SB atan(
const SB& x) {
return mc::atan(x); }
1425 static SB arh (
const SB& x,
const double k) {
return mc::exp(-k/x); }
1426 template <
typename X,
typename Y>
static SB pow(
const X& x,
const Y& y) {
return mc::pow(x,y); }
1427 static SB monomial (
const unsigned int n,
const T* x,
const int* k) {
return mc::monomial(n,x,k); }
1429 static bool eq(
const SB& x,
const SB& y) {
return x.
SI()==y.
SI() && x.
FI()==y.
FI(); }
1430 static bool ne(
const SB& x,
const SB& y) {
return x.
SI()!=y.
SI() || x.
FI()==y.
FI(); }
1431 static bool lt(
const SB& x,
const SB& y) {
return x.
SI()<y.
SI() && x.
FI()<y.
FI(); }
1432 static bool le(
const SB& x,
const SB& y) {
return x.
SI()<=y.
SI() && x.
FI()<=y.
FI(); }
1433 static bool gt(
const SB& x,
const SB& y) {
return x.
SI()>y.
SI() && x.
FI()>y.
FI(); }
1434 static bool ge(
const SB& x,
const SB& y) {
return x.
SI()>=y.
SI() && x.
FI()>=y.
FI(); }
const fadbad::F< T > & FI() const
Return function and gradient bounds.
Definition: specbnd.hpp:528
TYPE
Enumeration type for Specbnd exception handling.
Definition: specbnd.hpp:432
Failed to compute spectrum in Specbnd::spectrum.
Definition: specbnd.hpp:433
HESSBND_STRATEGY HESSBND
Method to bound eignevalues in interval Hessian matrix using mc::Specbnd::spectral_bound.
Definition: specbnd.hpp:424
std::string what()
Error description.
Definition: specbnd.hpp:443
unsigned int dep() const
Return number of independent variables.
Definition: specbnd.hpp:516
int ierr()
Inline function returning the error flag.
Definition: specbnd.hpp:441
HESSBND_STRATEGY
Strategy for computing spectral bounds in interval Hessian matrix.
Definition: specbnd.hpp:419
Exceptions(TYPE ierr)
Constructor for error ierr
Definition: specbnd.hpp:439
static std::pair< double, double > spectrum(const fadbad::B< fadbad::F< double > > *D2X)
Compute spectrum of Hessian matrix D2X of type fadbad::B< fadbad::F<double> >* using LAPACK function ...
Definition: specbnd.hpp:592
Feature not yet implemented in mc::Specbnd.
Definition: specbnd.hpp:436
static std::pair< double, double > spectral_bound(const fadbad::B< fadbad::F< T > > *D2X)
Compute spectral bound of interval Hessian matrix D2X of type fadbad::B< fadbad::F<T> >*...
Definition: specbnd.hpp:641
Options of mc::Specbnd.
Definition: specbnd.hpp:412
~Specbnd()
Destructor.
Definition: specbnd.hpp:490
Specbnd(const double c)
Constructor for real scalar c
Definition: specbnd.hpp:468
Exceptions of mc::Specbnd.
Definition: specbnd.hpp:428
Specbnd(const Specbnd< T > &x)
Copy constructor.
Definition: specbnd.hpp:485
Hertz & Rohn's method.
Definition: specbnd.hpp:421
Options()
Constructor.
Definition: specbnd.hpp:415
Operation between variables with different numbers of dependents.
Definition: specbnd.hpp:435
Specbnd< T > & dep(const unsigned int i, const unsigned int n)
Set the index of a variable (and total number of variables)
Definition: specbnd.hpp:507
const T & I() const
Return function bounds.
Definition: specbnd.hpp:522
Gershgorin circle's criterion.
Definition: specbnd.hpp:420
Specbnd< T > & set(const T &B, const unsigned int i, const unsigned int n)
Set variable with range B and index i of n independent variables.
Definition: specbnd.hpp:494
Specbnd()
Default constructor (needed to declare arrays of Specbnd)
Definition: specbnd.hpp:463
C++ template class computing spectral bounds for the Hessian matrix of a factorable function on a box...
Definition: specbnd.hpp:317
const T & SI() const
Return spectral bounds for Hessian matrix.
Definition: specbnd.hpp:534
Specbnd(const T &B)
Constructor for an interval B
Definition: specbnd.hpp:473
Specbnd(const T &B, const unsigned int i, const unsigned int n)
Constructor for a variable with range B and index i of n independent variables.
Definition: specbnd.hpp:478
Failed to compute spectral bound in Specbnd::spectral_bound.
Definition: specbnd.hpp:434
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