MoochoPack_CheckBasisFromCPy_Step.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 <typeinfo>
00030 
00031 #include "MoochoPack/src/std/CheckBasisFromCPy_Step.h"
00032 #include "MoochoPack_moocho_algo_conversion.hpp"
00033 #include "IterationPack_print_algorithm_step.hpp"
00034 #include "ConstrainedOptPack_DecompositionSystemVarReduct.hpp"
00035 #include "ConstrainedOptPack/src/VectorWithNorms.h"
00036 #include "AbstractLinAlgPack/src/AbstractLinAlgPack_MatrixOp.hpp"
00037 #include "DenseLinAlgPack_DVectorClass.hpp"
00038 #include "DenseLinAlgPack_DVectorOp.hpp"
00039 #include "DenseLinAlgPack_LinAlgOpPack.hpp"
00040 #include "Midynamic_cast_verbose.h"
00041 #include "MiWorkspacePack.h"
00042 
00043 namespace LinAlgOpPack {
00044   using AbstractLinAlgPack::Vp_StMtV;
00045 }
00046 
00047 namespace MoochoPack {
00048 
00049 CheckBasisFromCPy_Step::CheckBasisFromCPy_Step(
00050   const new_basis_strategy_ptr_t& new_basis_strategy
00051   ,value_type max_basis_cond_change_frac
00052   )
00053   : new_basis_strategy_(new_basis_strategy)
00054   , max_basis_cond_change_frac_( max_basis_cond_change_frac )
00055 {
00056   reset();
00057 }
00058 
00059 void CheckBasisFromCPy_Step::reset() {
00060   beta_min_ = std::numeric_limits<value_type>::max();
00061 }
00062 
00063 // Overridden
00064 
00065 bool CheckBasisFromCPy_Step::do_step( Algorithm& _algo, poss_type step_poss
00066   , IterationPack::EDoStepType type, poss_type assoc_step_poss )
00067 {
00068   using DenseLinAlgPack::norm_inf;
00069   using LinAlgOpPack::V_MtV;
00070   using LinAlgOpPack::Vp_V;
00071   using Teuchos::dyn_cast;
00072   using Teuchos::Workspace;
00073   Teuchos::WorkspaceStore* wss = Teuchos::get_default_workspace_store().get();
00074 
00075   NLPAlgo &algo = rsqp_algo(_algo);
00076   NLPAlgoState  &s    = algo.rsqp_state();
00077   Range1D   decomp  = s.con_decomp();
00078 
00079   EJournalOutputLevel olevel = algo.algo_cntr().journal_output_level();
00080   std::ostream& out = algo.track().journal_out();
00081 
00082   // print step header.
00083   if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
00084     using IterationPack::print_algorithm_step;
00085     print_algorithm_step( algo, step_poss, type, assoc_step_poss, out );
00086   }
00087 
00088   bool select_new_basis = false;
00089 
00090   // Compute: resid = (Gc(decomp)'*Y) * py + c(decomp)
00091   DVectorSlice                  py_k = s.py().get_k(0)();
00092   DVectorSlice                  c_k  = s.c().get_k(0)();
00093   DVectorSlice                  c_decomp_k  = c_k(decomp);
00094   Workspace<value_type>   resid_ws(wss,py_k.size());
00095   DVectorSlice                  resid(&resid_ws[0],resid_ws.size());
00096   {
00097 #ifdef _WINDOWS
00098     DecompositionSystemVarReduct
00099       &decomp_sys = dynamic_cast<DecompositionSystemVarReduct&>(s.decomp_sys());
00100 #else
00101     DecompositionSystemVarReduct
00102       &decomp_sys = dyn_cast<DecompositionSystemVarReduct>(s.decomp_sys());
00103 #endif
00104     V_MtV( &resid, decomp_sys.C(), BLAS_Cpp::no_trans, py_k );
00105     Vp_V( &resid, c_decomp_k );
00106   }
00107 
00108   const value_type
00109     small_num    = std::numeric_limits<value_type>::min(),
00110     nrm_resid    = norm_inf(resid),
00111     nrm_c_decomp = norm_inf(c_decomp_k),
00112     beta         = nrm_resid / (nrm_c_decomp+small_num);
00113 
00114   if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) {
00115     out << "\nbeta = ||(Gc(decomp)'*Y)*py_k + c_k(decomp)||inf / (||c_k(decomp)||inf + small_number)"
00116       << "\n     = "<<nrm_resid<<" / ("<<nrm_c_decomp<<" + "<<small_num<<")"
00117       << "\n     = " << beta << std::endl;
00118   }
00119   
00120   if( beta != 0.0 ) {
00121     if( beta < beta_min_ ) {
00122       beta_min_ = beta;
00123     }
00124     else {
00125       if( beta / beta_min_ > max_basis_cond_change_frac() ) {
00126         if( (int)olevel >= (int)PRINT_BASIC_ALGORITHM_INFO ) {
00127           out
00128             << "\nbeta_change = ( ||R*py+c||/||c|| = " << beta
00129             << " ) / ( min ||R*py+c||/||c|| = " << beta_min_ << " )\n"
00130             << "              = " << (beta/beta_min_) << " > max_basis_cond_change_frac = "
00131             << max_basis_cond_change_frac()
00132             << "\n\nSelecting a new basis"
00133             << " (k = " << algo.state().k() << ") ...\n";
00134         }
00135         select_new_basis = true;
00136       }
00137     }
00138     if(select_new_basis) {
00139       // reset memory
00140       // ToDo: You really need to check to see if someone else
00141       //     changed the basis also.
00142       beta_min_ = std::numeric_limits<value_type>::max();
00143       return new_basis_strategy().transistion_basis(algo,step_poss,type,assoc_step_poss);
00144     }
00145   }
00146   return true;
00147 }
00148 
00149 void CheckBasisFromCPy_Step::print_step( const Algorithm& algo, poss_type step_poss
00150   , IterationPack::EDoStepType type, poss_type assoc_step_poss
00151   , std::ostream& out, const std::string& L ) const
00152 {
00153   out
00154     << L << "*** Try to detect when the basis is becomming illconditioned\n"
00155     << L << "default: beta_min = inf\n"
00156     << L << "         max_basis_cond_change_frac = " << max_basis_cond_change_frac() << std::endl
00157     << L << "beta = norm_inf((Gc(decomp)'*Y)*py_k + c_k(decomp)) / (norm_inf(c_k(decomp))+small_number)\n"
00158     << L << "select_new_basis = false\n"
00159     << L << "if beta < beta_min then\n"
00160     << L << "    beta_min = beta\n"
00161     << L << "else\n"
00162     << L << "    if beta/ beta_min > max_basis_cond_change_frac then\n"
00163     << L << "        select_new_basis = true\n"
00164     << L << "    end\n"
00165     << L << "end\n"
00166     << L << "if select_new_basis == true then\n"
00167     << L << "    new basis selection : " << typeid(new_basis_strategy()).name() << std::endl;
00168   new_basis_strategy().print_method( rsqp_algo(algo),step_poss,type,assoc_step_poss,out
00169     ,L+"        " );
00170   out
00171     << L << "    end new basis selection\n"
00172     << L << "end\n";
00173 }
00174 
00175 } // end namespace MoochoPack 

Generated on Thu Sep 18 12:35:16 2008 for MOOCHO (Single Doxygen Collection) by doxygen 1.3.9.1