Actual source code: lmvmmat.c

  1: #include "lmvmmat.h"
  2: #include "tao_general.h"
  3: #include "tao_solver.h"
  4: #include "taovec.h"

  6: #define Scale_None                0
  7: #define Scale_Scalar                1
  8: #define Scale_Broyden                2
  9: #define Scale_Types             3

 11: #define Rescale_None                0
 12: #define Rescale_Scalar                1
 13: #define Rescale_GL                2
 14: #define Rescale_Types                  3

 16: #define Limit_None                0
 17: #define Limit_Average                1
 18: #define Limit_Relative                2
 19: #define Limit_Absolute                3
 20: #define Limit_Types                4

 22: #define TAO_ZER_SAFEGUARD        1e-8
 23: #define TAO_INF_SAFEGUARD        1e+8

 25: static const char *Scale_Table[64] = {
 26:   "none", "scalar", "broyden"
 27: };

 29: static const char *Rescale_Table[64] = {
 30:   "none", "scalar", "gl"
 31: };

 33: static const char *Limit_Table[64] = {
 34:   "none", "average", "relative", "absolute"
 35: };

 39: TaoLMVMMat::TaoLMVMMat(TaoVec *tv)
 40: {
 41:   // Set default values
 42:   lm = 5;
 43:   eps = 0.0;

 45:   limitType = Limit_None;
 46:   scaleType = Scale_Broyden;
 47:   rScaleType = Rescale_Scalar;

 49:   s_alpha = 1.0;
 50:   r_alpha = 1.0;
 51:   r_beta = 0.5;
 52:   mu = 1.0;
 53:   nu = 100.0;

 55:   phi = 0.125;                

 57:   scalar_history = 1;
 58:   rescale_history = 1;

 60:   delta_min = 1e-7;
 61:   delta_max = 100;

 63:   // Begin configuration
 64:   TaoOptionsHead("Limited memory matrix options");
 65:   TaoOptionInt("-tao_lmm_vectors", "vectors to use for approximation", "", lm, &lm, 0);
 66:   TaoOptionDouble("-tao_lmm_limit_mu", "mu limiting factor", "", mu, &mu, 0);
 67:   TaoOptionDouble("-tao_lmm_limit_nu", "nu limiting factor", "", nu, &nu, 0);
 68:   TaoOptionDouble("-tao_lmm_broyden_phi", "phi factor for Broyden scaling", "", phi, &phi, 0);
 69:   TaoOptionDouble("-tao_lmm_scalar_alpha", "alpha factor for scalar scaling", "", s_alpha, &s_alpha, 0);
 70:   TaoOptionDouble("-tao_lmm_rescale_alpha", "alpha factor for rescaling diagonal", "", r_alpha, &r_alpha, 0);
 71:   TaoOptionDouble("-tao_lmm_rescale_beta", "beta factor for rescaling diagonal", "", r_beta, &r_beta, 0);
 72:   TaoOptionInt("-tao_lmm_scalar_history", "amount of history for scalar scaling", "", scalar_history, &scalar_history, 0);
 73:   TaoOptionInt("-tao_lmm_rescale_history", "amount of history for rescaling diagonal", "", rescale_history, &rescale_history, 0);
 74:   TaoOptionDouble("-tao_lmm_eps", "rejection tolerance", "", eps, &eps, 0);
 75:   TaoOptionList("-tao_lmm_scale_type", "scale type", "", Scale_Table, Scale_Types, Scale_Table[scaleType], &scaleType, 0);
 76:   TaoOptionList("-tao_lmm_rescale_type", "rescale type", "", Rescale_Table, Rescale_Types, Rescale_Table[rScaleType], &rScaleType, 0);
 77:   TaoOptionList("-tao_lmm_limit_type", "limit type", "", Limit_Table, Limit_Types, Limit_Table[limitType], &limitType, 0);
 78:   TaoOptionDouble("-tao_lmm_delta_min", "minimum delta value", "", delta_min, &delta_min, 0);
 79:   TaoOptionDouble("-tao_lmm_delta_max", "maximum delta value", "", delta_max, &delta_max, 0);
 80:   TaoOptionsTail();

 82:   // Complete configuration
 83:   rescale_history = TaoMin(rescale_history, lm);

 85:   // Perform allocations
 86:   tv->CloneVecs(lm+1, &S);
 87:   tv->CloneVecs(lm+1, &Y);
 88:   tv->Clone(&D);
 89:   tv->Clone(&U);
 90:   tv->Clone(&V);
 91:   tv->Clone(&W);
 92:   tv->Clone(&P);
 93:   tv->Clone(&Q);

 95:   rho = new double[lm+1];
 96:   beta = new double[lm+1];

 98:   yy_history = new double[TaoMax(scalar_history, 1)];
 99:   ys_history = new double[TaoMax(scalar_history, 1)];
100:   ss_history = new double[TaoMax(scalar_history, 1)];

102:   yy_rhistory = new double[TaoMax(rescale_history, 1)];
103:   ys_rhistory = new double[TaoMax(rescale_history, 1)];
104:   ss_rhistory = new double[TaoMax(rescale_history, 1)];

106:   // Finish initializations
107:   lmnow = 0;
108:   iter = 0;
109:   updates = 0;
110:   rejects = 0;
111:   delta = 1.0;

113:   Gprev = 0;
114:   Xprev = 0;

116:   scale = 0;

118:   H0 = 0;
119:   H0default = TAO_TRUE;
120: }

