MoochoPack_BFGSUpdate_Strategy.cpp

Go to the documentation of this file.
00001 // @HEADER
00002 // ***********************************************************************
00003 // 
00004 // Moocho: Multi-functional Object-Oriented arCHitecture for Optimization
00005 //                  Copyright (2003) Sandia Corporation
00006 // 
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 // 
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //  
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //  
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov) 
00025 // 
00026 // ***********************************************************************
00027 // @HEADER
00028 
00029 #include <limits>
00030 
00031 #include "MoochoPack_BFGSUpdate_Strategy.hpp"
00032 #include "MoochoPack_Exceptions.hpp"
00033 #include "AbstractLinAlgPack_TestMatrixSymSecant.hpp"
00034 #include "AbstractLinAlgPack_MatrixSymSecant.hpp"
00035 #include "AbstractLinAlgPack_MatrixSymOp.hpp"
00036 #include "AbstractLinAlgPack_MatrixOpOut.hpp"
00037 #include "AbstractLinAlgPack_VectorSpace.hpp"
00038 #include "AbstractLinAlgPack_VectorMutable.hpp"
00039 #include "AbstractLinAlgPack_VectorOut.hpp"
00040 #include "AbstractLinAlgPack_VectorStdOps.hpp"
00041 #include "AbstractLinAlgPack_LinAlgOpPack.hpp"
00042 #include "Teuchos_dyn_cast.hpp"
00043 #include "Teuchos_TestForException.hpp"
00044 
00045 namespace MoochoPack {
00046 
00047 BFGSUpdate_Strategy::BFGSUpdate_Strategy(
00048   bool               rescale_init_identity
00049   ,bool              use_dampening
00050   ,ESecantTesting    secant_testing
00051   ,value_type        secant_warning_tol
00052   ,value_type        secant_error_tol
00053   )
00054   :rescale_init_identity_(rescale_init_identity)
00055   ,use_dampening_(use_dampening)
00056   ,secant_testing_(secant_testing)
00057   ,secant_warning_tol_(secant_warning_tol)
00058   ,secant_error_tol_(secant_error_tol)
00059 {}
00060 
00061 void BFGSUpdate_Strategy::perform_update(
00062   VectorMutable            *s_bfgs
00063   ,VectorMutable           *y_bfgs
00064   ,bool                    first_update
00065   ,std::ostream            &out
00066   ,EJournalOutputLevel     olevel
00067   ,bool                    check_results
00068   ,MatrixSymOp             *B
00069   ,QuasiNewtonStats        *quasi_newton_stats 
00070   )
00071 {
00072   namespace rcp = MemMngPack;
00073   using Teuchos::dyn_cast;
00074   using AbstractLinAlgPack::dot;
00075   using AbstractLinAlgPack::Vt_S;
00076   using AbstractLinAlgPack::Vp_StV;
00077   using LinAlgOpPack::Vp_V;
00078   using LinAlgOpPack::V_StV;
00079   using LinAlgOpPack::V_VpV;
00080   using LinAlgOpPack::V_VmV;
00081   using LinAlgOpPack::V_MtV;
00082 
00083   const value_type
00084     NOT_CALCULATED = std::numeric_limits<value_type>::max();
00085   value_type
00086     sTy = NOT_CALCULATED,
00087     yTy = NOT_CALCULATED;
00088   bool used_dampening = false;
00089 
00090   // Make sure the update is defined!
00091   if( sTy == NOT_CALCULATED )
00092     sTy = dot( *s_bfgs, *y_bfgs );
00093   if( sTy == 0 ) {
00094     if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) {
00095       out << "\n(y'*s) == 0, skipping the update ...\n";
00096     }
00097     quasi_newton_stats->set_updated_stats( QuasiNewtonStats::INDEF_SKIPED );
00098     return;
00099   }
00100 
00101   MatrixSymSecant
00102       &B_updatable = dyn_cast<MatrixSymSecant>(*B);
00103 
00104   // /////////////////////////////////////////////////////////////
00105   // Consider rescaling the initial identity hessian before
00106   // the update.
00107   // 
00108   // This was taken from Nocedal & Wright, 1999, p. 200 
00109   // 
00110   // Bo = (y'*y)/(y'*s) * I
00111   // 
00112   if( first_update && rescale_init_identity() ) {
00113     if( sTy == NOT_CALCULATED )
00114       sTy = dot( *s_bfgs, *y_bfgs );
00115     if( yTy == NOT_CALCULATED )
00116       yTy = dot( *y_bfgs, *y_bfgs );
00117     const value_type
00118       Iscale = yTy/sTy;
00119     const value_type
00120       Iscale_too_small = 1e-5;  // ToDo: Make this adjustable
00121     if( Iscale >= Iscale_too_small ) {  
00122       if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) {
00123         out << "\nRescaling the initial identity matrix before the update as:\n"
00124           << "Iscale = (y'*y)/(y'*s) = ("<<yTy<<")/("<<sTy<<") = "<<Iscale<<"\n"
00125           << "B =  Iscale * eye(n-r) ...\n";
00126       }
00127       B_updatable.init_identity( y_bfgs->space(), Iscale );
00128       if( (int)olevel >= (int)PRINT_ITERATION_QUANTITIES ) {
00129         out << "\nB after rescaling = \n" << *B;
00130       }
00131     }
00132     else {
00133       if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) {
00134         out << "\nSkipping rescaling of the initial identity matrix since:\n"
00135           << "Iscale = (y'*y)/(y'*s) = ("<<yTy<<")/("<<sTy<<") = "<<Iscale
00136           << " < Iscale_too_small = " << Iscale_too_small << std::endl;
00137       }
00138     }
00139   }
00140 
00141   // ////////////////////////////////////////////////////
00142   // Modify the s_bfgs and y_bfgs vectors for dampening?
00143   VectorSpace::vec_mut_ptr_t
00144     Bs = Teuchos::null;
00145   if( use_dampening() ) {
00146     if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS )
00147     {
00148       out
00149         << "\nConsidering the dampened update ...\n";
00150     }
00151     // Bs = Bm1 * s_bfgs
00152     Bs = y_bfgs->space().create_member();
00153     V_MtV( Bs.get(), *B, BLAS_Cpp::no_trans, *s_bfgs );
00154     // sTy = s_bfgs' * y_bfgs
00155     if( sTy == NOT_CALCULATED )
00156       sTy = dot( *s_bfgs, *y_bfgs );
00157     // sTBs = s_bfgs' * Bs
00158     const value_type sTBs = dot( *s_bfgs, *Bs );
00159     // Select dampening parameter theta
00160     const value_type
00161       theta = ( sTy >= 0.2 * sTBs )
00162       ? 1.0
00163       : (0.8 * sTBs ) / ( sTBs - sTy );
00164     if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS )
00165     {
00166       out
00167         << "\ns_bfgs'*y_bfgs = " << sTy
00168         << ( theta == 1.0 ? " >= " : " < " )
00169         << " s_bfgs' * B * s_bfgs = " << sTBs << std::endl;
00170     }
00171     if( theta == 1.0 ) {
00172       if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS )
00173       {
00174         out
00175           << "Perform the undamped update ...\n";
00176       }
00177     }
00178     else {
00179       // y_bfgs = theta*y_bfgs + (1-theta)*B*s_bfgs
00180       Vt_S( y_bfgs, theta );
00181       Vp_StV( y_bfgs, (1.0-theta), *Bs );
00182 
00183       if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS )
00184       {
00185         out
00186           << "Dampen the update ...\n"
00187           << "theta = " << theta << std::endl
00188           << "y_bfgs = theta*y_bfgs + (1-theta)*B*s_bfgs ...\n"
00189           << "||y_bfgs||inf = " << y_bfgs->norm_inf() << std::endl;
00190       }
00191 
00192       if( (int)olevel >= (int)PRINT_VECTORS )
00193       {
00194         out
00195           << "y_bfgs =\n" << *y_bfgs;
00196       }
00197 
00198       used_dampening = true;
00199     }
00200   }
00201 
00202   // Perform the update if it is defined (s_bfgs' * y_bfgs > 0.0)
00203         
00204   VectorSpace::vec_mut_ptr_t
00205     s_bfgs_save = Teuchos::null,
00206     y_bfgs_save = Teuchos::null;
00207   if( check_results ) {
00208     // Save s and y since they may be overwritten in the update.
00209     s_bfgs_save = s_bfgs->clone();
00210     y_bfgs_save = y_bfgs->clone();
00211   }
00212   try {
00213     B_updatable.secant_update(
00214       s_bfgs
00215       ,y_bfgs
00216       ,Bs.get()
00217       );
00218   }
00219   catch( const MatrixSymSecant::UpdateFailedException& excpt ) {
00220     if( (int)olevel >= (int)PRINT_BASIC_ALGORITHM_INFO ) {
00221       out << "\nThe factorization of B failed in a major way!  Rethrow!\n";
00222     }
00223     throw;
00224   }         
00225   catch( const MatrixSymSecant::UpdateSkippedException& excpt ) {
00226     if( (int)olevel >= (int)PRINT_BASIC_ALGORITHM_INFO ) {
00227       out << excpt.what() << std::endl
00228         << "\nSkipping BFGS update.  B = B ...\n";
00229     }
00230     quasi_newton_stats->set_updated_stats(
00231       QuasiNewtonStats::INDEF_SKIPED );
00232     return;
00233   }         
00234     
00235   if( (int)olevel >= (int)PRINT_ITERATION_QUANTITIES ) {
00236     out << "\nB after the BFGS update = \n" << *B;
00237   }
00238   
00239   if( ( check_results && secant_testing() == SECANT_TEST_DEFAULT )
00240     || secant_testing() == SECANT_TEST_ALWAYS )
00241   {
00242     const bool result =
00243       AbstractLinAlgPack::TestMatrixSymSecant(
00244         *B, *s_bfgs_save, *y_bfgs_save
00245         , secant_warning_tol(), secant_error_tol()
00246         , (int)olevel >= (int)PRINT_VECTORS
00247         , (int)olevel >  (int)PRINT_NOTHING ? &out : NULL
00248         , (int)olevel >= (int)PRINT_ALGORITHM_STEPS
00249         );
00250     if( !result ) {
00251       const char
00252         msg[] = "Error, the secant property for the BFGS update failed\n"
00253         "Stopping the algorithm ...\n";
00254       out << msg;
00255       TEST_FOR_EXCEPTION(
00256         true, TestFailed
00257         ," BFGSUpdate_Strategy::perform_update(...) : " << msg );
00258     }
00259   }
00260   
00261   quasi_newton_stats->set_updated_stats(
00262     used_dampening
00263     ? QuasiNewtonStats::DAMPENED_UPDATED
00264     : QuasiNewtonStats::UPDATED );
00265 }
00266 
00267 void BFGSUpdate_Strategy::print_step( std::ostream& out, const std::string& L ) const
00268 {
00269   out
00270     << L << "if use_dampening == true then\n"
00271     << L << "    if s'*y >= 0.2 * s'*B*s then\n"
00272     << L << "        theta = 1.0\n"
00273     << L << "    else\n"
00274     << L << "        theta = 0.8*s'*B*s / (s'*B*s - s'*y)\n"
00275     << L << "    end\n"
00276     << L << "    y = theta*y + (1-theta)*B*s\n"
00277     << L << "end\n"
00278     << L << "if first_update && rescale_init_identity and y'*s is sufficently positive then\n"
00279     << L << "    B = |(y'*y)/(y'*s)| * eye(size(s))\n"
00280     << L << "end\n"
00281     << L << "if s'*y is sufficently positive then\n"
00282     << L << "    *** Peform BFGS update\n"
00283     << L << "    (B, s, y ) -> B (through MatrixSymSecantUpdate interface)\n"
00284     << L << "    if ( check_results && secant_testing == SECANT_TEST_DEFAULT )\n"
00285     << L << "    or ( secant_testing == SECANT_TEST_ALWAYS ) then\n"
00286     << L << "        if B*s != y then\n"
00287     << L << "            *** The secant condition does not check out\n"
00288     << L << "            Throw TestFailed!\n"
00289     << L << "        end\n"
00290     << L << "    end\n"
00291     << L << "end\n"
00292     ;
00293 }
00294 
00295 }  // end namespace MoochoPack

Generated on Wed May 12 21:52:30 2010 for MOOCHO (Single Doxygen Collection) by  doxygen 1.4.7