|
MoochoPack : Framework for Large-Scale Optimization Algorithms Version of the Day
|
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 <typeinfo> 00030 00031 #include "MoochoPack_CheckDecompositionFromRPy_Step.hpp" 00032 #include "MoochoPack_moocho_algo_conversion.hpp" 00033 #include "IterationPack_print_algorithm_step.hpp" 00034 #include "AbstractLinAlgPack_MatrixOpNonsing.hpp" 00035 #include "AbstractLinAlgPack_Vector.hpp" 00036 #include "AbstractLinAlgPack_LinAlgOpPack.hpp" 00037 00038 namespace MoochoPack { 00039 00040 CheckDecompositionFromRPy_Step::CheckDecompositionFromRPy_Step( 00041 const new_decomp_strategy_ptr_t &new_decomp_strategy 00042 ,value_type max_decomposition_cond_change_frac 00043 ) 00044 :new_decomp_strategy_(new_decomp_strategy) 00045 ,max_decomposition_cond_change_frac_(max_decomposition_cond_change_frac) 00046 { 00047 reset(); 00048 } 00049 00050 void CheckDecompositionFromRPy_Step::reset() { 00051 beta_min_ = std::numeric_limits<value_type>::max(); 00052 } 00053 00054 // Overridden 00055 00056 bool CheckDecompositionFromRPy_Step::do_step( Algorithm& _algo, poss_type step_poss 00057 , IterationPack::EDoStepType type, poss_type assoc_step_poss ) 00058 { 00059 NLPAlgo &algo = rsqp_algo(_algo); 00060 NLPAlgoState &s = algo.rsqp_state(); 00061 const Range1D equ_decomp = s.equ_decomp(); 00062 EJournalOutputLevel olevel = algo.algo_cntr().journal_output_level(); 00063 std::ostream &out = algo.track().journal_out(); 00064 00065 // print step header. 00066 if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) { 00067 using IterationPack::print_algorithm_step; 00068 print_algorithm_step( algo, step_poss, type, assoc_step_poss, out ); 00069 } 00070 00071 bool select_new_decomposition = false; 00072 00073 // Compute: resid = (Gc(decomp)'*Y) * py + c(decomp) 00074 const Vector &py_k = s.py().get_k(0); 00075 const Vector &c_k = s.c().get_k(0); 00076 Vector::vec_ptr_t c_decomp_k = c_k.sub_view(equ_decomp); 00077 VectorMutable::vec_mut_ptr_t resid = c_decomp_k->space().create_member(); 00078 00079 // resid = R*py + c(equ_decomp) 00080 LinAlgOpPack::V_MtV( resid.get(), s.R().get_k(0), BLAS_Cpp::no_trans, py_k ); 00081 LinAlgOpPack::Vp_V( resid.get(), *c_decomp_k ); 00082 00083 const value_type 00084 small_num = std::numeric_limits<value_type>::min(), 00085 epsilon = std::numeric_limits<value_type>::epsilon(), 00086 nrm_resid = resid->norm_inf(), 00087 nrm_c_decomp = c_decomp_k->norm_inf(), 00088 beta = nrm_resid / (nrm_c_decomp+small_num); 00089 00090 if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) { 00091 out << "\nbeta = ||R*py_k + c_k(decomp)||inf / (||c_k(decomp)||inf + small_number)" 00092 << "\n = "<<nrm_resid<<" / ("<<nrm_c_decomp<<" + "<<small_num<<")" 00093 << "\n = " << beta << std::endl; 00094 } 00095 00096 // Check to see if a new basis was selected or not 00097 IterQuantityAccess<index_type> 00098 &num_basis_iq = s.num_basis(); 00099 if( num_basis_iq.updated_k(0) ) { 00100 if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) 00101 out << "\nnum_basis_k was updated so the basis changed so we will skip this check\n" 00102 << " reset min ||R*py+c||/||c|| to current value + epsilon(" << epsilon << ")\n"; 00103 beta_min_ = beta + epsilon; 00104 return true; 00105 } 00106 00107 if( beta != 0.0 ) { 00108 if( beta < beta_min_ ) { 00109 beta_min_ = beta; 00110 } 00111 else { 00112 if( beta / beta_min_ > max_decomposition_cond_change_frac() ) { 00113 if( (int)olevel >= (int)PRINT_BASIC_ALGORITHM_INFO ) { 00114 out 00115 << "\nbeta_change = ( ||R*py+c||/||c|| = " << beta 00116 << " ) / ( min ||R*py+c||/||c|| = " << beta_min_ << " )\n" 00117 << " = " << (beta/beta_min_) << " > max_decomposition_cond_change_frac = " 00118 << max_decomposition_cond_change_frac() 00119 << "\n\nSelecting a new decomposition" 00120 << " (k = " << algo.state().k() << ") ...\n"; 00121 } 00122 select_new_decomposition = true; 00123 } 00124 } 00125 if(select_new_decomposition) { 00126 reset(); 00127 return new_decomp_strategy().new_decomposition(algo,step_poss,type,assoc_step_poss); 00128 } 00129 } 00130 return true; 00131 } 00132 00133 void CheckDecompositionFromRPy_Step::print_step( const Algorithm& algo, poss_type step_poss 00134 , IterationPack::EDoStepType type, poss_type assoc_step_poss 00135 , std::ostream& out, const std::string& L ) const 00136 { 00137 out 00138 << L << "*** Try to detect when the decomposition is becomming illconditioned\n" 00139 << L << "default: beta_min = inf\n" 00140 << L << " max_decomposition_cond_change_frac = " << max_decomposition_cond_change_frac() << std::endl 00141 << L << "beta = norm_inf(R*py_k + c_k(decomp)) / (norm_inf(c_k(decomp))+small_number)\n" 00142 << L << "select_new_decomposition = false\n" 00143 << L << "if num_basis_k is updated then\n" 00144 << L << " beta_min = beta\n" 00145 << L << "end\n" 00146 << L << "if beta < beta_min then\n" 00147 << L << " beta_min = beta\n" 00148 << L << "else\n" 00149 << L << " if beta/ beta_min > max_decomposition_cond_change_frac then\n" 00150 << L << " select_new_decomposition = true\n" 00151 << L << " end\n" 00152 << L << "end\n" 00153 << L << "if select_new_decomposition == true then\n" 00154 << L << " new decomposition selection : " << typeName(new_decomp_strategy()) << std::endl 00155 ; 00156 new_decomp_strategy().print_new_decomposition( 00157 rsqp_algo(algo),step_poss,type,assoc_step_poss,out, L + " " ); 00158 out 00159 << L << " end new decomposition selection\n" 00160 << L << "end\n" 00161 ; 00162 } 00163 00164 } // end namespace MoochoPack
1.7.4