124: TaoLMVMMat::~TaoLMVMMat()
125: {
126:   int info;

128:   info = S[0]->DestroyVecs(lm+1, S);
129:   info = Y[0]->DestroyVecs(lm+1, Y);
130:   info = TaoVecDestroy(D);
131:   info = TaoVecDestroy(U);
132:   info = TaoVecDestroy(V);
133:   info = TaoVecDestroy(W);
134:   info = TaoVecDestroy(P);
135:   info = TaoVecDestroy(Q);

137:   delete[] rho;
138:   delete[] beta;
139:   if (info) rho=0;

141:   delete[] yy_history;
142:   delete[] ys_history;
143:   delete[] ss_history;

145:   delete[] yy_rhistory;
146:   delete[] ys_rhistory;
147:   delete[] ss_rhistory;
148: }

152: int TaoLMVMMat::Refine(TaoLMVMMat *coarse, TaoMat *op, TaoVec *fineX, TaoVec *fineG)
153: {
154:   double rhotemp, rhotol;
155:   double y0temp, s0temp;
156:   int i, info;

158:   TaoFunctionBegin;

160:   info = Reset(); CHKERRQ(info);

162:   for (i = 0; i < coarse->lmnow; ++i) {
163:     // Refine S[i] and Y[i]
164:     info = op->Multiply(coarse->S[i], S[lmnow]); CHKERRQ(info);
165:     info = op->Multiply(coarse->Y[i], Y[lmnow]); CHKERRQ(info);

167:     // Check to see if refined vectors are fine
168:     info = Y[lmnow]->Dot(S[lmnow], &rhotemp); CHKERRQ(info);
169:     info = Y[lmnow]->Dot(Y[lmnow], &y0temp); CHKERRQ(info);

171:     rhotol = eps * y0temp;
172:     if (rhotemp > rhotol) {
173:       rho[lmnow] = 1.0 / rhotemp;

175:       // Add information to the scaling
176:       switch(scaleType) {
177:       case Scale_None:
178:         break;
179:         
180:       case Scale_Scalar:
181:         // Compute s^T s 
182:         info = S[lmnow]->Dot(S[lmnow], &s0temp); CHKERRQ(info);

184:         // Save information for scalar scaling
185:         yy_history[updates % scalar_history] = y0temp;
186:         ys_history[updates % scalar_history] = rhotemp;
187:         ss_history[updates % scalar_history] = s0temp;

189:         sigma = coarse->sigma;
190:         break;

192:       case Scale_Broyden:
193:         switch(rScaleType) {
194:         case Rescale_None:
195:           break;
196:           
197:         case Rescale_Scalar:
198:         case Rescale_GL:
199:           // Compute s^T s 
200:           info = S[lmnow]->Dot(S[lmnow], &s0temp); CHKERRQ(info);

202:           // Save information for special cases of scalar rescaling
203:           yy_rhistory[updates % rescale_history] = y0temp;
204:           ys_rhistory[updates % rescale_history] = rhotemp;
205:           ss_rhistory[updates % rescale_history] = s0temp;
206:           break;
207:         }

209:         // Obtain diagonal scaling and ensure positive definite
210:         info = op->Multiply(coarse->D, D); CHKERRQ(info);
211:         info = D->AbsoluteValue(); CHKERRQ(info);
212:         break;
213:       }

215:       // Successful update
216:       ++updates;
217:       ++lmnow;
218:     }
219:     else {
220:       ++rejects;
221:     }
222:   }

224:   // Save the number of iterations and previous values for x and g
225:   // Need to use the actual value for the gradient here and not
226:   // the refined version of the coarse gradient in order for the
227:   // Wolfe conditions to guarantee the update is acceptable.
228:   iter = coarse->iter;
229:   info = Xprev->CopyFrom(fineX); CHKERRQ(info);
230:   info = Gprev->CopyFrom(fineG); CHKERRQ(info);
231:   TaoFunctionReturn(0);
232: }

