00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, Andrew W. Steiner 00005 00006 This file is part of O2scl. 00007 00008 Some of the source code in this file has been adapted from OOL 00009 library written by Sergio Drumond Ventura, Luis Alberto D'Afonseca, 00010 and Ricardo Biloti available from http://ool.sourceforge.net. 00011 00012 O2scl is free software; you can redistribute it and/or modify 00013 it under the terms of the GNU General Public License as published by 00014 the Free Software Foundation; either version 3 of the License, or 00015 (at your option) any later version. 00016 00017 O2scl is distributed in the hope that it will be useful, 00018 but WITHOUT ANY WARRANTY; without even the implied warranty of 00019 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00020 GNU General Public License for more details. 00021 00022 You should have received a copy of the GNU General Public License 00023 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00024 00025 ------------------------------------------------------------------- 00026 */ 00027 #ifndef O2SCL_OOL_MMIN_GENCAN_H 00028 #define O2SCL_OOL_MMIN_GENCAN_H 00029 00030 #include <o2scl/text_file.h> 00031 #include <o2scl/ovector_tlate.h> 00032 #include <o2scl/multi_funct.h> 00033 #include <o2scl/ool_constr_min.h> 00034 #include <gsl/gsl_math.h> 00035 00036 #ifndef DOXYGENP 00037 namespace o2scl { 00038 #endif 00039 00040 /** 00041 \brief Constrained minimization by the "GENCAN" method (OOL) 00042 00043 \note Not working yet 00044 */ 00045 template<class param_t, class func_t, class dfunc_t=func_t, 00046 class hfunc_t=func_t, class vec_t=ovector_view, class alloc_vec_t=ovector, 00047 class alloc_t=ovector_alloc> class ool_mmin_gencan : 00048 public ool_constr_min<param_t,func_t,dfunc_t,hfunc_t,vec_t, 00049 alloc_vec_t,alloc_t> { 00050 00051 #ifndef DOXYGEN_INTERNAL 00052 00053 protected: 00054 00055 /// Desc (default 1.0) 00056 double cg_src; 00057 00058 /// Temporary vector 00059 alloc_vec_t S; 00060 /// Temporary vector 00061 alloc_vec_t Y; 00062 /// Temporary vector 00063 alloc_vec_t D; 00064 /// Temporary vector 00065 alloc_vec_t cg_W; 00066 /// Temporary vector 00067 alloc_vec_t cg_R; 00068 /// Temporary vector 00069 alloc_vec_t cg_D; 00070 /// Temporary vector 00071 alloc_vec_t cg_Sprev; 00072 /// Temporary vector 00073 alloc_vec_t Xtrial; 00074 /// Temporary vector 00075 alloc_vec_t tnls_Xtemp; 00076 /// Temporary vector 00077 alloc_vec_t near_l; 00078 /// Temporary vector 00079 alloc_vec_t near_u; 00080 00081 /// Desc 00082 int *Ind; 00083 00084 #ifdef NEVER_DEFINED 00085 00086 /// Desc 00087 int spgls() { 00088 00089 gsl_vector *X = M->x; 00090 gsl_vector *gradient = M->gradient; 00091 00092 size_t nn = X->size; 00093 00094 /* Direct access to vector data */ 00095 double *l = st->L->data; 00096 double *u = st->U->data; 00097 double *d = st->D->data; 00098 double *x = X->data; 00099 00100 double *xtrial = st->Xtrial->data; 00101 00102 /* Internal variables */ 00103 size_t interp; 00104 size_t imax; 00105 00106 double alpha; 00107 double dinfn; 00108 double gtd; 00109 double ftrial; 00110 00111 /* Compute first trial point, spectral projected gradient direction, 00112 * and directional derivative <g,d> */ 00113 alpha = 1; 00114 00115 /* Xtrial = min{ U, max[ L, ( X-lambda G ) ] } */ 00116 gsl_vector_memcpy( st->Xtrial, X ); 00117 gsl_blas_daxpy( -(st->lambda), gradient, st->Xtrial ); 00118 conmin_vector_minofmax( st->n, xtrial, u, l, xtrial ); 00119 00120 /* D = Xtrial - X */ 00121 gsl_vector_memcpy( st->D, st->Xtrial ); 00122 gsl_vector_sub( st->D, X ); 00123 00124 /* Inifite norm of D and < G, D > */ 00125 imax = gsl_blas_idamax( st->D ); 00126 dinfn = fabs( gsl_vector_get( st->D, imax ) ); 00127 gsl_blas_ddot( gradient, st->D, >d ); 00128 00129 /* Evaluate objective function */ 00130 OOL_CONMIN_EVAL_F( M, st->Xtrial, ftrial ); 00131 00132 interp = 0; 00133 00134 /* While Armijo isn't satizfied repeat */ 00135 while (ftrial > M->f + P->gamma*alpha*gtd ) { 00136 00137 /* Test if the trial point has a function value lower than fmin */ 00138 if (ftrial < M->f ) { 00139 00140 M->f = ftrial; 00141 gsl_vector_memcpy( X, st->Xtrial ); 00142 00143 return OOL_UNBOUNDEDF; 00144 } 00145 00146 interp++; 00147 00148 if (alpha < P->sigma1 ) { 00149 alpha /= P->nint; 00150 } else { 00151 /* quadratic model */ 00152 double atemp = ( -gtd*alpha*alpha )/ 00153 ( 2.0*(ftrial-M->f-alpha*gtd) ); 00154 00155 if (atemp < P->sigma1 || atemp > P->sigma2*alpha ) { 00156 alpha /= P->nint; 00157 } else { 00158 alpha = atemp; 00159 } 00160 } 00161 00162 /* Compute new trial point 00163 * Xtrial = X + alpha D */ 00164 gsl_vector_memcpy( st->Xtrial, X ); 00165 gsl_blas_daxpy( alpha, st->D, st->Xtrial ); 00166 00167 /* Evaluate objective function */ 00168 OOL_CONMIN_EVAL_F( M, st->Xtrial, ftrial ); 00169 00170 /* Test whether at least mininterp interpolations were made 00171 * and two consecutive iterates are close enough */ 00172 if( interp > P->mininterp && 00173 are_close( nn, alpha, d, x, P->epsrel, P->epsabs )) { 00174 00175 M->f = ftrial; 00176 gsl_vector_memcpy( X, st->Xtrial ); 00177 00178 return OOL_FLSEARCH; 00179 } 00180 } 00181 00182 /* Set new values of f and X */ 00183 M->f = ftrial; 00184 gsl_vector_memcpy( X, st->Xtrial ); 00185 00186 return OOL_SUCCESS; 00187 } 00188 00189 /// Desc 00190 int tnls() { 00191 00192 gsl_vector *X = M->x; 00193 gsl_vector *gradient = M->gradient; 00194 gsl_vector *Xplus = st->Xtrial; 00195 00196 /* Direct access to vector data */ 00197 double *x = X->data; 00198 double *g = gradient->data; 00199 double *d = st->D->data; 00200 double *xplus = Xplus->data; 00201 00202 /* Constant values */ 00203 const size_t nind = st->nind; 00204 00205 /* Internal variables */ 00206 double fplus; 00207 double gtd; 00208 double alpha; 00209 double gptd; 00210 00211 /* Compute directional derivative */ 00212 gtd = cblas_ddot( (int)nind, g, 1, d, 1 ); 00213 00214 /* Compute first trial */ 00215 alpha = GSL_MIN( 1.0, st->tnls_amax ); 00216 00217 /* Xplus = X + alpha D */ 00218 conmin_vector_memcpy( nind, xplus, x ); 00219 cblas_daxpy( (int)nind, alpha, d, 1, xplus, 1 ); 00220 00221 /* Evaluate objective function */ 00222 fplus = conmin_calc_f( M, nind, st->Ind, Xplus, X ); 00223 00224 /* Test Armijo and beta-condition and decide for: 00225 * 1 - accepting the trial point, 00226 * 2 - interpolate or 00227 * 3 - extrapolate. */ 00228 if( st->tnls_amax > 1.0 ) { 00229 00230 /* X+D belongs to the interior of the feasible set (amax > 1) */ 00231 00232 /* Verify Armijo */ 00233 if( fplus <= M->f + P->gamma * alpha * gtd ) { 00234 00235 /* Armijo condition holds */ 00236 00237 /* Evaluate the gradient of objective function */ 00238 conmin_calc_g( M, nind, st->Ind, Xplus, X, gradient ); 00239 /* Eval gptd = < g, d > */ 00240 gptd = cblas_ddot( (int)nind, g, 1, d, 1 ); 00241 00242 /* Verify directional derivative (beta condition) */ 00243 if ( gptd >= P->beta * gtd ) { 00244 00245 /* Step = 1 was ok, finish the line search */ 00246 00247 M->f = fplus; 00248 conmin_vector_memcpy( nind, x, xplus ); 00249 00250 return OOL_SUCCESS; 00251 } else { 00252 return tnls_extrapolation( M, st, P, alpha, fplus ); 00253 } 00254 } else { 00255 return tnls_interpolation(M, st, P, alpha, fplus, gtd); 00256 } 00257 } else { 00258 /* x + d does not belong to the feasible set (amax <= 1) */ 00259 if( fplus < M->f ) { 00260 return tnls_extrapolation( M, st, P, alpha, fplus ); 00261 } else { 00262 return tnls_interpolation(M, st, P, alpha, fplus, gtd); 00263 } 00264 } 00265 } 00266 00267 /// Desc 00268 int tnls_extrapolation(double alpha, double fplus) { 00269 00270 gsl_vector *X = M->x; 00271 gsl_vector *gradient = M->gradient; 00272 gsl_vector *Xplus = st->Xtrial; 00273 00274 /* Direct access to vector data */ 00275 double *x = X->data; 00276 double *d = st->D->data; 00277 double *l = st->L->data; 00278 double *u = st->U->data; 00279 00280 double *xplus = Xplus->data; 00281 double *xtemp = st->tnls_Xtemp->data; 00282 00283 /* Constant values */ 00284 const size_t nind = st->nind; 00285 00286 /* Internal variables */ 00287 double atemp; 00288 double ftemp; 00289 00290 size_t ii, extrap; 00291 short same; 00292 00293 /* Iterations */ 00294 extrap = 0; 00295 do { 00296 00297 extrap++; 00298 00299 /* Test if maximum number of extrapolation was exceeded */ 00300 if ( extrap > P->maxextrap ) { 00301 00302 M->f = fplus; 00303 conmin_vector_memcpy( nind, x, xplus ); 00304 00305 if (extrap > 0 || st->tnls_amax < 1){ 00306 conmin_calc_g( M, nind, st->Ind, Xplus, X, gradient ); 00307 } 00308 return TNLS_MAXEXTRAP; 00309 } 00310 00311 /* Chose new step */ 00312 if (alpha < st->tnls_amax && st->tnls_amax < P->next*alpha ) { 00313 atemp = st->tnls_amax; 00314 } else { 00315 atemp = P->next * alpha; 00316 } 00317 00318 /* Compute new trial point. Xtemp = X + atemp*D */ 00319 conmin_vector_memcpy( nind, xtemp, x ); 00320 cblas_daxpy( (int)nind, atemp, d, 1, xtemp, 1 ); 00321 00322 /* Project */ 00323 if (atemp > st->tnls_amax ) { 00324 conmin_vector_maxofmin( nind, xtemp, l, xtemp, u ); 00325 } 00326 00327 /* Test if this is not the same point as the previous one. 00328 * This test is performed only when alpha >= amax. */ 00329 if( alpha > st->tnls_amax ) { 00330 00331 same = 1; 00332 for (ii = 0; ii<nind && same; ii++) { 00333 00334 double aux; 00335 00336 aux = P->epsrel * fabs( xplus[ii] ); 00337 00338 if ( fabs(xtemp[ii]-xplus[ii]) > GSL_MAX(aux,P->epsabs)) { 00339 same = 0; 00340 } 00341 } 00342 00343 if (same) { 00344 00345 /* Finish the extrapolation with the current point */ 00346 M->f = fplus; 00347 00348 conmin_vector_memcpy( nind, x, xplus ); 00349 00350 if (extrap > 0 || st->tnls_amax < 1){ 00351 conmin_calc_g( M, nind, st->Ind, Xplus, X, gradient ); 00352 } 00353 return OOL_SUCCESS; 00354 } 00355 } 00356 00357 ftemp = conmin_calc_f( M, nind, st->Ind, st->tnls_Xtemp, X ); 00358 00359 if (ftemp < fplus) { 00360 00361 /* If the functional value decreases then set the current 00362 * point and continue the extrapolation */ 00363 00364 alpha = atemp; 00365 fplus = ftemp; 00366 conmin_vector_memcpy( nind, xplus, xtemp ); 00367 00368 continue; 00369 00370 } else { 00371 00372 /* If the functional value does not decrease then discard the 00373 * last trial and finish the extrapolation with the previous 00374 * point */ 00375 00376 M->f = fplus; 00377 00378 conmin_vector_memcpy( nind, x, xplus ); 00379 if (extrap > 0 || st->tnls_amax < 1) { 00380 conmin_calc_g( M, nind, st->Ind, X, X, gradient ); 00381 } 00382 00383 return OOL_SUCCESS; 00384 } 00385 00386 } while (1); 00387 00388 /* Just to make gcc happy */ 00389 return OOL_SUCCESS; 00390 } 00391 00392 /// Desc 00393 int tnls_interpolation(double alpha, double fplus, double gtd) { 00394 00395 gsl_vector *X = M->x; 00396 gsl_vector *gradient = M->gradient; 00397 gsl_vector *Xplus = st->Xtrial; 00398 00399 /* Direct access to vector data */ 00400 double *x = X->data; 00401 double *d = st->D->data; 00402 double *xplus = Xplus->data; 00403 00404 /* Constant values */ 00405 const size_t nind = st->nind; 00406 00407 /* Internal variables */ 00408 size_t interp; 00409 double atemp; 00410 00411 /* Initialization */ 00412 interp = 0; 00413 00414 /* Iterations */ 00415 do { 00416 interp++; 00417 00418 /* Test Armijo condition */ 00419 if (fplus <= M->f + P->gamma * alpha * gtd ) { 00420 00421 /* Finish the line search */ 00422 M->f = fplus; 00423 00424 /* X = Xplus */ 00425 conmin_vector_memcpy( nind, x, xplus ); 00426 00427 /* Evaluate objective function gradient */ 00428 conmin_calc_g( M, nind, st->Ind, X, X, gradient ); 00429 00430 return OOL_SUCCESS; 00431 } 00432 /* Compute new step */ 00433 if (alpha < P->sigma1 ) { 00434 alpha= alpha / P->nint; 00435 } else { 00436 /* quadratic model */ 00437 atemp = -gtd * alpha*alpha / 00438 (2 * (fplus - M->f - alpha * gtd)); 00439 00440 if( atemp < P->sigma1 || 00441 atemp > P->sigma2*alpha ) { 00442 alpha = alpha / P->nint; 00443 } else { 00444 alpha = atemp; 00445 } 00446 } 00447 00448 /* Compute new trial point: xplus = x + alpha d */ 00449 conmin_vector_memcpy( nind, xplus, x ); 00450 cblas_daxpy( (int)nind, alpha, d, 1, xplus, 1 ); 00451 00452 /* Evaluate objective function */ 00453 fplus = conmin_calc_f( M, nind, st->Ind, Xplus, X ); 00454 00455 /* Test whether at least mininterp interpolations were made 00456 * and the steplength is soo small */ 00457 if ( are_close( nind, alpha, d, x, P->epsrel, P->epsabs ) && 00458 interp > P->mininterp ){ 00459 return OOL_FLSEARCH; 00460 } 00461 } 00462 while( 1 ); 00463 00464 /* Just to make gcc happy */ 00465 return OOL_SUCCESS; 00466 } 00467 00468 /** Truncated Newton maximum step*/ 00469 double tnls_maximum_step() { 00470 00471 /* Direct access to vector data */ 00472 double *x = M->x->data; 00473 double *l = st->L->data; 00474 double *u = st->U->data; 00475 double *d = st->D->data; 00476 00477 double step = P->infabs; 00478 size_t ii; 00479 00480 for( ii = 0; ii < st->nind; ii++ ) { 00481 00482 if( *d > 0 ) { 00483 double aux = ( *u - *x ) / *d; 00484 step = GSL_MIN( step, aux ); 00485 } else if( *d < 0 ) { 00486 double aux = ( *l - *x ) / *d; 00487 step = GSL_MIN( step, aux ); 00488 } 00489 00490 x++; 00491 l++; 00492 u++; 00493 d++; 00494 } 00495 00496 return step; 00497 } 00498 00499 /** Spectral step length */ 00500 void spg_steplength() { 00501 00502 if (st->sty <= 0.0) { 00503 st->lambda = GSL_MAX( 1.0, st->xeucn ) / sqrt( st->gpeucn2 ); 00504 } else { 00505 double aux; 00506 double ss = st->sts / st->sty; 00507 00508 aux = GSL_MAX( P->lspgmi, ss ); 00509 st->lambda = GSL_MIN( P->lspgma, aux ); 00510 } 00511 } 00512 00513 /** Iterate */ 00514 int actual_iterate() { 00515 00516 /* Direct access to vector data */ 00517 double *x = M->x->data; 00518 double *l = st->L->data; 00519 double *u = st->U->data; 00520 /* double *d = st->D->data; */ 00521 00522 /* Status of internal iterations */ 00523 int lsflag; 00524 int cgflag; 00525 00526 /* Internal variables */ 00527 size_t ii, imax; 00528 00529 /* Saving previous values */ 00530 gsl_vector_memcpy( st->S, M->x ); 00531 gsl_vector_memcpy( st->Y, M->gradient ); 00532 00533 /* The actual iteration */ 00534 if ( st->gieucn2 <= st->ometa2 * st->gpeucn2 ) { 00535 /* Compute the new iterate using an SPG iteration */ 00536 00537 /* Perform a line search with spectral continuous 00538 projected gradient */ 00539 lsflag = spgls( M, st, P ); 00540 00541 /* Compute the gradient for the new iterate */ 00542 OOL_CONMIN_EVAL_DF( M, M->x, M->gradient ); 00543 00544 } else { 00545 00546 /* The new iterate will belong to the closure of the actual face */ 00547 00548 /* Shrink the point, its gradient and the bounds */ 00549 conmin_shrink( st->nind, st->Ind, M->x ); 00550 conmin_shrink( st->nind, st->Ind, M->gradient ); 00551 conmin_shrink( st->nind, st->Ind, st->L ); 00552 conmin_shrink( st->nind, st->Ind, st->U ); 00553 00554 /* Compute the descent direction solving the newtonian system */ 00555 cgflag = cg( M, st, P ); 00556 00557 /* Compute maximum step for truncated newton line search */ 00558 if ( cgflag == CG_BOUNDARY ) { 00559 st->tnls_amax = 1.0; 00560 } else { 00561 st->tnls_amax = tnls_maximum_step( M, st, P ); 00562 } 00563 00564 /* Perform the line search */ 00565 lsflag = tnls( M, st, P ); 00566 00567 /* Expand the point, its gradient and the bounds */ 00568 conmin_expand( st->nind, st->Ind, M->x ); 00569 conmin_expand( st->nind, st->Ind, M->gradient ); 00570 conmin_expand( st->nind, st->Ind, st->L ); 00571 conmin_expand( st->nind, st->Ind, st->U ); 00572 00573 /* If the line search in the Truncated Newton direction 00574 stopped due to a very small step discard this iteration 00575 and force a SPG iteration */ 00576 if ( lsflag == OOL_FLSEARCH ) { 00577 00578 /* Perform a line search with spectral projected gradient */ 00579 lsflag = spgls( M, st, P ); 00580 00581 /* Compute the gradient for the new iterate */ 00582 OOL_CONMIN_EVAL_DF( M, M->x, M->gradient ); 00583 } 00584 } 00585 00586 /* Prepare for the next iteration */ 00587 /* Adjustment */ 00588 for( ii = 0; ii < st->n; ii++ ) { 00589 /* In principle, the current point could be slightly changed 00590 * here, requiring a new function and gradient 00591 * evaluation. But, according to the algorithms authors, this 00592 * is done just to account for points that are "numerically¨ 00593 * at faces already. Thus, no additional evaluations are 00594 * performed. (May 11th, 2005). 00595 */ 00596 if ( x[ii] <= st->near_l[ii] ) x[ii] = l[ii]; 00597 else if( x[ii] >= st->near_u[ii] ) x[ii] = u[ii]; 00598 } 00599 00600 /* Compute infinite and Euclidian-norm of X */ 00601 imax = gsl_blas_idamax( M->x ); 00602 st->xsupn = fabs( gsl_vector_get( M->x, imax ) ); 00603 st->xeucn = gsl_blas_dnrm2 ( M->x ); 00604 00605 /* Until now S = X_prev, now S = X - X_prev 00606 * Compute s = x_{k+1} - x_k = X - S 00607 * and y = g_{k+1} - g_k = G - Y */ 00608 gsl_vector_sub ( st->S, M->x ); /* S = S - X */ 00609 gsl_vector_scale( st->S, -1.0 ); /* S = -S = X - S_prev */ 00610 gsl_vector_sub ( st->Y, M->gradient ); /* Y = Y - G */ 00611 gsl_vector_scale( st->Y, -1.0 ); /* Y = -Y = G - Y_prev */ 00612 00613 /* Compute sts = s dot s 00614 * sty = s dot y 00615 * and sinf = |s|_inf */ 00616 gsl_blas_ddot( st->S, st->S, &(st->sts) ); 00617 gsl_blas_ddot( st->S, st->Y, &(st->sty) ); 00618 imax = gsl_blas_idamax( st->S ); 00619 st->sinf = fabs( gsl_vector_get( st->S, imax ) ); 00620 00621 /* Compute continuous project gradient */ 00622 projected_gradient( st, M->x, M->gradient ); 00623 00624 /* Update spectral steplength */ 00625 spg_steplength( st, P ); 00626 00627 /* Update trust-region radius */ 00628 if ( P->trtype ) st->cg_delta = GSL_MAX( P->delmin, 10*sqrt( st->sts ) ); 00629 else st->cg_delta = GSL_MAX( P->delmin, 10* ( st->sinf) ); 00630 00631 return lsflag; 00632 } 00633 00634 #endif 00635 00636 #endif 00637 00638 public: 00639 00640 ool_mmin_gencan() { 00641 epsgpen=1.0e-5; 00642 epsgpsn=1.0e-5; 00643 fmin=-1.0e99; 00644 udelta0=-1.0; 00645 ucgmia=-1.0; 00646 ucgmib=-1.0; 00647 cg_src=1.0; 00648 cg_scre=1.0; 00649 cg_gpnf=1.0e-5; 00650 cg_epsi=1.0e-1; 00651 cg_epsf=1.0e-5; 00652 cg_epsnqmp=1.0e-4; 00653 cg_maxitnqmp=5; 00654 nearlyq=0; 00655 nint=2.0; 00656 next=2.0; 00657 mininterp=4; 00658 maxextrap=100; 00659 trtype=0; 00660 eta=0.9; 00661 delmin=0.1; 00662 lspgmi=1.0e-10; 00663 lspgma=1.0e10; 00664 theta=1.0e-6; 00665 gamma=1.0e-4; 00666 beta=0.5; 00667 sigma1=0.1; 00668 sigma2=0.9; 00669 epsrel=1.0e-7; 00670 epsabs=1.0e-10; 00671 infrel=1.0e20; 00672 infabs=1.0e99; 00673 } 00674 00675 /// Tolerance on Euclidean norm of projected gradient (default 1.0e-5) 00676 double epsgpen; 00677 /// Tolerance on infinite norm of projected gradient (default 1.0e-5) 00678 double epsgpsn; 00679 /** 00680 \brief Minimum function value (default \f$ 10^{-99} \f$) 00681 00682 If the function value is below this value, then the algorithm 00683 assumes that the function is not bounded and exits. 00684 */ 00685 double fmin; 00686 /// Trust-region radius (default -1.0) 00687 double udelta0; 00688 /// Maximum interations per variable (default -1.0) 00689 double ucgmia; 00690 /// Extra maximum iterations (default -1.0) 00691 double ucgmib; 00692 /// Conjugate gradient condition type (default 1) 00693 int cg_scre; 00694 /// Projected gradient norm (default 1.0e-5) 00695 double cg_gpnf; 00696 /// Desc (default 1.0e-1) 00697 double cg_epsi; 00698 /// Desc (default 1.0e-5) 00699 double cg_epsf; 00700 /// Stopping fractional tolerance for conjugate gradient (default 1.0e-4) 00701 double cg_epsnqmp; 00702 /// Maximum iterations for conjugate gradient (default 5) 00703 int cg_maxitnqmp; 00704 /// Set to 1 if the function is nearly quadratic (default 0) 00705 int nearlyq; 00706 /// Interpolation constant (default 2.0) 00707 double nint; 00708 /// Extrapolation constant (default 2.0) 00709 double next; 00710 /// Minimum interpolation size (default 4) 00711 int mininterp; 00712 /// Maximum extrapolations in truncated Newton direction (default 100) 00713 int maxextrap; 00714 /// Type of trust region (default 0) 00715 int trtype; 00716 /// Threshold for abandoning current face (default 0.9) 00717 double eta; 00718 /// Minimum trust region for truncated Newton direction (default 0.1) 00719 double delmin; 00720 /// Minimum spectral steplength (default 1.0e-10) 00721 double lspgmi; 00722 /// Maximum spectral steplength (default 1.0e10) 00723 double lspgma; 00724 /// Constant for the angle condition (default 1.0e-6) 00725 double theta; 00726 /// Constant for Armijo condition (default 1.0e-4) 00727 double gamma; 00728 /// Constant for beta condition (default 0.5) 00729 double beta; 00730 /// Lower bound to the step length reduction (default 0.1) 00731 double sigma1; 00732 /// Upper bound to the step length reduction (default 0.9) 00733 double sigma2; 00734 /// Relative small number (default 1.0e-7) 00735 double epsrel; 00736 /// Absolute small number (default 1.0e-10) 00737 double epsabs; 00738 /// Relative infinite number (default 1.0e20) 00739 double infrel; 00740 /// Absolute infinite number (default 1.0e99) 00741 double infabs; 00742 00743 /// Allocate memory 00744 virtual int alloc(const size_t n) { 00745 if (this->dim>0) free(); 00746 this->ao.allocate(xx,n); 00747 this->ao.allocate(d,n); 00748 this->ao.allocate(s,n); 00749 this->ao.allocate(y,n); 00750 return ool_constr_min<param_t,func_t,dfunc_t,hfunc_t,vec_t,alloc_vec_t, 00751 alloc_t>::alloc(n); 00752 } 00753 00754 /// Free previously allocated memory 00755 virtual int free() { 00756 if (this->dim>0) this->ao.free(xx); 00757 return ool_constr_min<param_t,func_t,dfunc_t,hfunc_t,vec_t,alloc_vec_t, 00758 alloc_t>::free(); 00759 } 00760 00761 /// Set the function, the initial guess, and the parameters 00762 virtual int set(func_t &fn, dfunc_t &dfn, hfunc_t &hfn, 00763 vec_t &init, param_t &par) { 00764 00765 int ret=ool_constr_min<param_t,func_t,dfunc_t,hfunc_t,vec_t,alloc_vec_t, 00766 alloc_t>::set(fn,dfn,hfn,init,par); 00767 00768 #ifdef NEVER_DEFINED 00769 // Internal variables 00770 size_t nn = M->x->size; 00771 00772 // Checking dimensions 00773 if( nn != st->n || nn != M->fdf->n || nn != M->con->n ) 00774 { 00775 return OOL_EINVAL; 00776 } 00777 00778 // Copy boundary vectors 00779 gsl_vector_memcpy( st->L, M->con->L ); 00780 gsl_vector_memcpy( st->U, M->con->U ); 00781 00782 #endif 00783 00784 prepare_iteration(); 00785 00786 return 0; 00787 } 00788 00789 #ifdef NEVER_DEFINED 00790 00791 int prepare_iteration { 00792 00793 /* Direct access to vector data */ 00794 double *x = M->x->data; 00795 double *l = st->L->data; 00796 double *u = st->U->data; 00797 00798 /* Internal variables */ 00799 size_t nn = M->x->size; 00800 size_t ii, imax; 00801 00802 /* Impose factibility */ 00803 conmin_vector_maxofmin( nn, x, l, u, x ); 00804 00805 /* Eval Euclidean and infinity norms of X */ 00806 st->xeucn = gsl_blas_dnrm2 ( M->x ); 00807 imax = gsl_blas_idamax( M->x ); 00808 st->xsupn = fabs( gsl_vector_get( M->x, imax ) ); 00809 00810 /* Evaluate objective function and its gradient */ 00811 OOL_CONMIN_EVAL_FDF( M, M->x, &(M->f), M->gradient ); 00812 00813 /* Define near_l and near_u vector */ 00814 for (ii=0; ii < nn; ii++){ 00815 st->near_l[ii] = l[ii] + GSL_MAX( P->epsrel*fabs( l[ii] ), P->epsabs ); 00816 st->near_u[ii] = u[ii] - GSL_MAX( P->epsrel*fabs( u[ii] ), P->epsabs ); 00817 } 00818 00819 /* Setting constant parameters */ 00820 st->ometa2 = gsl_pow_2( 1.0 - P->eta ); 00821 st->epsgpen2 = gsl_pow_2( P->epsgpen ); 00822 00823 /* Compute continuous project gradient */ 00824 projected_gradient( st, M->x, M->gradient ); 00825 00826 /* Compute a linear relation between gpeucn2 and cgeps2, i.e., 00827 * scalars a and b such that 00828 * 00829 * a * log10(||g_P(x_0)||_2^2) + b = log10(cgeps_0^2) and 00830 * 00831 * a * log10(||g_P(x_f)||_2^2) + b = log10(cgeps_f^2), 00832 * 00833 * where cgeps_0 and cgeps_f are provided. Note that if 00834 * cgeps_0 is equal to cgeps_f then cgeps will be always 00835 * equal to cgeps_0 and cgeps_f. 00836 * 00837 * We introduce now a linear relation between gpsupn and cgeps also. 00838 */ 00839 if (P->cg_scre == 1 ) { 00840 st->acgeps = 2 *( log10( P->cg_epsf / P->cg_epsi ) / 00841 log10( P->cg_gpnf * P->cg_gpnf / st->gpeucn2 )); 00842 00843 st->bcgeps = 2 * log10( P->cg_epsi ) - 00844 st->acgeps * log10( st->gpeucn2 ); 00845 } else { 00846 st->acgeps = ( log10( P->cg_epsf / P->cg_epsi ) / 00847 log10( P->cg_gpnf / st->gpsupn ) ); 00848 st->bcgeps = log10( P->cg_epsi ) - st->acgeps * log10( st->gpsupn ); 00849 } 00850 00851 /* And it will be used for the linear relation of cgmaxit */ 00852 st->gpsupn0 = st->gpsupn; 00853 st->gpeucn20 = st->gpeucn2; 00854 00855 /* Initialize the spectral steplength */ 00856 if ( st->gpeucn2 != 0.0 ) { 00857 st->lambda = GSL_MAX( 1.0, st->xeucn ) / sqrt( st->gpeucn2 ); 00858 } 00859 00860 /* Initialize the trust-region radius */ 00861 if (P->udelta0 < 0.0 ) { 00862 00863 double aux; 00864 if ( P->trtype ) { 00865 aux = 0.1 * GSL_MAX( 1.0, st->xeucn ); 00866 } else { 00867 aux = 0.1 * GSL_MAX( 1.0, st->xsupn ); 00868 } 00869 00870 st->cg_delta = GSL_MAX( P->delmin, aux ); 00871 00872 } else { 00873 st->cg_delta = GSL_MAX( P->delmin, P->udelta0 ); 00874 } 00875 00876 /* Export size */ 00877 M->size = st->gpsupn; 00878 00879 return OOL_SUCCESS; 00880 } 00881 00882 #endif 00883 00884 /// Restart the minimizer 00885 virtual int restart() { 00886 00887 /* 00888 // Restarting dx 00889 gsl_vector_set_zero( M->dx ); 00890 00891 return prepare_iteration( M ); 00892 */ 00893 00894 return 0; 00895 } 00896 00897 /// Perform an iteration 00898 virtual int iterate() { 00899 00900 #ifdef NEVER_DEFINED 00901 00902 int status; 00903 00904 status = actual_iterate( M, st, P ); 00905 00906 /* Export size and dx variables */ 00907 M->size = st->gpsupn; 00908 00909 /* In the next version does dx replace st->S ? */ 00910 gsl_vector_memcpy( M->dx, st->S ); 00911 00912 return status; 00913 00914 #endif 00915 00916 return 0; 00917 } 00918 00919 /// See if we're finished 00920 virtual int is_optimal() { 00921 00922 //return (( st->gpeucn2 <= st->epsgpen2 || 00923 //st->gpsupn <= P->epsgpsn || 00924 //M->f <= P->fmin )? OOL_SUCCESS : OOL_CONTINUE ); 00925 00926 } 00927 00928 /// Return string denoting type ("ool_mmin_gencan") 00929 const char *type() { return "ool_mmin_gencan"; } 00930 00931 }; 00932 00933 #ifndef DOXYGENP 00934 } 00935 #endif 00936 00937 #endif 00938
Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.
Project hosting provided by
,
O2scl Sourceforge Project Page