MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MoochoPack_InitFinDiffReducedHessian_Step.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Moocho: Multi-functional Object-Oriented arCHitecture for Optimization
5 // Copyright (2003) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #include <math.h>
43 
44 #include <ostream>
45 
53 //#include "AbstractLinAlgPack_SpVectorClass.hpp"
54 //#include "AbstractLinAlgPack/src/max_near_feas_step.hpp"
59 #include "Teuchos_dyn_cast.hpp"
60 
61 namespace {
62 template< class T >
63 inline
64 T my_max( const T& v1, const T& v2 ) { return v1 > v2 ? v1 : v2; }
65 } // end namespace
66 
67 namespace MoochoPack {
68 
70  EInitializationMethod initialization_method
71  ,value_type max_cond
72  ,value_type min_diag
73  ,value_type step_scale
74  )
75  :initialization_method_(initialization_method)
76  ,max_cond_(max_cond)
77  ,min_diag_(min_diag)
78  ,step_scale_(step_scale)
79 {}
80 
82  Algorithm& _algo, poss_type step_poss, IterationPack::EDoStepType type
83  ,poss_type assoc_step_poss
84  )
85 {
86  using Teuchos::dyn_cast;
87  using LinAlgOpPack::Vt_S;
89  using LinAlgOpPack::V_MtV;
91 
92  NLPAlgo &algo = rsqp_algo(_algo);
93  NLPAlgoState &s = algo.rsqp_state();
94  NLPObjGrad &nlp = dyn_cast<NLPObjGrad>(algo.nlp());
95 
96  EJournalOutputLevel olevel = algo.algo_cntr().journal_output_level();
97  std::ostream& out = algo.track().journal_out();
98 
99  // print step header.
100  if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
102  print_algorithm_step( algo, step_poss, type, assoc_step_poss, out );
103  }
104 
105  // Get the iteration quantity container objects
106  IterQuantityAccess<index_type>
107  &num_basis_iq = s.num_basis();
108 
109  const bool new_basis = num_basis_iq.updated_k(0);
110  const int k_last_offset = s.rHL().last_updated();
111 
112  // If the basis has changed or there is no previous matrix to use
113  // then reinitialize.
114 
115  if( new_basis || k_last_offset == IterQuantity::NONE_UPDATED ) {
116 
117  // Compute a finite difference along the null space of the
118  // constraints
119 
120  IterQuantityAccess<VectorMutable>
121  &x_iq = s.x(),
122  &rGf_iq = s.rGf();
123  IterQuantityAccess<MatrixOp>
124  &Z_iq = s.Z();
125  IterQuantityAccess<MatrixSymOp>
126  &rHL_iq = s.rHL();
127 
128  if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) {
129  out << "\nReinitializing the reduced Hessain using a finite difference\n";
130  }
131 
132  MatrixSymInitDiag &rHL_diag = dyn_cast<MatrixSymInitDiag>(rHL_iq.set_k(0));
133  const MatrixOp &Z_k = Z_iq.get_k(0);
134  const Vector &x_k = x_iq.get_k(0);
135  const Vector &rGf_k = rGf_iq.get_k(0);
136 
137  // one vector
138  VectorSpace::vec_mut_ptr_t e = Z_k.space_rows().create_member(1.0);
139 
140  // Ze
141  VectorSpace::vec_mut_ptr_t Ze = x_k.space().create_member();
142  V_MtV( Ze.get(), Z_k, BLAS_Cpp::no_trans, *e );
143 
144  // This does not have to be an accurate finite difference so lets just
145  // take step_scale/||Ze|| as the step size unless it is out of bounds
146  // If we assume that our variables are scaled near
147  // one (step_scale == 1?) then this will give us an appreciable step. Beside we
148  // should not be near the solution so the reduced gradient should not
149  // be near zero.
150 
151  const value_type nrm_Ze = Ze->norm_inf();
152  value_type u = step_scale() / nrm_Ze;
153 
154  if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) {
155  out << "\n||Ze||inf = " << nrm_Ze << std::endl;
156  }
157 
158  if( (int)olevel >= (int)PRINT_VECTORS ) {
159  out << "\nZe =\n" << *Ze;
160  }
161 
162  if( algo.nlp().num_bounded_x() ) {
163 
164  // Find the maximum step u
165  // we can take along x_fd = x_k + u*Ze
166  // that don't violate variable bounds by too much.
167  // If a positive step can't be found then this may be a negative step.
168 
169  std::pair<value_type,value_type>
170  u_steps = max_near_feas_step(
171  x_k, *Ze
172  ,nlp.xl(), nlp.xu()
173  ,nlp.max_var_bounds_viol()
174  );
175 
176  if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) {
177  out << "\nMaximum steps ( possitive, negative ) in bounds u = ("
178  << u_steps.first << "," << u_steps.second << ")\n";
179  }
180 
181  if( u_steps.first < u )
182  u = u_steps.first;
183  if( ::fabs(u_steps.second) < u )
184  u = u_steps.second;
185  }
186 
187  if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) {
188  out << "\nFinite difference step length u = " << u << std::endl;
189  }
190 
191  // Take the finite difference from x_k to x_fd = x_k + u * Ze
192  //
193  // rGf_fd = ( Z_k'*Gf(x_k + u*Ze) - rGf_k ) / u
194  //
195 
196  VectorSpace::vec_mut_ptr_t x_fd = x_k.space().create_member();
197  Vp_StV( x_fd.get(), u, *Ze );
198 
199  // Gf_fd = Gf(x_fd)
200  VectorSpace::vec_mut_ptr_t Gf_fd = x_k.space().create_member();
201  nlp.unset_quantities();
202  nlp.set_Gf( Gf_fd.get() );
203  nlp.calc_Gf( *x_fd );
204 
205  if( (int)olevel >= (int)PRINT_VECTORS ) {
206  out << "\nGf_fd =\n" << *Gf_fd;
207  }
208 
209  // rGf_fd = Z'*Gf_fd
210  VectorSpace::vec_mut_ptr_t rGf_fd = Z_k.space_rows().create_member();
211  V_MtV( rGf_fd.get(), Z_k, BLAS_Cpp::trans, *Gf_fd );
212 
213  // rGf_fd = rGf_fd - rGf_k
214  Vp_StV( rGf_fd.get(), -1.0, rGf_k );
215 
216  // rGf_fd = rGf_fd / u
217  Vt_S( rGf_fd.get(), 1.0 / u );
218 
219  const value_type
220  nrm_rGf_fd = rGf_fd->norm_inf();
221 
222  if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) {
223  out << "\n||(rGf_fd - rGf_k)/u||inf = " << nrm_rGf_fd << std::endl;
224  }
225  if( (int)olevel >= (int)PRINT_VECTORS ) {
226  out << "\n(rGf_fd - rGf_k)/u =\n" << *rGf_fd;
227  }
228 
229  if( nrm_rGf_fd <= min_diag() ) {
230  if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) {
231  out << "\n||(rGf_fd - rGf_k)/u||inf = " << nrm_rGf_fd
232  << " < min_diag = " << min_diag() << std::endl
233  << "\nScale by min_diag ... \n";
234  }
235  rHL_diag.init_identity(Z_k.space_rows(),min_diag());
236  }
237  else {
238  switch( initialization_method() ) {
239  case SCALE_IDENTITY: {
240  if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) {
241  out << "\nScale the identity matrix by ||(rGf_fd - rGf_k)/u||inf ... \n";
242  }
243  rHL_diag.init_identity(Z_k.space_rows(),nrm_rGf_fd);
244  break;
245  }
246  case SCALE_DIAGONAL:
247  case SCALE_DIAGONAL_ABS:
248  {
249  if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) {
250  out << "\nScale diagonal by modified finite difference ... \n";
251  }
252  // In order to keep the initial reduced Hessian well conditioned we
253  // will not let any diagonal element drop below
254  // ||rGf_fd||inf / max_cond
255 
256  const value_type
257  min_ele = my_max( nrm_rGf_fd / max_cond(), min_diag() );
258 
259  if( initialization_method() == SCALE_DIAGONAL )
260  AbstractLinAlgPack::max_vec_scalar( min_ele, rGf_fd.get() );
261  else
262  AbstractLinAlgPack::max_abs_vec_scalar( min_ele, rGf_fd.get() );
263 
264  if( (int)olevel >= (int)PRINT_ALGORITHM_STEPS ) {
265  out << "\n||diag||inf = " << rGf_fd->norm_inf() << std::endl;
266  }
267  if( (int)olevel >= (int)PRINT_VECTORS ) {
268  out << "\ndiag =\n" << *rGf_fd;
269  }
270  rHL_diag.init_diagonal(*rGf_fd);
271  break;
272  }
273  default:
274  TEUCHOS_TEST_FOR_EXCEPT(true); // only local programming error?
275  }
276  }
277  nlp.unset_quantities();
278 
279  quasi_newton_stats_(s).set_k(0).set_updated_stats(
281 
282  if( (int)olevel >= (int)PRINT_ITERATION_QUANTITIES ) {
283  out << "\nrHL_k =\n" << rHL_iq.get_k(0);
284  }
285 
286  }
287 
288  return true;
289 }
290 
292  const Algorithm& algo
293  ,poss_type step_poss, IterationPack::EDoStepType type, poss_type assoc_step_poss
294  ,std::ostream& out, const std::string& L
295  ) const
296 {
297  out
298  << L << "*** Initialize the reduced Hessian using a single finite difference.\n"
299  << L << "*** Where the nlp must support the NLPObjGrad interface and\n"
300  << L << "*** rHL_k must support the MatrixSymInitDiag interface or exceptions\n"
301  << L << "*** will be thrown.\n"
302  << L << "default: num_basis_remembered = NO_BASIS_UPDATED_YET\n"
303  << L << " initialization_method = SCALE_DIAGONAL\n"
304  << L << " max_cond = " << max_cond() << std::endl
305  << L << " min_diag = " << min_diag() << std::endl
306  << L << " step_scale = " << step_scale() << std::endl
307  << L << "if num_basis_k is updated then\n"
308  << L << " new_basis = true\n"
309  << L << "else\n"
310  << L << " new_basis = false\n"
311  << L << "end\n"
312  << L << "if new_basis == true or no past rHL as been updated then\n"
313  << L << " *** Reinitialize the reduced Hessian using finite differences\n"
314  << L << " Ze = Z * e\n"
315  << L << " u = step_scale / norm(Ze,inf)\n"
316  << L << " if there are bounds on the problem then\n"
317  << L << " Find the largest (in magnitude) positive (u_pos) and\n"
318  << L << " negative (u_neg) step u where the slightly relaxed variable\n"
319  << L << " bounds:\n"
320  << L << " xl - delta <= x_k + u * Ze <= xu + delta\n"
321  << L << " are strictly satisfied (where delta = max_var_bounds_viol).\n"
322  << L << " if u_pos < u then\n"
323  << L << " u = u_pos\n"
324  << L << " end\n"
325  << L << " if abs(u_neg) < u then\n"
326  << L << " u = u_neg\n"
327  << L << " end\n"
328  << L << " end\n"
329  << L << " x_fd = x_k + u * Ze\n"
330  << L << " rGf_fd = ( Z_k' * Gf(x_fd) - rGf_k ) / u\n"
331  << L << " if norm(rGf_fd,inf) <= min_diag then\n"
332  << L << " rHL_k = min_diag * eye(n-r)\n"
333  << L << " else\n"
334  << L << " if initialization_method == SCALE_IDENTITY then\n"
335  << L << " rHL_k = norm(rGf_fd,inf) * eye(n-r)\n"
336  << L << " else if initialization_method == SCALE_DIAGONAL or SCALE_DIAGONAL_ABS then\n"
337  << L << " *** Make sure that all of the diagonal elements are\n"
338  << L << " *** positive and that the smallest element is\n"
339  << L << " *** no smaller than norm(rGf_fd,inf) / max_cond\n"
340  << L << " *** So that rHL_k will be positive definite an\n"
341  << L << " *** well conditioned\n"
342  << L << " min_ele = max( norm(rGf_fd,inf)/max_cond, min_diag )\n"
343  << L << " if initialization_method == SCALE_DIAGONAL then\n"
344  << L << " for i = 1 ... n-r\n"
345  << L << " diag(i) = max( rGf_fd(i), min_ele )\n"
346  << L << " end\n"
347  << L << " else *** SCALE_DIAGONAL_ABS\n"
348  << L << " for i = 1 ... n-r\n"
349  << L << " diag(i) = max( abs(rGf_fd(i)), min_ele )\n"
350  << L << " end\n"
351  << L << " end\n"
352  << L << " rHL_k = diag(diag)\n"
353  << L << " end\n"
354  << L << " end\n"
355  << L << "end\n";
356 }
357 
358 } // end namespace MoochoPack
void max_vec_scalar(value_type min_ele, VectorMutable *y)
Take the maximum value of the vector elements and a scalar.
void Vt_S(VectorMutable *v_lhs, const value_type &alpha)
v_lhs *= alpha
bool do_step(Algorithm &algo, poss_type step_poss, IterationPack::EDoStepType type, poss_type assoc_step_poss)
void Vp_StV(VectorMutable *v_lhs, const value_type &alpha, const Vector &v_rhs)
v_lhs = alpha * v_rhs + v_lhs
void max_abs_vec_scalar(value_type min_ele, VectorMutable *y)
Take the maximum value of the absolute vector elements and a scalar.
rSQP Algorithm control class.
Transposed.
Not transposed.
virtual std::ostream & journal_out() const
Return a reference to a std::ostream to be used to output debug information and the like...
EJournalOutputLevel
enum for journal output.
T_To & dyn_cast(T_From &from)
void print_step(const Algorithm &algo, poss_type step_poss, IterationPack::EDoStepType type, poss_type assoc_step_poss, std::ostream &out, const std::string &leading_str) const
Reduced space SQP state encapsulation interface.
std::ostream * out
void print_algorithm_step(const Algorithm &algo, Algorithm::poss_type step_poss, EDoStepType type, Algorithm::poss_type assoc_step_poss, std::ostream &out)
Prints to 'out' the algorithm step.
AlgorithmTracker & track()
void V_MtV(VectorMutable *v_lhs, const MatrixOp &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const V &V_rhs2)
v_lhs = op(M_rhs1) * V_rhs2.
std::pair< value_type, value_type > max_near_feas_step(const Vector &x, const Vector &d, const Vector &xl, const Vector &xu, value_type max_bnd_viol)
Computes the maximum positive and negative step that can be taken that are within the relaxed bounds...
AbstractLinAlgPack::value_type value_type
Acts as the central hub for an iterative algorithm.
NLPAlgoState & rsqp_state()
<<std aggr>="">> members for algo_cntr
InitFinDiffReducedHessian_Step(EInitializationMethod initialization_method=SCALE_IDENTITY, value_type max_cond=1e+1, value_type min_diag=1e-8, value_type step_scale=1e-1)
const f_int f_dbl_prec const f_int f_int const f_int f_int const f_dbl_prec & u
NLPAlgo & rsqp_algo(Algorithm &algo)
Convert from a Algorithm to a NLPAlgo.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)