236: int TaoLMVMMat::Reset()
237: {
238:   int i;

240:   TaoFunctionBegin;
241:   Gprev = Y[lm];
242:   Xprev = S[lm];
243:   for (i = 0; i < lm; ++i) {
244:     rho[i] = 0.0;
245:   }
246:   rho[0] = 1.0;

248:   // Set the scaling and diagonal scaling matrix
249:   switch(scaleType) {
250:   case Scale_None:
251:     sigma = 1.0;
252:     break;

254:   case Scale_Scalar:
255:     sigma = delta;
256:     break;

258:   case Scale_Broyden:
259:     D->SetToConstant(delta);
260:     break;
261:   }

263:   iter = 0;
264:   updates = 0;
265:   lmnow = 0;
266:   TaoFunctionReturn(0);
267: }

271: int TaoLMVMMat::Presolve()
272: {
273:   TaoFunctionBegin;
274:   TaoFunctionReturn(0);
275: }

279: int TaoLMVMMat::Solve(TaoVec *tb, TaoVec *dx, TaoTruth *tt)
280: {
281:   double   sq, yq, dd;
282:   int      ll, info;
283:   TaoTruth scaled;

285:   TaoFunctionBegin;
286:   if (lmnow < 1) {
287:     rho[0] = 1.0;
288:   }

290:   info = dx->CopyFrom(tb); CHKERRQ(info);
291:   for (ll = 0; ll < lmnow; ++ll) {
292:     info = dx->Dot(S[ll], &sq); CHKERRQ(info);
293:     beta[ll] = sq * rho[ll];
294:     info = dx->Axpy(-beta[ll], Y[ll]); CHKERRQ(info);
295:   }

297:   scaled = TAO_FALSE;
298:   if (!scaled && H0) {
299:     info = H0->Solve(dx, U, tt); CHKERRQ(info);
300:     if (*tt) {
301:       info = dx->Dot(U, &dd); CHKERRQ(info);
302:       if ((dd > 0.0) && !TaoInfOrNaN(dd)) {
303:         // Accept Hessian solve
304:         info = dx->CopyFrom(U); CHKERRQ(info);
305:         scaled = TAO_TRUE;
306:       }
307:     }
308:   }

310:   if (!scaled && scale) {
311:     info = U->PointwiseMultiply(dx, scale); CHKERRQ(info);
312:     info = dx->Dot(U, &dd); CHKERRQ(info);
313:     if ((dd > 0.0) && !TaoInfOrNaN(dd)) {
314:       // Accept scaling
315:       info = dx->CopyFrom(U); CHKERRQ(info);
316:       scaled = TAO_TRUE;
317:     }
318:   }
319:   
320:   if (!scaled) {
321:     switch(scaleType) {
322:     case Scale_None:
323:       break;

325:     case Scale_Scalar:
326:       info = dx->Scale(sigma); CHKERRQ(info);
327:       break;
328:   
329:     case Scale_Broyden:
330:       info = dx->PointwiseMultiply(dx, D); CHKERRQ(info);
331:       break;
332:     }
333:   } 

335:   for (ll = lmnow-1; ll >= 0; --ll) {
336:     info = dx->Dot(Y[ll], &yq); CHKERRQ(info);
337:     info = dx->Axpy(beta[ll] - yq * rho[ll], S[ll]); CHKERRQ(info);
338:   }

340:   *tt = TAO_TRUE;
341:   TaoFunctionReturn(0);
342: }

