MoochoPack_CalcLambdaIndepStd_AddedStep.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 <ostream>
00030 
00031 #include "MoochoPack_CalcLambdaIndepStd_AddedStep.hpp"
00032 #include "MoochoPack_moocho_algo_conversion.hpp"
00033 #include "IterationPack_print_algorithm_step.hpp"
00034 #include "ConstrainedOptPack_ComputeMinMult.hpp"
00035 #include "ConstrainedOptPack/src/VectorWithNorms.h"
00036 #include "AbstractLinAlgPack_SpVectorOp.hpp"
00037 #include "AbstractLinAlgPack/src/AbstractLinAlgPack_MatrixOp.hpp"
00038 #include "DenseLinAlgPack_LinAlgOpPack.hpp"
00039 #include "DenseLinAlgPack_DVectorOp.hpp"
00040 #include "DenseLinAlgPack_DVectorOut.hpp"
00041 
00042 namespace LinAlgOpPack {
00043   using AbstractLinAlgPack::Vp_StV;
00044   using AbstractLinAlgPack::Vp_StMtV;
00045 }
00046 
00047 bool MoochoPack::CalcLambdaIndepStd_AddedStep::do_step(Algorithm& _algo
00048   , poss_type step_poss, IterationPack::EDoStepType type, poss_type assoc_step_poss)
00049 {
00050 
00051   using BLAS_Cpp::no_trans;
00052   using BLAS_Cpp::trans;
00053 
00054   using DenseLinAlgPack::Vt_S;
00055   using DenseLinAlgPack::norm_inf;
00056 
00057   using AbstractLinAlgPack::Vp_StMtV;
00058   
00059   using ConstrainedOptPack::min_abs;
00060 
00061   using LinAlgOpPack::Vp_V;
00062   using LinAlgOpPack::V_StMtV;
00063   using LinAlgOpPack::Vp_MtV;
00064 
00065   NLPAlgo &algo = rsqp_algo(_algo);
00066   NLPAlgoState  &s    = algo.rsqp_state();
00067   Range1D   decomp  = s.equ_decomp();
00068   
00069   EJournalOutputLevel olevel = algo.algo_cntr().journal_output_level();
00070   std::ostream& out = algo.track().journal_out();
00071 
00072   // print step header.
00073   if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
00074     using IterationPack::print_algorithm_step;
00075     print_algorithm_step( algo, step_poss, type, assoc_step_poss, out );
00076   }
00077 
00078   // Compute: lambda(decomp) = inv(Gc(decomp)'* Y)' * ( - Y' * (Gf + nu) - U' * lambda(undecomp) )
00079   // where U = Gc(undecomp)' * Y
00080 
00081   // Must resize lambda here explicitly since we will only be updating a region of it.
00082   // If lambda(undecomp) has already been updated then lambda will have been resized
00083   // already but lambda(decomp) will not be initialized yet.
00084   if( !s.lambda().updated_k(0) ) s.lambda().set_k(0).v().resize( algo.nlp().m() );
00085   
00086   DVectorSlice lambda_decomp = s.lambda().get_k(0).v()(decomp);
00087   
00088   // lambda_decomp_tmp1 = - Y' * (Gf + nu)
00089   if( algo.nlp().has_bounds() ) {
00090     // _tmp = Gf + nu
00091     DVector _tmp = s.Gf().get_k(0)();
00092     DVectorSlice _vs_tmp = _tmp;  // only create this DVectorSlice once
00093     Vp_V( &_vs_tmp, s.nu().get_k(0)() );
00094     // lambda_decomp_tmp1 = - Y' * _tmp
00095     V_StMtV( &lambda_decomp, -1.0, s.Y().get_k(0), trans, _vs_tmp );
00096   }
00097   else {
00098     // lambda_decomp__tmp1 = - Y' * Gf
00099     V_StMtV( &lambda_decomp, -1.0, s.Y().get_k(0), trans, s.Gf().get_k(0)() );
00100   }
00101 
00102   // lambda_decomp_tmp2 = lambda_decomp_tmp1 - U' * lambda(undecomp)
00103   if( algo.nlp().r() < algo.nlp().m() ) {
00104     Range1D undecomp = s.equ_undecomp();
00105     Vp_StMtV( &lambda_decomp, -1.0, s.U().get_k(0), trans, s.lambda().get_k(0).v()(undecomp) );
00106   }
00107   // else lambda(decomp)_tmp2 = lambda(decomp)_tmp1
00108 
00109   // lambda_decomp = inv(Gc(decomp)'* Y)' * lambda_decomp_tmp2
00110   s.decomp_sys().solve_transAtY( lambda_decomp, trans, &lambda_decomp );
00111 
00112   if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
00113     out << "\nmax(|lambda_k(equ_decomp)(i)|) = " << norm_inf(lambda_decomp)
00114       << "\nmin(|lambda_k(equ_decomp)(i)|) = " << min_abs(lambda_decomp)  << std::endl;
00115   }
00116 
00117   if( static_cast<int>(olevel) >= static_cast<int>(PRINT_VECTORS) ) {
00118     out << "\nlambda_k(equ_decomp) = \n" << lambda_decomp;
00119   }
00120 
00121   return true;
00122 }
00123 
00124 void MoochoPack::CalcLambdaIndepStd_AddedStep::print_step( const Algorithm& algo
00125   , poss_type step_poss, IterationPack::EDoStepType type, poss_type assoc_step_poss
00126   , std::ostream& out, const std::string& L ) const
00127 {
00128   out
00129     << L << "*** Calculate the Lagrange multipliers for the decomposed constraints\n"
00130     << L << "lambda_k(equ_decomp) = - inv(Gc_k(:,equ_decomp)'*Y_k)\n"
00131     << L << "                        * (Y_k'*(Gf_k + nu_k) + U_k'*lambda_k(equ_undecomp))\n";
00132 }

Generated on Tue Oct 20 12:51:46 2009 for MOOCHO (Single Doxygen Collection) by doxygen 1.4.7