MC++
specbnd.hpp
1 // Copyright (C) 2009-2013 Benoit Chachuat, Imperial College London.
2 // All Rights Reserved.
3 
293 #ifndef MC__SPECBND_H
294 #define MC__SPECBND_H
295 
296 #include <iostream>
297 #include <cmath>
298 
299 #include "mclapack.hpp"
300 #include "mcop.hpp"
301 #include "fadiff.h"
302 #include "badiff.h"
303 
304 #undef MC__SPECBND_DEBUG_SPECTRUM
305 #undef MC__SPECBND_DEBUG_HESSBND
306 
307 namespace mc
308 {
316 template <typename T>
317 class Specbnd
319 {
320  template <typename U> friend class Specbnd;
321 
322  template <class U> friend Specbnd<U> operator+(const Specbnd<U> &x );
323  template <class U> friend Specbnd<U> operator+(const Specbnd<U> &x, const Specbnd<U> &y);
324  template <class U> friend Specbnd<U> operator+(const Specbnd<U> &y, const double c);
325  template <class U> friend Specbnd<U> operator+(const double c, const Specbnd<U> &y);
326  template <class U> friend Specbnd<U> operator+(const int c, const Specbnd<U> &y);
327  template <class U> friend Specbnd<U> operator-(const Specbnd<U> &x );
328  template <class U> friend Specbnd<U> operator-(const Specbnd<U> &x, const Specbnd<U> &y);
329  template <class U> friend Specbnd<U> operator-(const Specbnd<U> &y, const double c);
330  template <class U> friend Specbnd<U> operator-(const double c, const Specbnd<U> &y);
331  template <class U> friend Specbnd<U> operator*(const Specbnd<U> &x, const Specbnd<U> &y);
332  template <class U> friend Specbnd<U> operator*(const double c, const Specbnd<U> &y);
333  template <class U> friend Specbnd<U> operator*(const Specbnd<U> &y, const double c);
334  template <class U> friend Specbnd<U> pow(const Specbnd<U> &x, const int m);
335  template <class U> friend Specbnd<U> pow(const Specbnd<U> &x, const double c);
336  template <class U> friend Specbnd<U> pow(const Specbnd<U> &x, const Specbnd<U> &y);
337  template <class U> friend Specbnd<U> pow(const double c, const Specbnd<U> &y);
338  template <class U> friend Specbnd<U> inv(const Specbnd<U> &y);
339  template <class U> friend Specbnd<U> operator/(const Specbnd<U> &x, const Specbnd<U> &y);
340  template <class U> friend Specbnd<U> operator/(const Specbnd<U> &y, const double c);
341  template <class U> friend Specbnd<U> operator/(const double c, const Specbnd<U> &y);
342  template <class U> friend Specbnd<U> sqr(const Specbnd<U> &y);
343  template <class U> friend Specbnd<U> sqrt(const Specbnd<U> &y);
344  template <class U> friend Specbnd<U> exp(const Specbnd<U> &y);
345  template <class U> friend Specbnd<U> log(const Specbnd<U> &y);
346  template <class U> friend Specbnd<U> xlog(const Specbnd<U> &y);
347  template <class U> friend Specbnd<U> cos(const Specbnd<U> &y);
348  template <class U> friend Specbnd<U> sin(const Specbnd<U> &y);
349  template <class U> friend Specbnd<U> acos(const Specbnd<U> &y);
350  template <class U> friend Specbnd<U> asin(const Specbnd<U> &y);
351  template <class U> friend Specbnd<U> tan(const Specbnd<U> &y);
352  template <class U> friend Specbnd<U> atan(const Specbnd<U> &y);
353  template <class U> friend Specbnd<U> erf(const Specbnd<U> &y);
354  template <class U> friend Specbnd<U> erfc(const Specbnd<U> &y);
355  template <class U> friend std::ostream& operator<<(std::ostream&, const Specbnd<U>&);
356 
357 private:
359  unsigned int _n;
361  fadbad::F<T> _FI;
363  T _spec;
364 
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 );
369 
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 );
382 
383 public:
384  // other operator overloadings
385  Specbnd<T>& operator+=
386  ( const Specbnd<T>& );
387  Specbnd<T>& operator+=
388  ( const double );
389  Specbnd<T>& operator-=
390  ( const Specbnd<T>& );
391  Specbnd<T>& operator-=
392  ( const double );
393  Specbnd<T>& operator*=
394  ( const Specbnd<T>& );
395  Specbnd<T>& operator*=
396  ( const double );
397  Specbnd<T>& operator/=
398  ( const Specbnd<T>& );
399  Specbnd<T>& operator/=
400  ( const double );
401  Specbnd<T> & operator=
402  ( const Specbnd<T> &x );
403  Specbnd<T> & operator=
404  ( const double c );
405  Specbnd<T> & operator=
406  ( const T &c );
407 
411  static struct Options
413  {
417  {}
422  };
425  } options;
426 
429  {
430  public:
432  enum TYPE{
433  SPECTR=1,
435  SIZE=-1,
436  UNDEF=-33
437  };
439  Exceptions( TYPE ierr ) : _ierr( ierr ){}
441  int ierr(){ return _ierr; }
443  std::string what(){
444  switch( _ierr ){
445  case SPECTR:
446  return "mc::Specbnd\t Computation of Hessian spectrum failed";
447  case HESSBND:
448  return "mc::Specbnd\t Computation of Hessian spectral bound failed";
449  case SIZE:
450  return "mc::Specbnd\t Operation between variables with different numbers of dependents";
451  case UNDEF:
452  return "mc::Specbnd\t Feature not yet implemented in mc::Specbnd class";
453  default:
454  return "mc::Specbnd\t Undocumented error";
455  }
456  }
457 
458  private:
459  TYPE _ierr;
460  };
461 
464  _n(0), _FI(0), _spec(0.)
465  {}
466 
468  Specbnd( const double c ):
469  _n(0), _FI(c), _spec(0.)
470  {}
471 
473  Specbnd( const T &B ):
474  _n(0), _FI(B), _spec(0.)
475  {}
476 
478  Specbnd( const T &B, const unsigned int i, const unsigned int n ):
479  _n(n), _FI(B), _spec(0.)
480  {
481  _FI.diff(i,_n);
482  }
483 
485  Specbnd(const Specbnd<T> &x):
486  _n(x._n), _FI(x._FI), _spec(x._spec)
487  {}
488 
491  {}
492 
494  Specbnd<T>& set( const T &B, const unsigned int i, const unsigned int n )
495  {
496  _n = n;
497  _FI = B;
498  _FI.diff(i,_n);
499  _spec = 0.;
500  return *this;
501  }
502 
504  Specbnd<T>& set( const fadbad::F<T>&FB, const T&SB );
505 
507  Specbnd<T>& dep( const unsigned int i, const unsigned int n )
508  {
509  _n = n;
510  _FI.diff( i, n );
511  _spec = 0.;
512  return *this;
513  }
514 
516  unsigned int dep() const
517  {
518  return _n;
519  }
520 
522  const T& I() const
523  {
524  return _FI.val();
525  }
526 
528  const fadbad::F<T>& FI() const
529  {
530  return _FI;
531  }
532 
534  const T& SI() const
535  {
536  return _spec;
537  }
538 
540  static std::pair<double,double> spectrum
541  ( const fadbad::B< fadbad::F< double > >* D2X );
542 
544  static std::pair<double,double> spectrum
545  ( const fadbad::F< fadbad::F< double > >& D2X );
546 
548  static std::pair<double,double> spectral_bound
549  ( const fadbad::B< fadbad::F< T > >* D2X );
550 
552  static std::pair<double,double> spectral_bound
553  ( const fadbad::F< fadbad::F< T > >& D2X );
555 };
556 
558 
559 template <typename T> typename Specbnd<T>::Options Specbnd<T>::options;
560 
561 template <class T> inline T
562 Specbnd<T>::_LambdaS
563 ( const fadbad::F<T> &a, unsigned int n )
564 {
565  if( !n ) return 0.;
566  if( n == 1 ) return a.size()? Op<T>::sqr( a[0] ): 0.;
567  double upbnd=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;
571 }
572 
573 template <class T> inline T
574 Specbnd<T>::_LambdaT
575 ( const fadbad::F<T> &a, const fadbad::F<T> &b, unsigned int n )
576 {
577  if( !n ) return 0.;
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] ) );
583  }
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];
587  return lamb;
588 }
589 
590 template <typename T> inline std::pair<double,double>
592 ( const fadbad::B< fadbad::F< double > >* D2X )
593 {
594  if( !D2X ) throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::SPECTR );
595  const unsigned int N = D2X->val().size();
596  if( !N ) throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::SPECTR );
597  double*H = new double[N*N];
598 
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 );
605 #endif
606 
607  double*D = mc::dsyev_wrapper( N, H );
608  delete[] H;
609  if( !D ) throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::SPECTR );
610  std::pair<double,double> spbnd = std::make_pair( min(N,D), max(N,D) );
611  delete[] D;
612  return spbnd;
613 }
614 
615 template <typename T> inline std::pair<double,double>
617 ( const fadbad::F< fadbad::F< double > >& D2F )
618 {
619  const unsigned int N = D2F.size();
620  if( !N ) throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::SPECTR );
621  double*H = new double[N*N];
622 
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 );
629 #endif
630 
631  double*D = mc::dsyev_wrapper( N, H );
632  delete[] H;
633  if( !D ) throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::SPECTR );
634  std::pair<double,double> spbnd = std::make_pair( min(N,D), max(N,D) );
635  delete[] D;
636  return spbnd;
637 }
638 
639 template <typename T> inline std::pair<double,double>
641 ( const fadbad::B< fadbad::F< T > >* D2X )
642 {
643  if( !D2X ) throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::HESSBND );
644  switch( options.HESSBND ){
645  case Options::GERSHGORIN:
646  return _gershgorin_bound( D2X );
647  case Options::HERTZROHN: default:
648  return _hertzrohn_bound( D2X );
649  }
650 }
651 
652 template <typename T> inline std::pair<double,double>
654 ( const fadbad::B< fadbad::F< T > >* D2X )
655 {
656  const unsigned int N = D2X->val().size();
657  if( !N ) throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::HESSBND );
658 
659  std::pair<double,double> spbnd;
660  for( unsigned int i=0; i<N; i++ ){
661  double ri = 0.;
662  for( unsigned int j=0; j<N; j++ ){
663  if( j == i ) continue;
664  T D2Fij;
665  if( !inter( D2Fij, D2X[i].deriv(0).deriv(j), D2X[i].deriv(0).deriv(j) ) )
666  D2Fij = D2X[i].deriv(0).deriv(j);
667  ri += Op<T>::abs( D2Fij );
668  }
669 
670  if( !i ){
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;
673  }
674  else{
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 );
677  }
678  }
679  return spbnd;
680 }
681 
682 template <typename T> inline std::pair<double,double>
683 Specbnd<T>::_hertzrohn_bound
684 ( const fadbad::B< fadbad::F< T > >* D2X )
685 {
686  const unsigned int N = D2X->val().size();
687  if( !N ) throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::HESSBND );
688 
689  // Form matrix S^{(N)} recursively
690  unsigned int col_S = 2;
691  short *S = new short[col_S];
692  S[0] = 1; S[1] = -1;
693 #ifdef MC__SPECBND_DEBUG_HESSBND
694  mc::display( 1, col_S, S, 1, "\nMatrix S1", std::cout );
695 #endif
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];
705  }
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;
708  }
709  delete[] Sprev;
710 #ifdef MC__SPECBND_DEBUG_HESSBND
711  mc::display( k+1, 2*col_S, S, k+1, "\nMatrix Sk", std::cout );
712 #endif
713  }
714 
715  // Compute lower and upper bound on spectral radius
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++ ){
719 
720  for( unsigned int j=0; j<N; j++ )
721  for( unsigned int i=j; i<N; i++ ){
722  if( i == j ){
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) );
725  continue;
726  }
727  T D2Fij;
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 );
732  }
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 );
736 #endif
737 
738  double*DLk = mc::dsyev_wrapper( N, Lk );
739  if( !DLk ){
740  delete[] Lk; delete[] Uk; delete[] S;
741  throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::HESSBND );
742  }
743 #ifdef MC__SPECBND_DEBUG_HESSBND
744  mc::display( 1, N, DLk, 1, "\nMatrix DLk", std::cout );
745 #endif
746  spbnd.first = k? std::min( spbnd.first, min(N,DLk) ): min(N,DLk);
747  delete[] DLk;
748 
749  double*DUk = mc::dsyev_wrapper( N, Uk );
750  if( !DUk ){
751  delete[] Lk; delete[] Uk; delete[] S;
752  throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::HESSBND );
753  }
754 #ifdef MC__SPECBND_DEBUG_HESSBND
755  mc::display( 1, N, DUk, 1, "\nMatrix DUk", std::cout );
756 #endif
757  spbnd.second = k? std::max( spbnd.second, max(N,DUk) ): max(N,DUk);
758  delete[] DUk;
759 
760  }
761  delete[] Lk; delete[] Uk; delete[] S;
762  return spbnd;
763 }
764 
765 template <typename T> inline std::pair<double,double>
767 ( const fadbad::F< fadbad::F< T > >& D2F )
768 {
769  switch( options.HESSBND ){
770  case Options::GERSHGORIN:
771  return _gershgorin_bound( D2F );
772  case Options::HERTZROHN: default:
773  return _hertzrohn_bound( D2F );
774  }
775 }
776 
777 template <typename T> inline std::pair<double,double>
779 ( const fadbad::F< fadbad::F< T > >& D2F )
780 {
781  const unsigned int N = D2F.size();
782  if( !N ) throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::HESSBND );
783 
784  std::pair<double,double> spbnd;
785  for( unsigned int i=0; i<N; i++ ){
786  double ri = 0.;
787  for( unsigned int j=0; j<N; j++ ){
788  if( j == i ) continue;
789  T D2Fij;
790  if( !inter( D2Fij, D2F.deriv(i).deriv(j), D2F.deriv(i).deriv(j) ) )
791  D2Fij = D2F.deriv(i).deriv(j);
792  ri += Op<T>::abs( D2Fij );
793  }
794 
795  if( !i ){
796  spbnd.first = Op<T>::l(D2F.deriv(i).deriv(i)) - ri;
797  spbnd.second = Op<T>::u(D2F.deriv(i).deriv(i)) + ri;
798  }
799  else{
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 );
802  }
803  }
804  return spbnd;
805 }
806 
807 template <typename T> inline std::pair<double,double>
808 Specbnd<T>::_hertzrohn_bound
809 ( const fadbad::F< fadbad::F< T > >& D2F )
810 {
811  const unsigned int N = D2F.size();
812  if( !N ) throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::HESSBND );
813 
814  // Form matrix S^{(N)} recursively
815  unsigned int col_S = 2;
816  short *S = new short[col_S];
817  S[0] = 1; S[1] = -1;
818 #ifdef MC__SPECBND_DEBUG_HESSBND
819  mc::display( 1, col_S, S, 1, "\nMatrix S1", std::cout );
820 #endif
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];
830  }
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;
833  }
834  delete[] Sprev;
835 #ifdef MC__SPECBND_DEBUG_HESSBND
836  mc::display( k+1, 2*col_S, S, k+1, "\nMatrix Sk", std::cout );
837 #endif
838  }
839 
840 #ifdef MC__SPECBND_DEBUG_HESSBND
841  T *H = new T[N*N];
842  for( unsigned int j=0; j<N; j++ )
843  for( unsigned int i=j; i<N; i++ ){
844  if( i == j ){
845  H[i*N+i] = D2F.deriv(i).deriv(i);
846  continue;
847  }
848  T D2Fij;
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;
852  }
853  mc::display( N, N, H, N, "\nMatrix H", std::cout );
854  delete[] H;
855 #endif
856 
857  // Compute lower and upper bound on spectral radius
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++ ){
861 
862  for( unsigned int j=0; j<N; j++ )
863  for( unsigned int i=j; i<N; i++ ){
864  if( i == j ){
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) );
867  continue;
868  }
869  T D2Fij;
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 );
874  }
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 );
878 #endif
879 
880  double*DLk = mc::dsyev_wrapper( N, Lk );
881  if( !DLk ){
882  delete[] Lk; delete[] Uk; delete[] S;
883  throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::HESSBND );
884  }
885 #ifdef MC__SPECBND_DEBUG_HESSBND
886  mc::display( 1, N, DLk, 1, "\nMatrix DLk", std::cout );
887 #endif
888  spbnd.first = k? std::min( spbnd.first, min(N,DLk) ): min(N,DLk);
889  delete[] DLk;
890 
891  double*DUk = mc::dsyev_wrapper( N, Uk );
892  if( !DUk ){
893  delete[] Lk; delete[] Uk; delete[] S;
894  throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::HESSBND );
895  }
896 #ifdef MC__SPECBND_DEBUG_HESSBND
897  mc::display( 1, N, DUk, 1, "\nMatrix DUk", std::cout );
898 #endif
899  spbnd.second = k? std::max( spbnd.second, max(N,DUk) ): max(N,DUk);
900  delete[] DUk;
901 
902  }
903  delete[] Lk; delete[] Uk; delete[] S;
904  return spbnd;
905 }
906 
907 template <class T> inline std::ostream&
908 operator<<
909 ( std::ostream &out, const Specbnd<T> &y )
910 {
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;
915  return out;
916 }
917 
918 template <class T> inline Specbnd<T>&
919 Specbnd<T>::operator=
920 ( const double c )
921 {
922  _n = 0;
923  _FI = c;
924  _spec = 0.;
925  return *this;
926 }
927 
928 template <class T> inline Specbnd<T>&
929 Specbnd<T>::operator=
930 ( const T&I )
931 {
932  _n = 0;
933  _FI = I;
934  _spec = 0.;
935  return *this;
936 }
937 
938 template <class T> inline Specbnd<T>&
940 ( const fadbad::F<T>&FB, const T&SB )
941 {
942  _n = FB.size();
943  _FI = FB;
944  _spec = SB;
945  return *this;
946 }
947 
948 template <class T> inline Specbnd<T>&
950 ( const Specbnd<T> &x )
951 {
952  _n = x._n;
953  _FI = x._FI;
954  _spec = x._spec;
955  return *this;
956 }
957 
958 template <class T> inline Specbnd<T>&
959 Specbnd<T>::operator+=
960 ( const double c )
961 {
962  _FI += c;
963  return *this;
964 }
965 
966 template <class T> inline Specbnd<T>&
967 Specbnd<T>::operator+=
968 ( const Specbnd<T> &y )
969 {
970  if( _n && y._n && _n != y._n )
971  throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::SIZE );
972  _FI += y._FI;
973  _n = _FI.size();
974  _spec += y.spec;
975  return *this;
976 }
977 
978 template <class T> inline Specbnd<T>
979 operator+
980 ( const Specbnd<T> &y )
981 {
982  return y;
983 }
984 
985 template <class T> Specbnd<T>
986 operator+
987 ( const Specbnd<T> &x, const Specbnd<T> &y )
988 {
989  if( x._n && y._n && x._n != y._n )
990  throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::SIZE );
991  Specbnd<T> z;
992  z._FI = x._FI + y._FI;
993  z._n = z._FI.size();
994  z._spec = x._spec + y._spec;
995  return z;
996 }
997 
998 template <class T> Specbnd<T>
999 operator+
1000 ( const Specbnd<T> &y, const double c )
1001 {
1002  Specbnd<T> z;
1003  z._n = y._n;
1004  z._FI = y._FI + c;
1005  z._spec = y._spec;
1006  return z;
1007 }
1008 
1009 template <class T> Specbnd<T>
1010 operator+
1011 ( const double c, const Specbnd<T> &y )
1012 {
1013  return y + c;
1014 }
1015 
1016 template <class T> inline Specbnd<T>&
1017 Specbnd<T>::operator-=
1018 ( const double c )
1019 {
1020  _FI -= c;
1021  return *this;
1022 }
1023 
1024 template <class T> inline Specbnd<T>&
1025 Specbnd<T>::operator-=
1026 ( const Specbnd<T> &y )
1027 {
1028  if( _n && y._n && _n != y._n )
1029  throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::SIZE );
1030  _FI -= y._FI;
1031  _n = _FI.size();
1032  _spec -= y.spec;
1033  return *this;
1034 }
1035 
1036 template <class T> inline Specbnd<T>
1037 operator-
1038 ( const Specbnd<T> &y )
1039 {
1040  Specbnd<T> z;
1041  z._n = y._n;
1042  z._FI = -y._FI;
1043  z._spec = -y._spec;
1044  return z;
1045 }
1046 
1047 template <class T> inline Specbnd<T>
1048 operator-
1049 ( const Specbnd<T> &x, const Specbnd<T> &y )
1050 {
1051  if( x._n && y._n && x._n != y._n )
1052  throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::SIZE );
1053  Specbnd<T> z;
1054  z._FI = x._FI - y._FI;
1055  z._n = z._FI.size();
1056  z._spec = x._spec - y._spec;
1057  return z;
1058 }
1059 
1060 template <class T> inline Specbnd<T>
1061 operator-
1062 ( const Specbnd<T> &y, const double c )
1063 {
1064  return y + (-c);
1065 }
1066 
1067 template <class T> inline Specbnd<T>
1068 operator-
1069 ( const double c, const Specbnd<T> &y )
1070 {
1071  return (-y) + c;
1072 }
1073 
1074 template <typename T> inline Specbnd<T>&
1075 Specbnd<T>::operator*=
1076 ( const double c )
1077 {
1078  Specbnd<T> z = c * (*this);
1079  *this = z;
1080  return *this;
1081 }
1082 
1083 template <typename T> inline Specbnd<T>&
1084 Specbnd<T>::operator*=
1085 ( const Specbnd<T> &x )
1086 {
1087  Specbnd<T> z = x * (*this);
1088  *this = z;
1089  return *this;
1090 }
1091 
1092 template <class T> inline Specbnd<T>
1093 operator*
1094 ( const Specbnd<T> &x, const Specbnd<T> &y )
1095 {
1096  if( x._n && y._n && x._n != y._n )
1097  throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::SIZE );
1098  Specbnd<T> z;
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 );
1103  return z;
1104 }
1105 
1106 template <class T> inline Specbnd<T>
1107 operator*
1108 ( const double c, const Specbnd<T> &y )
1109 {
1110  Specbnd<T> z;
1111  z._n = y._n;
1112  z._FI = c * y._FI;
1113  z._spec = c * y._spec;
1114  return z;
1115 }
1116 
1117 template <class T> inline Specbnd<T>
1118 operator*
1119 ( const Specbnd<T> &y, const double c )
1120 {
1121  return c * y;
1122 }
1123 
1124 template <class T> inline Specbnd<T>
1125 sqr
1126 ( const Specbnd<T> &y )
1127 {
1128  Specbnd<T> z;
1129  z._n = y._n;
1130  z._FI = fadbad::sqr( y._FI );
1131  z._spec = 2.*( y._FI.val() * y._spec + Specbnd<T>::_LambdaS( y._FI, z._n ) );
1132  return z;
1133 }
1134 
1135 template <class T> inline Specbnd<T>
1136 pow
1137 ( const Specbnd<T> &y, const int m )
1138 {
1139  if( !m ) return 1.;
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 ) );
1144  Specbnd<T> z;
1145  z._n = y._n;
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 ) );
1149  return z;
1150 }
1151 
1152 template <class T> inline Specbnd<T>
1153 inv
1154 ( const Specbnd<T> &y )
1155 {
1156  Specbnd<T> z;
1157  z._n = y._n;
1158  z._FI = 1./y._FI;
1159  z._spec = Op<T>::sqr( z._FI.val() )
1160  * ( 2. * z._FI.val() * Specbnd<T>::_LambdaS( y._FI, z._n ) - y._spec );
1161  return z;
1162 }
1163 
1164 template <typename T> inline Specbnd<T>&
1165 Specbnd<T>::operator/=
1166 ( const double c )
1167 {
1168  return *this / c;
1169 }
1170 
1171 template <typename T> inline Specbnd<T>&
1172 Specbnd<T>::operator/=
1173 ( const Specbnd<T> &x )
1174 {
1175  return *this / x;
1176 }
1177 
1178 template <class T> inline Specbnd<T>
1179 operator/
1180 ( const Specbnd<T> &x, const Specbnd<T> &y )
1181 {
1182  return x * inv(y);
1183 }
1184 
1185 template <class T> inline Specbnd<T>
1186 operator/
1187 ( const Specbnd<T> &y, const double c )
1188 {
1189  return y * (1./c);
1190 }
1191 
1192 template <class T> inline Specbnd<T>
1193 operator/
1194 ( const double c, const Specbnd<T> &y )
1195 {
1196  return c * inv(y);
1197 }
1198 
1199 template <class T> inline Specbnd<T>
1200 sqrt
1201 ( const Specbnd<T> &y )
1202 {
1203  Specbnd<T> z;
1204  z._n = y._n;
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()) );
1208  return z;
1209 }
1210 
1211 template <class T> inline Specbnd<T>
1212 exp
1213 ( const Specbnd<T> &y )
1214 {
1215  Specbnd<T> z;
1216  z._n = y._n;
1217  z._FI = fadbad::exp( y._FI );
1218  z._spec = z._FI.val() * ( y._spec + Specbnd<T>::_LambdaS( y._FI, z._n ) );
1219  return z;
1220 }
1221 
1222 template <class T> inline Specbnd<T>
1223 log
1224 ( const Specbnd<T> &y )
1225 {
1226  Specbnd<T> z;
1227  z._n = y._n;
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() );
1230  return z;
1231 }
1232 
1233 template <class T> inline Specbnd<T>
1234 xlog
1235 ( const Specbnd<T> &y )
1236 {
1237  Specbnd<T> z;
1238  z._n = y._n;
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();
1242  return z;
1243 }
1244 
1245 template <class T> inline Specbnd<T>
1246 pow
1247 ( const Specbnd<T> &x, const Specbnd<T> &y )
1248 {
1249  return exp( y * log( x ) );
1250 }
1251 
1252 template <class T> inline Specbnd<T>
1253 pow
1254 ( const double c, const Specbnd<T> &y )
1255 {
1256  return exp( y * std::log( c ) );
1257 }
1258 
1259 template <class T> inline Specbnd<T>
1260 pow
1261 ( const Specbnd<T> &x, const double c )
1262 {
1263  return exp( c * log( x ) );
1264 }
1265 
1266 template <class T> inline Specbnd<T>
1267 monomial
1268 ( const unsigned int n, const Specbnd<T>*x, const int*k )
1269 {
1270  if( n == 0 ){
1271  return 1.;
1272  }
1273  if( n == 1 ){
1274  return pow( x[0], k[0] );
1275  }
1276  return pow( x[0], k[0] ) * monomial( n-1, x+1, k+1 );
1277 }
1278 
1279 template <class T> inline Specbnd<T>
1280 fabs
1281 ( const Specbnd<T> &y )
1282 {
1283  throw typename Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::UNDEF );
1284 }
1285 
1286 template <class T> inline Specbnd<T>
1287 cos
1288 ( const Specbnd<T> &y )
1289 {
1290  Specbnd<T> z;
1291  z._n = y._n;
1292  z._FI = fadbad::cos( y._FI );
1293  z._spec = - z._FI.val() * ( y._spec + Specbnd<T>::_LambdaS( y._FI, z._n ) );
1294  return z;
1295 }
1296 
1297 template <class T> inline Specbnd<T>
1298 sin
1299 ( const Specbnd<T> &y )
1300 {
1301  Specbnd<T> z;
1302  z._n = y._n;
1303  z._FI = fadbad::sin( y._FI );
1304  z._spec = z._FI.val() * ( y._spec - Specbnd<T>::_LambdaS( y._FI, z._n ) );
1305  return z;
1306 }
1307 
1308 template <class T> inline Specbnd<T>
1309 tan
1310 ( const Specbnd<T> &y )
1311 {
1312  Specbnd<T> z;
1313  z._n = y._n;
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 ) );
1317  return z;
1318 }
1319 
1320 template <class T> inline Specbnd<T>
1321 acos
1322 ( const Specbnd<T> &y )
1323 {
1324  Specbnd<T> z;
1325  z._n = y._n;
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 ) );
1330  return z;
1331 }
1332 
1333 template <class T> inline Specbnd<T>
1334 asin
1335 ( const Specbnd<T> &y )
1336 {
1337  Specbnd<T> z;
1338  z._n = y._n;
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 ) );
1343  return z;
1344 }
1345 
1346 template <class T> inline Specbnd<T>
1347 atan
1348 ( const Specbnd<T> &y )
1349 {
1350  Specbnd<T> z;
1351  z._n = y._n;
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 ) );
1356  return z;
1357 }
1358 
1359 template <class T> inline Specbnd<T>
1360 erf
1361 ( const Specbnd<T> &y )
1362 {
1364 // Specbnd<T> z;
1365 // z._n = y._n;
1366 // z._FI = fadbad::erf( y._FI );
1367 // z._spec = 2./std::sqrt(PI)*Op<T>::exp(-Op<T>::sqr(y._FI.val())) * ( y._spec
1368 // - 2.*y._FI.val()*Specbnd<T>::_LambdaS( y._FI, z._n ) );
1369 // return z;
1370 }
1371 
1372 template <class T> inline Specbnd<T>
1373 erfc
1374 ( const Specbnd<T> &y )
1375 {
1377 // Specbnd<T> z;
1378 // z._n = y._n;
1379 // z._FI = fadbad::erfc( y._FI );
1380 // z._spec = - 2./std::sqrt(PI)*Op<T>::exp(-Op<T>::sqr(y._FI.val())) * ( y._spec
1381 // - 2.*y._FI.val()*Specbnd<T>::_LambdaS( y._FI, z._n ) );
1382 // return z;
1383 }
1384 
1385 } // namespace mc
1386 
1387 
1388 #include "mcop.hpp"
1389 
1390 namespace mc
1391 {
1392 
1394 template <> template<typename T> struct Op< mc::Specbnd<T> >
1395 {
1396  typedef mc::Specbnd<T> SB;
1397  static SB point( const double c ) { return SB(c); }
1398  static SB zeroone() { return SB( mc::Op<T>::zeroone() ); }
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()); }
1402  static double abs (const SB& x) { return mc::Op<T>::abs(x.I()); }
1403  static double mid (const SB& x) { return mc::Op<T>::mid(x.I()); }
1404  static double diam(const SB& x) { return mc::Op<T>::diam(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); }
1418  static SB erf (const SB& x) { throw typename mc::Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::UNDEF ); }
1419  static SB erfc(const SB& x) { throw typename mc::Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::UNDEF ); }
1420  static SB fstep(const SB& x) { throw typename mc::Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::UNDEF ); }
1421  static SB bstep(const SB& x) { throw typename mc::Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::UNDEF ); }
1422  static SB hull(const SB& x, const SB& y) { throw typename mc::Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::UNDEF ); }
1423  static SB min (const SB& x, const SB& y) { throw typename mc::Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::UNDEF ); }
1424  static SB max (const SB& x, const SB& y) { throw typename mc::Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::UNDEF ); }
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); }
1428  static bool inter(SB& xIy, const SB& x, const SB& y) { throw typename mc::Specbnd<T>::Exceptions( Specbnd<T>::Exceptions::UNDEF ); }
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(); }
1435 };
1436 
1437 } // namespace mc
1438 
1439 #endif
1440 
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