346: int TaoLMVMMat::Update(TaoVec *x, TaoVec *g)
347: {
348:   double rhotemp, rhotol;
349:   double y0temp, s0temp;
350:   double yDy, yDs, sDs;
351:   double sigmanew, denom;
352:   int i, info;

354:   double yy_sum=0, ys_sum=0, ss_sum=0;

356:   TaoFunctionBegin;

358:   if (0 == iter) {
359:     info = Reset(); CHKERRQ(info);
360:   } 
361:   else {
362:     info = Gprev->Aypx(-1.0, g); CHKERRQ(info);
363:     info = Xprev->Aypx(-1.0, x); CHKERRQ(info);

365:     info = Gprev->Dot(Xprev, &rhotemp); CHKERRQ(info);
366:     info = Gprev->Dot(Gprev, &y0temp); CHKERRQ(info);

368:     rhotol = eps * y0temp;
369:     if (rhotemp > rhotol) {
370:       ++updates;

372:       lmnow = TaoMin(lmnow+1, lm);
373:       for (i = lm-1; i >= 0; --i) {
374:         S[i+1] = S[i];
375:         Y[i+1] = Y[i];
376:         rho[i+1] = rho[i];
377:       }
378:       S[0] = Xprev;
379:       Y[0] = Gprev;
380:       rho[0] = 1.0 / rhotemp;

382:       // Compute the scaling
383:       switch(scaleType) {
384:       case Scale_None:
385:         break;

387:       case Scale_Scalar:
388:         // Compute s^T s 
389:         info = Xprev->Dot(Xprev, &s0temp); CHKERRQ(info);

391:         // Scalar is positive; safeguards are not required.

393:         // Save information for scalar scaling
394:         yy_history[(updates - 1) % scalar_history] = y0temp;
395:         ys_history[(updates - 1) % scalar_history] = rhotemp;
396:         ss_history[(updates - 1) % scalar_history] = s0temp;

398:         // Compute summations for scalar scaling
399:         yy_sum = 0;        // No safeguard required; y^T y > 0
400:         ys_sum = 0;        // No safeguard required; y^T s > 0
401:         ss_sum = 0;        // No safeguard required; s^T s > 0
402:         for (i = 0; i < TaoMin(updates, scalar_history); ++i) {
403:           yy_sum += yy_history[i];
404:           ys_sum += ys_history[i];
405:           ss_sum += ss_history[i];
406:         }

408:         if (0.0 == s_alpha) {
409:           // Safeguard ys_sum 
410:           if (0.0 == ys_sum) {
411:             ys_sum = TAO_ZER_SAFEGUARD;
412:           }

414:           sigmanew = ss_sum / ys_sum;
415:         }
416:         else if (1.0 == s_alpha) {
417:           // Safeguard yy_sum 
418:           if (0.0 == yy_sum) {
419:             yy_sum = TAO_ZER_SAFEGUARD;
420:           }

422:           sigmanew = ys_sum / yy_sum;
423:         }
424:         else {
425:           denom = 2*s_alpha*yy_sum;

427:           // Safeguard denom
428:           if (0.0 == denom) {
429:             denom = TAO_ZER_SAFEGUARD;
430:           }

432:           sigmanew = ((2*s_alpha-1)*ys_sum + 
433:                       sqrt((2*s_alpha-1)*(2*s_alpha-1)*ys_sum*ys_sum - 
434:                            4*s_alpha*(s_alpha-1)*yy_sum*ss_sum)) / denom;
435:         }

437:         switch(limitType) {
438:         case Limit_Average:
439:           if (1.0 == mu) {
440:             sigma = sigmanew;
441:           }
442:           else if (mu) {
443:             sigma = mu * sigmanew + (1.0 - mu) * sigma;
444:           }
445:           break;

447:         case Limit_Relative:
448:           if (mu) {
449:             sigma = TaoMid((1.0 - mu) * sigma, sigmanew, (1.0 + mu) * sigma);
450:           }
451:           break;

453:         case Limit_Absolute:
454:           if (nu) {
455:             sigma = TaoMid(sigma - nu, sigmanew, sigma + nu);
456:           }
457:           break;

459:         default:
460:           sigma = sigmanew;
461:           break;
462:         }
463:         break;

465:       case Scale_Broyden:
466:         // Original version
467:         // Combine DFP and BFGS

469:         // This code appears to be numerically unstable.  We use the
470:         // original version because this was used to generate all of
471:         // the data and because it may be the least unstable of the
472:         // bunch.

474:         // P = Q = inv(D);
475:         info = P->CopyFrom(D); CHKERRQ(info);
476:         info = P->Reciprocal(); CHKERRQ(info);
477:         info = Q->CopyFrom(P); CHKERRQ(info);

479:         // V = y*y
480:         info = V->PointwiseMultiply(Gprev, Gprev); CHKERRQ(info);

482:         // W = inv(D)*s
483:         info = W->PointwiseMultiply(Xprev, P); CHKERRQ(info);
484:         info = W->Dot(Xprev, &sDs); CHKERRQ(info);

486:         // Safeguard rhotemp and sDs
487:         if (0.0 == rhotemp) {
488:           rhotemp = TAO_ZER_SAFEGUARD;
489:         }

491:         if (0.0 == sDs) {
492:           sDs = TAO_ZER_SAFEGUARD;
493:         }

495:         if (1.0 != phi) {
496:           // BFGS portion of the update
497:           // U = (inv(D)*s)*(inv(D)*s)
498:           info = U->PointwiseMultiply(W, W); CHKERRQ(info);

500:           // Assemble
501:           info = P->Axpy(1.0 / rhotemp, V); CHKERRQ(info);
502:           info = P->Axpy(-1.0 / sDs, U); CHKERRQ(info);
503:         }

505:         if (0.0 != phi) {
506:           // DFP portion of the update
507:           // U = inv(D)*s*y
508:           info = U->PointwiseMultiply(W, Gprev); CHKERRQ(info);

510:           // Assemble
511:           info = Q->Axpy(1.0 / rhotemp + sDs / (rhotemp * rhotemp), V); CHKERRQ(info);
512:           info = Q->Axpy(-2.0 / rhotemp, U); CHKERRQ(info);
513:         }

515:         if (0.0 == phi) {
516:           info = U->CopyFrom(P); CHKERRQ(info);
517:         }
518:         else if (1.0 == phi) {
519:           info = U->CopyFrom(Q); CHKERRQ(info);
520:         }
521:         else {
522:           // Broyden update
523:           info = U->Waxpby(1.0 - phi, P, phi, Q); CHKERRQ(info);
524:         }

526:         // Obtain inverse and ensure positive definite
527:         info = U->Reciprocal(); CHKERRQ(info);
528:         info = U->AbsoluteValue(); CHKERRQ(info);

530:         // Checking the diagonal scaling for not a number and infinity
531:         // should not be necessary for the Broyden update
532:         // info = U->Dot(U, &sDs); CHKERRQ(info);
533:         // if (sDs != sDs) {
534:         //   // not a number
535:         //   info = U->SetToConstant(TAO_ZER_SAFEGUARD); CHKERRQ(info);
536:         // }
537:         // else if ((sDs - sDs) != 0.0)) {
538:         //   // infinity
539:         //   info = U->SetToConstant(TAO_INF_SAFEGUARD); CHKERRQ(info);
540:         // }

542:         switch(rScaleType) {
543:         case Rescale_None:
544:           break;

546:         case Rescale_Scalar:
547:         case Rescale_GL:
548:           if (rScaleType == Rescale_GL) {
549:             // Gilbert and Lemarachal use the old diagonal
550:             info = P->CopyFrom(D); CHKERRQ(info);
551:           }
552:           else {
553:             // The default version uses the current diagonal
554:             info = P->CopyFrom(U); CHKERRQ(info);
555:           }

557:           // Compute s^T s 
558:           info = Xprev->Dot(Xprev, &s0temp); CHKERRQ(info);

560:           // Save information for special cases of scalar rescaling
561:           yy_rhistory[(updates - 1) % rescale_history] = y0temp;
562:           ys_rhistory[(updates - 1) % rescale_history] = rhotemp;
563:           ss_rhistory[(updates - 1) % rescale_history] = s0temp;

565:           if (0.5 == r_beta) {
566:             if (1 == TaoMin(updates, rescale_history)) {
567:               info = V->PointwiseMultiply(Y[0], P); CHKERRQ(info);
568:               info = V->Dot(Y[0], &yy_sum); CHKERRQ(info);

570:               info = W->PointwiseDivide(S[0], P); CHKERRQ(info);
571:               info = W->Dot(S[0], &ss_sum); CHKERRQ(info);

573:               ys_sum = ys_rhistory[0];
574:             }
575:             else {
576:               info = Q->CopyFrom(P); CHKERRQ(info);
577:               info = Q->Reciprocal(); CHKERRQ(info);

579:               // Compute summations for scalar scaling
580:               yy_sum = 0;        // No safeguard required
581:               ys_sum = 0;        // No safeguard required
582:               ss_sum = 0;        // No safeguard required
583:               for (i = 0; i < TaoMin(updates, rescale_history); ++i) {
584:                 info = V->PointwiseMultiply(Y[i], P); CHKERRQ(info);
585:                 info = V->Dot(Y[i], &yDy); CHKERRQ(info);
586:                 yy_sum += yDy;

588:                 info = W->PointwiseMultiply(S[i], Q); CHKERRQ(info);
589:                 info = W->Dot(S[i], &sDs); CHKERRQ(info);
590:                 ss_sum += sDs;

592:                 ys_sum += ys_rhistory[i];
593:               }
594:             }
595:           }
596:           else if (0.0 == r_beta) {
597:             if (1 == TaoMin(updates, rescale_history)) {
598:               // Compute summations for scalar scaling
599:               info = W->PointwiseDivide(S[0], P); CHKERRQ(info);
600:   
601:               info = W->Dot(Y[0], &ys_sum); CHKERRQ(info);
602:               info = W->Dot(W, &ss_sum); CHKERRQ(info);
603:   
604:               yy_sum += yy_rhistory[0];
605:             }
606:             else {
607:               info = Q->CopyFrom(P); CHKERRQ(info);
608:               info = Q->Reciprocal(); CHKERRQ(info);

610:               // Compute summations for scalar scaling
611:               yy_sum = 0;        // No safeguard required
612:               ys_sum = 0;        // No safeguard required
613:               ss_sum = 0;        // No safeguard required
614:               for (i = 0; i < TaoMin(updates, rescale_history); ++i) {
615:                 info = W->PointwiseMultiply(S[i], Q); CHKERRQ(info);
616:                 info = W->Dot(Y[i], &yDs); CHKERRQ(info);
617:                 ys_sum += yDs;

619:                 info = W->Dot(W, &sDs); CHKERRQ(info);
620:                 ss_sum += sDs;
621:   
622:                 yy_sum += yy_rhistory[i];
623:               }
624:             }
625:           }
626:           else if (1.0 == r_beta) {
627:             // Compute summations for scalar scaling
628:             yy_sum = 0;        // No safeguard required
629:             ys_sum = 0;        // No safeguard required
630:             ss_sum = 0;        // No safeguard required
631:             for (i = 0; i < TaoMin(updates, rescale_history); ++i) {
632:               info = V->PointwiseMultiply(Y[i], P); CHKERRQ(info);
633:               info = V->Dot(S[i], &yDs); CHKERRQ(info);
634:               ys_sum += yDs;

636:               info = V->Dot(V, &yDy); CHKERRQ(info);
637:               yy_sum += yDy;

639:               ss_sum += ss_rhistory[i];
640:             }
641:           }
642:           else {
643:             info = Q->CopyFrom(P); CHKERRQ(info);

645:             info = P->Pow(r_beta); CHKERRQ(info);
646:             info = Q->PointwiseDivide(P, Q); CHKERRQ(info);

648:             // Compute summations for scalar scaling
649:             yy_sum = 0;        // No safeguard required
650:             ys_sum = 0;        // No safeguard required
651:             ss_sum = 0;        // No safeguard required
652:             for (i = 0; i < TaoMin(updates, rescale_history); ++i) {
653:               info = V->PointwiseMultiply(P, Y[i]); CHKERRQ(info);
654:               info = W->PointwiseMultiply(Q, S[i]); CHKERRQ(info);

656:               info = V->Dot(V, &yDy); CHKERRQ(info);
657:               info = V->Dot(W, &yDs); CHKERRQ(info);
658:               info = W->Dot(W, &sDs); CHKERRQ(info);

660:               yy_sum += yDy;
661:               ys_sum += yDs;
662:               ss_sum += sDs;
663:             }
664:           }

666:           if (0.0 == r_alpha) {
667:             // Safeguard ys_sum 
668:             if (0.0 == ys_sum) {
669:               ys_sum = TAO_ZER_SAFEGUARD;
670:             }

672:             sigmanew = ss_sum / ys_sum;
673:           }
674:           else if (1.0 == r_alpha) {
675:             // Safeguard yy_sum 
676:             if (0.0 == yy_sum) {
677:               ys_sum = TAO_ZER_SAFEGUARD;
678:             }

680:             sigmanew = ys_sum / yy_sum;
681:           }
682:           else {
683:             denom = 2*r_alpha*yy_sum;

685:             // Safeguard denom
686:             if (0.0 == denom) {
687:               denom = TAO_ZER_SAFEGUARD;
688:             }

690:             sigmanew = ((2*r_alpha-1)*ys_sum +
691:                         sqrt((2*r_alpha-1)*(2*r_alpha-1)*ys_sum*ys_sum -
692:                              4*r_alpha*(r_alpha-1)*yy_sum*ss_sum)) / denom;
693:           }

695:           // If Q has small values, then Q^(r_beta - 1)
696:           // can have very large values.  Hence, ys_sum
697:           // and ss_sum can be infinity.  In this case,
698:           // sigmanew can either be not-a-number or infinity.

700:           if (TaoInfOrNaN(sigmanew)) {
701:             // sigmanew is not-a-number; skip rescaling
702:           }
703:           else if (!sigmanew) {
704:             // sigmanew is zero; this is a bad case; skip rescaling
705:           }
706:           else {
707:             // sigmanew is positive
708:             info = U->Scale(sigmanew); CHKERRQ(info);
709:           }
710:           break;
711:         }

713:         // Modify for previous information
714:         switch(limitType) {
715:         case Limit_Average:
716:           if (1.0 == mu) {
717:             info = D->CopyFrom(U); CHKERRQ(info);
718:           }
719:           else if (mu) {
720:             info = D->Axpby(mu, U, 1.0 - mu); CHKERRQ(info);
721:           }
722:           break;
723:  
724:         case Limit_Relative:
725:           if (mu) {
726:             info = P->ScaleCopyFrom(1.0 - mu, D); CHKERRQ(info);
727:             info = Q->ScaleCopyFrom(1.0 + mu, D); CHKERRQ(info);
728:             info = D->Median(P, U, Q); CHKERRQ(info);
729:           }
730:           break;

732:         case Limit_Absolute:
733:           if (nu) {
734:             info = P->CopyFrom(D); CHKERRQ(info);
735:             info = P->AddConstant(-nu); CHKERRQ(info);
736:                info = Q->CopyFrom(D); CHKERRQ(info);
737:             info = Q->AddConstant(nu); CHKERRQ(info);
738:             info = D->Median(P, U, Q); CHKERRQ(info);
739:           }
740:           break;

742:         default:
743:           info = D->CopyFrom(U); CHKERRQ(info);
744:           break;
745:         } 
746:         break;
747:       }

749:       Xprev = S[lm]; Gprev = Y[lm];
750:     } 
751:     else { 
752:       ++rejects;
753:     }
754:   }
755:   
756:   ++iter;
757:   info = Xprev->CopyFrom(x); CHKERRQ(info);
758:   info = Gprev->CopyFrom(g); CHKERRQ(info);
759:   TaoFunctionReturn(0);
760: }

764: int TaoLMVMMat::View()
765: {
766:   TaoFunctionBegin;
767:   TaoFunctionReturn(0);
768: }

772: int TaoLMVMMat::SetDelta(double d)
773: {
774:   TaoFunctionBegin;
775:   delta = TaoAbsDouble(d);
776:   delta = TaoMax(delta_min, delta);
777:   delta = TaoMin(delta_max, delta);
778:   TaoFunctionReturn(0);
779: }

783: int TaoLMVMMat::SetScale(TaoVec *s)
784: {
785:   TaoFunctionBegin;
786:   scale = s;
787:   TaoFunctionReturn(0);
788: }

792: int TaoLMVMMat::SetH0(TaoMat *HH0)
793: {
794:   TaoFunctionBegin;
795:   if (H0) { 
796:     H0 = 0;
797:   }
798:   H0default = TAO_TRUE;

800:   if (HH0) {
801:     H0 = HH0;
802:     H0default = TAO_FALSE;
803:   }
804:   TaoFunctionReturn(0);
805: }

809: int TaoLMVMMat::GetX0(TaoVec *x)
810: {
811:   int i,info;

813:   TaoFunctionBegin;
814:   info = x->CopyFrom(Xprev); CHKERRQ(info);
815:   for (i = 0; i < lmnow; ++i) {
816:     info = x->Axpy(-1.0, S[i]); CHKERRQ(info);
817:   }
818:   TaoFunctionReturn(0);
819: }

823: int TaoLMVMMat::InitialApproximation(TaoVec *x)
824: {
825:   TaoFunctionBegin;
826:   TaoFunctionReturn(0);
827: }