Xyce  6.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
N_NLS_DampedNewton.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------------
2 // Copyright Notice
3 //
4 // Copyright 2002 Sandia Corporation. Under the terms
5 // of Contract DE-AC04-94AL85000 with Sandia Corporation, the U.S.
6 // Government retains certain rights in this software.
7 //
8 // Xyce(TM) Parallel Electrical Simulator
9 // Copyright (C) 2002-2014 Sandia Corporation
10 //
11 // This program is free software: you can redistribute it and/or modify
12 // it under the terms of the GNU General Public License as published by
13 // the Free Software Foundation, either version 3 of the License, or
14 // (at your option) any later version.
15 //
16 // This program is distributed in the hope that it will be useful,
17 // but WITHOUT ANY WARRANTY; without even the implied warranty of
18 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 // GNU General Public License for more details.
20 //
21 // You should have received a copy of the GNU General Public License
22 // along with this program. If not, see <http://www.gnu.org/licenses/>.
23 //-----------------------------------------------------------------------------
24 
25 //-------------------------------------------------------------------------
26 // Filename : $RCSfile: N_NLS_DampedNewton.C,v $
27 //
28 // Purpose : Body file for the implemenation of the Newton trust-region
29 // related methods.
30 //
31 // Special Notes :
32 //
33 // Creator : Eric R. Keiter, SNL, Parallel Computational Sciences
34 //
35 // Creation Date : 05/28/00
36 //
37 // Revision Information:
38 // ---------------------
39 //
40 // Revision Number: $Revision: 1.226 $
41 //
42 // Revision Date : $Date: 2014/08/07 23:08:54 $
43 //
44 // Current Owner : $Author: dgbaur $
45 //-------------------------------------------------------------------------
46 
47 #include <Xyce_config.h>
48 
49 
50 // ---------- Standard Includes ----------
51 
52 #include <N_UTL_Misc.h>
53 
54 #ifdef HAVE_CSTDIO
55 #include <cstdio>
56 #else
57 #include <stdio.h>
58 #endif
59 
60 #include <vector>
61 
62 // ---------- Other Includes ----------
63 
64 // ---------- Xyce Includes ----------
65 
66 #include <N_NLS_ConstraintBT.h>
67 #include <N_NLS_DampedNewton.h>
68 
69 #include <N_NLS_TwoLevelNewton.h>
70 #include <N_NLS_ParamMgr.h>
71 
72 #include <N_LOA_Loader.h>
73 
74 #include <N_ANP_AnalysisManager.h>
75 
76 #include <N_LAS_System.h>
77 #include <N_LAS_Vector.h>
78 #include <N_LAS_Matrix.h>
79 #include <N_LAS_Solver.h>
80 #include <N_LAS_Problem.h>
81 #include <N_LAS_Builder.h>
82 #include <N_NLS_NonLinearSolver.h>
83 #include <N_ERH_ErrorMgr.h>
84 #include <N_UTL_Param.h>
85 #include <N_UTL_OptionBlock.h>
86 #include <N_UTL_Xyce.h>
87 
88 #include <N_IO_CmdParse.h>
89 
90 // ---------- Static Initializations ----------
91 
92 namespace Xyce {
93 namespace Nonlinear {
94 
95 //-----------------------------------------------------------------------------
96 // Function : DampedNewton::DampedNewton
97 // Purpose : constructor
98 // Special Notes :
99 // Scope : public
100 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
101 // Creation Date : 5/01/00
102 //-----------------------------------------------------------------------------
103 DampedNewton::DampedNewton(N_IO_CmdParse & cp):
104  NonLinearSolver (cp),
105  nlParams(DC_OP,cp),
106  iNumCalls_(0),
107  loadJacobianFlag_(true),
108  nlConstraintPtr_(0),
109  searchDirectionPtr_(0),
110  tmpVectorPtr_(0),
111  delta_(pow(N_UTL_MachineDependentParams::MachineEpsilon(), 0.67)),
112  // set the measures:
113  normRHS_(0.0),
114  maxNormRHS_(0.0),
115  maxNormRHSindex_(-1),
116  normRHS_init_(0.0),
117  normSoln_(0.0),
118  normDX_(0.0),
119  wtNormDX_(0.0),
120  stepLength_(1.0),
121  constraintFactor_(1.0),
122  nlStep_(0),
123  newtonStep_(0),
124  modNewtonStep_(0),
125  descentStep_(0),
126  searchStep_(0),
127  firstTime(true),
128  initialDeltaXTol(0.0),
129  etaOld(0.1),
130  nlResNormOld(0.0),
131  tmpConvRate(0.0),
132  count (0)
133 {
136 
138 }
139 
140 //-----------------------------------------------------------------------------
141 // Function : DampedNewton::~DampedNewton
142 // Purpose : destructor
143 // Special Notes :
144 // Scope : public
145 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
146 // Creation Date : 5/01/00
147 //-----------------------------------------------------------------------------
149 {
150  delete tmpVectorPtr_;
151  tmpVectorPtr_ = 0;
152 
153  delete nlConstraintPtr_;
154  nlConstraintPtr_ = 0;
155 
156  delete searchDirectionPtr_;
158 
159  delete nlpMgrPtr_;
160  nlpMgrPtr_ = 0;
161 
162 #ifdef Xyce_DEBUG_VOLTLIM
163  delete jdxVLVectorPtr_;
164  jdxVLVectorPtr_ = 0;
165  delete fdxVLVectorPtr_;
166  fdxVLVectorPtr_ = 0;
167  delete qdxVLVectorPtr_;
168  qdxVLVectorPtr_ = 0;
169 #endif
170 
171 }
172 
173 //-----------------------------------------------------------------------------
174 // Function : DampedNewton::setOptions
175 // Purpose : Sets the nonlinear solver options.
176 // Special Notes :
177 // Scope : public
178 // Creator : Robert Hoekstra, SNL, Parallel Computational Sciences
179 // Creation Date : 9/29/00
180 //-----------------------------------------------------------------------------
181 bool DampedNewton::setOptions(const N_UTL_OptionBlock & OB)
182 {
183  bool bsuccess = nlParams.setOptions(OB);
184 
187 
188  return bsuccess;
189 }
190 
191 //-----------------------------------------------------------------------------
192 // Function : DampedNewton::setTranOptions
193 // Purpose : Sets the nonlinear solver options.
194 // Special Notes :
195 // Scope : public
196 // Creator : Eric R. Keiter, SNL, Parallel Computational Sciences
197 // Creation Date : 9/05/01
198 //-----------------------------------------------------------------------------
199 bool DampedNewton::setTranOptions(const N_UTL_OptionBlock & OB)
200 {
201  NLParams nlTranParams(TRANSIENT, commandLine_);
202  bool bsuccess = nlTranParams.setOptions(OB);
203 
204  nlpMgrPtr_->addParameterSet(TRANSIENT, nlTranParams);
205  return bsuccess;
206 }
207 
208 //-----------------------------------------------------------------------------
209 // Function : DampedNewton::setHBOptions
210 // Purpose : Sets the nonlinear solver options.
211 // Special Notes :
212 // Scope : public
213 // Creator : Ting Mei, SNL
214 // Creation Date : 02/03/2009
215 //-----------------------------------------------------------------------------
216 bool DampedNewton::setHBOptions(const N_UTL_OptionBlock & OB)
217 {
218  NLParams nlHBParams(HB_MODE, commandLine_);
219  bool bsuccess = nlHBParams.setOptions(OB);
220 
221  nlpMgrPtr_->addParameterSet(HB_MODE, nlHBParams);
222  return bsuccess;
223 }
224 
225 //-----------------------------------------------------------------------------
226 // Function : DampedNewton::initializeAll
227 // Purpose : Called after all register and set functions.
228 // Once the various registrations have taken place,
229 // this function sets the remaining pointers.
230 // Special Notes :
231 // Scope : public
232 // Creator : Tamara G. Kolda, SNL, Compuational Sciences and
233 // Creation Date : 1/31/02
234 //-----------------------------------------------------------------------------
236 {
237  bool bsuccess = NonLinearSolver::initializeAll();
238 
239  tmpVectorPtr_ = lasSysPtr_->builder().createVector();
240 
241  // make sure the current nlParams is correct.
243 
244  // create and initialize constraint backtracking vectors:
246 
247  searchDirectionPtr_ = lasSysPtr_->builder().createVector();
248 
249 #ifdef Xyce_DEBUG_VOLTLIM
250  // get the dx vector (due to voltage limiting)
251  dxVoltlimVectorPtr_ = lasSysPtr_->getDxVoltlimVector();
252 
253  // this vector is for debug purposes:
254  jdxVLVectorPtr_ = lasSysPtr_->builder().createVector();
255  fdxVLVectorPtr_ = lasSysPtr_->builder().createVector();
256  qdxVLVectorPtr_ = lasSysPtr_->builder().createVector();
257 
258  // get the test matrices
259  jacTestMatrixPtr_ = lasSysPtr_->getJacTestMatrix(); // old-DAE
260  dFdxTestMatrixPtr_ = lasSysPtr_->getdFdxTestMatrix(); // new-DAE
261  dQdxTestMatrixPtr_ = lasSysPtr_->getdQdxTestMatrix(); // new-DAE
262 #endif
263 
264  bsuccess = bsuccess && (jacobianMatrixPtr_ != 0);
265 
266  return bsuccess;
267 }
268 
269 //-----------------------------------------------------------------------------
270 // Function : DampedNewton::printHeader_
271 // Purpose : Print out header for step information.
272 // Special Notes :
273 // Scope : private
274 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
275 // Creation Date : 5/03/00
276 //-----------------------------------------------------------------------------
277 void DampedNewton::printHeader_(std::ostream &os)
278 {
279  os << std::endl
280  << " Iter Step Wt DX Inf-Norm 2-Norm (rel)\n"
281  << " -------------------------------------------------------------------\n";
282 }
283 
284 //-----------------------------------------------------------------------------
285 // Function : DampedNewton::printFooter_
286 // Purpose :
287 // Special Notes :
288 // Scope : private
289 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
290 // Creation Date : 5/03/11
291 //-----------------------------------------------------------------------------
292 void DampedNewton::printFooter_(std::ostream &os)
293 {
294  os << Xyce::section_divider << std::endl;
295 }
296 
297 //-----------------------------------------------------------------------------
298 // Function : DampedNewton::printStepInfo_
299 // Purpose : Print out Newton step information
300 // Special Notes :
301 // Scope : public
302 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
303 // Creation Date : 5/03/00
304 //-----------------------------------------------------------------------------
305 void DampedNewton::printStepInfo_(std::ostream &os, int step)
306 {
307  os << "Niter: " << std::fixed << step
308  << " " << std::setprecision(4) << std::scientific << stepLength_
309  << " " << std::setprecision(4) << std::scientific << wtNormDX_
310  << " " << std::setprecision(4) << std::scientific << maxNormRHS_
311  << " " << std::setprecision(4) << std::scientific << normRHS_rel_
312  << std::endl;
313 }
314 
315 //-----------------------------------------------------------------------------
316 // Function : DampedNewton::updateWeights_
317 // Purpose : Updates the error weight vector using the previous solution
318 // values as "typical" and the absolute and relative
319 // tolerances.
320 // Special Notes :
321 // Scope : private
322 // Creator : Scott A. Hutchinson, SNL, Computational Sciences
323 // Creation Date : 01/17/01
324 //-----------------------------------------------------------------------------
326 {
327  // ***** Weighting Vector for Solution *****
328 
329  // On the first call to the nonlinear solver,
330  // weigh based on the minimums.
331  double tmp = 0.0;
332  (*nextSolVectorPtrPtr_)->maxValue(&tmp);
333  double solnNorm = fabs(tmp);
334 
335  if ((iNumCalls_ == 0) && (solnNorm <= N_UTL_MachineDependentParams::DoubleMin()))
336  {
338  }
339  else
340  {
341  double newSoln, oldSoln;
342 
343  int length = (*nextSolVectorPtrPtr_)->localLength();
344  for (int i = 0; i < length; ++i)
345  {
346  newSoln = (*(*nextSolVectorPtrPtr_))[i];
347  oldSoln = (*(*currSolVectorPtrPtr_))[i];
348 
349  (*(solWtVectorPtr_))[i] =
350  nlParams.getRelTol() * Xycemax(fabs(oldSoln), fabs(newSoln)) +
352  }
353  }
354 }
355 
356 //-----------------------------------------------------------------------------
357 // Function : DampedNewton::solve
358 // Purpose : solve implements a damped Newton method nonlinear equations
359 // solver.
360 // Special Notes : The implementation is based upon the original damped Newton
361 // method implemented in Xyce(TM) by Eric Keiter.
362 //
363 // Selected variables and methods:
364 //
365 // nextSolVector the approximate solution to the nonlinear equation
366 // The rhsVector and jacobianMatrix load nextSolVector.
367 //
368 // equateTmpVectors_ tmpSolVector := nextSolVector
369 // updateTmpSol_ tmpSolVector := nextVector + searchDirection*stepLength_
370 //
371 // switchTmpVectors swap &tmpSolVector and &nextSolVector, makes it possible
372 // for the nonlinear solver to preserve the entire
373 // circuit state associated with nextSolVector.
374 //
375 // Scope : public
376 // Creator : Scott A. Hutchinson, SNL, Computational Sciences
377 // Creation Date : 01/09/01
378 //-----------------------------------------------------------------------------
380 {
381  static const char *trace = "DampedNewton::solve: ";
382 
383 #ifdef Xyce_DEBUG_NONLINEAR
384  setDebugFlags ();
385 
386  if (debugTimeFlag_ && nlParams.getDebugLevel() > 0 ) Xyce::dout() << retCodes_;
387 #endif
388 
390 
391  // Change the nlparams so that they are appropriate for the current
392  // time integration mode, if neccessary.
394 
395  // Get the current analysis mode
397 
398  if (mode1 == TRANSIENT)
399  {
400  // Set the tolerance
401  lasSolverPtr_->setTolerance(1.0e-09);
402  }
403 
404 #ifndef Xyce_SPICE_NORMS
405  if (firstTime)
406  {
407  // First time through we want to tighten up the convergence tolerance a la
408  // Petzold, et al. The DeltaXTol parameter is used by converge_().
411  firstTime = false;
412  }
413 #endif
414 
415 #ifdef Xyce_VERBOSE_NONLINEAR
416  // Output the nonlinear solver information header:
417  printHeader_(Xyce::lout());
418 #endif
419 
420 #ifdef Xyce_VERBOSE_NONLINEAR
421  // Print warning about using the gradient-only strategy
423  Xyce::Report::UserWarning0() << "Use of the gradient strategy (1) is discouraged";
424 #endif
425 
426  // For the initial RHS load, the step number needs to be zero. The Xyce
427  // device package needs to know this - either in the event that it might want
428  // to set initial conditions.
430 
431  // Recall that prior to this solver being called, the new solution has been
432  // predicted and resides in nextSolVectorPtr.
433  rhs_();
434 
435 #ifdef Xyce_VERBOSE_NONLINEAR
436  // Set some values and dummy values so that we can call printStepInfo_().
437 
438  // Relative RHS norm (normRHS_rel_). On the first evaluation this is one
439  // since normRHS_init = normRHS.
440  normRHS_rel_ = 1.0;
441 
442  // Weighted norm of the change (wtNormDX_) is set to zero initially.
443  wtNormDX_ = 0.0;
444 
445  // Max RHS norm (maxNormRHS_).
446  rhsVectorPtr_->infNorm(&maxNormRHS_);
447  std::vector<int> index(1, -1);
448  rhsVectorPtr_->infNormIndex( &index[0] );
449  maxNormRHSindex_ = index[0];
450 
451  // Print out the starting point information.
452  printStepInfo_(Xyce::lout(), nlStep_);
453 #endif
454 
455  // Check to ensure that we haven't been passed a converged solution already
456  // (i.e., normRHS < absTol)
457 
458  if (normRHS_ < N_UTL_MachineDependentParams::MachineEpsilon())
459  {
460  return true;
461  }
462 
463  // Used to calculate resConvRate which is used by converged_(). This changes
464  // each nonlinear iteration to be the norm of the RHS from the previous
465  // iteration.
466  double normRHS_old = normRHS_;
467 
468  // Used to calculate normRHS_rel_, which is output by printStepInfo_().
470 
471  // Set the flag that controls the loading and solution of the Jacobian
472  // matrix.
473  if (nlParams.getLinearOpt())
474  loadJacobianFlag_ = !(loaderPtr_->getLinearSystemFlag());
475  else
476  loadJacobianFlag_ = true;
477 
478  // If we're doing a modified-Newton iteration, initialize the solver flag.
480  {
481  // Reset these flags.
482 // lasSolverPtr_->setUpdatedJacobian(true);
483 // directSolverPtr_->setUpdatedJacobian(true);
484  }
485 
486  // Update the error weighting vector for use with the weighted norms.
487  if (mode1 == TRANSIENT)
488  updateWeights_();
489 
490  // Nonlinear loop. Loop continues until convergedStatus is nonzero. A
491  // positive value indicates that the method has converged. A negative value
492  // that the method has failed. This is updated using the subroutine
493  // converged_() at the end of each nonlinear solver iteration.
494  int convergedStatus = 0;
495  while (convergedStatus == 0)
496  {
497  // Increment step counters. This needs to be updated *before* the rhs_
498  // call, which is inside of computeStepLength_.
499  nlStep_++;
500 
501  // Calculate the Jacobian for the current iterate.
502  if (loadJacobianFlag_)
503  jacobian_();
504 
505 #ifdef Xyce_DEBUG_NONLINEAR
506  debugOutput1( *(lasSysPtr_->getJacobianMatrix()), *(lasSysPtr_->getRHSVector()));
507 #endif
508 
509  // Calculate a direction.
510  // Note: It is illegal to modify tmpVector between direction_ and
511  // computeStepLength_
512  direction_();
513 
514  // Perform update constraining...
517 
518  setX0_();
519 
520  // Calculate a step length (damping or backtracking for Newton, line search
521  // for others...) and take that step and recalulate the RHS.
523 
524 #ifdef Xyce_DEBUG_NONLINEAR
525  debugOutput3 ((**nextSolVectorPtrPtr_), *searchDirectionPtr_ );
526 #endif
527 
528  // Test for convergence based on the weighted norms of the RHS and DX.
529  // Also update convergence rate and counters.
530 
531  // Update the error weighting vector for use with the weighted norms.
532  if (mode1 != TRANSIENT)
533  updateWeights_();
534 
535  // Scale the search direction so that we can calculate some norms on it and
536  // use them in the convergence tests and output.
538 
539  // Weighted norm of the change (wtNormDX_) is used in converged_() and
540  // output by printStepInfo_().
541 
542 #ifdef Xyce_SPICE_NORMS
543  // TGKOLDA: Why is tmpVectorPtr_ used here and where was it set??
545  *tmpVectorPtr_,
546  &wtNormDX_);
547 #else
549 #endif
550 
551  // Relative RHS norm (normRHS_rel_) is output by printStepInfo_().
553 
554  // Max RHS norm (maxNormRHS_) is used in converged_() and output by
555  // printStepInfo_().
556  rhsVectorPtr_->infNorm(&maxNormRHS_);
557 
558  // Residual convergence rate (resConvRate) used by converged_(),
559  // also printed out below.
560  resConvRate_ = normRHS_ / normRHS_old;
561 
562  // Reset old norm for next nonlinear iteration.
563  normRHS_old = normRHS_;
564 
565 #ifdef Xyce_VERBOSE_NONLINEAR
566  printStepInfo_(Xyce::lout(), nlStep_);
567 #endif
568 
569  // Increment diagnostic step counters. These need to be updated after
570  // direction_, in which the current direction (NEWTON_DIR,
571  // MOD_NEWTON_DIR or GRADIENT_DIR) is set.
573  newtonStep_++;
574  else if (nlParams.getDirection() == MOD_NEWTON_DIR)
575  modNewtonStep_++;
576  else if (nlParams.getDirection() == GRADIENT_DIR)
577  descentStep_++;
578 
579 // #ifdef Xyce_VERBOSE_NONLINEAR
580 // Xyce::lout() << "Residual Convergence Rate: " << resConvRate_ << "\n" << std::endl;
581 // #endif
582 
583  // Evaluate the need for a fresh Jacobian and/or preconditioner...
585  evalModNewton_();
586 
587  // Check our iteration status. Returns positive if done, negative if error,
588  // zero otherwise.
589  convergedStatus = converged_();
590 
591  // If the norm is an NaN, then exit.
592  if (maxNormRHS_ != 0.0 &&
593  !(maxNormRHS_ < 0.0) &&
594  !(maxNormRHS_ > 0.0))
595  {
596  convergedStatus = retCodes_.nanFail; // default = -6
597  //break;
598  }
599 
600  } // while (convergedStatus == 0)
601 
602 #ifndef Xyce_SPICE_NORMS
603  // Reset the tolerance which is used in converged_().
605 #endif
606 
607  // Reset the linear solver tolerance.
608  lasSolverPtr_->setDefaultOption( "AZ_tol" );
609 
610  // Increment the number of calls to the nonlinear solver.
611  iNumCalls_++;
612 
613 #ifdef Xyce_VERBOSE_NONLINEAR
614  printFooter_(Xyce::lout());
615 #endif
616 
617  // A positive value of convergedStatus indicates a successful nonlinear
618  // solutions. Otherwise, there was an error.
619  return convergedStatus;
620 }
621 
622 //-----------------------------------------------------------------------------
623 // Function : DampedNewton::takeFirstSolveStep
624 //
625 // Purpose : This is the same as the function "solve", except that
626 // it only does the various initializations at the top
627 // of the "::solve" function, and only takes the first
628 // NL step.
629 //
630 // Special Notes :
631 // Scope : public
632 // Creator : Eric R. Keiter, SNL, Computational Sciences
633 // Creation Date : 01/09/01
634 //-----------------------------------------------------------------------------
636 {
637  static const char *trace = "DampedNewton::takeFirstSolveStep: ";
638 
639 #ifdef Xyce_DEBUG_NONLINEAR
640  setDebugFlags ();
641 #endif
642 
643  // resolve the counters issue later. It may be neccessary to save the
644  // counters for later. Or not.
646 
647  // Change the nlparams so that they are appropriate for the current
648  // time integration mode, if neccessary.
650 
652  lasSolverPtr_->setTolerance(1.0e-09);
653 
654 
655 #ifndef Xyce_SPICE_NORMS
656 #if 0
657  static bool firstTime = true;
658  static double initialDeltaXTol;
659 
660  if (firstTime)
661  {
662  // First time through we want to tighten up the convergence tolerance a la
663  // Petzold, et al. The DeltaXTol parameter is used by converge_().
664  initialDeltaXTol = nlParams.getDeltaXTol();
665  nlParams.setDeltaXTol(initialDeltaXTol * 0.01);
666  firstTime = false;
667  }
668 #endif // figure this out later. Mainly, understand what is really meant
669  // by "firstTime", and make sure it works correctly for this function.
670 #endif
671 
672 #ifdef Xyce_VERBOSE_NONLINEAR
673  // Output the nonlinear solver information header:
674  printHeader_(Xyce::lout());
675 
676 #endif
677 
678 #ifdef Xyce_VERBOSE_NONLINEAR
679  // Print warning about using the gradient-only strategy
681  Xyce::Report::UserWarning0() << "Use of the gradient strategy (1) is discouraged";
682 #endif
683 
684  // For the initial RHS load, the step number needs to be zero. The Xyce
685  // device package needs to know this - either in the event that it might want
686  // to set initial conditions.
688 
689  // Recall that prior to this solver being called, the new solution has been
690  // predicted and resides in nextSolVectorPtr.
691  rhs_();
692 
693 #ifdef Xyce_VERBOSE_NONLINEAR
694  // Set some values and dummy values so that we can call printStepInfo_().
695 
696  // Relative RHS norm (normRHS_rel_). On the first evaluation this is one
697  // since normRHS_init = normRHS.
698  normRHS_rel_ = 1.0;
699 
700  // Weighted norm of the change (wtNormDX_) is set to zero initially.
701  wtNormDX_ = 0.0;
702 
703  // Max RHS norm (maxNormRHS_).
704  rhsVectorPtr_->infNorm(&maxNormRHS_);
705 
706  // Print out the starting point information.
707  printStepInfo_(Xyce::lout(), nlStep_);
708 #endif
709 
710  // Check to ensure that we haven't been passed a converged solution already
711  // (i.e., normRHS < absTol)
712  if (normRHS_ < N_UTL_MachineDependentParams::MachineEpsilon()) return 1;
713 
714  // Used to calculate resConvRate which is used by converged_(). This changes
715  // each nonlinear iteration to be the norm of the RHS from the previous
716  // iteration.
717  double normRHS_old = normRHS_;
718 
719  // Used to calculate normRHS_rel_, which is output by printStepInfo_().
721 
722  // Set the flag that controls the loading and solution of the Jacobian
723  // matrix.
724  if (nlParams.getLinearOpt())
725  loadJacobianFlag_ = !(loaderPtr_->getLinearSystemFlag());
726  else
727  loadJacobianFlag_ = true;
728 
729  // The variable convergedStatus: positive value indicates that
730  // the method has converged. A negative value
731  // that the method has failed. This is updated using the subroutine
732  // converged_() at the end of each nonlinear solver iteration.
733  int convergedStatus = 0;
734 
735  // Increment step counters. This needs to be updated *before* the rhs_
736  // call, which is inside of computeStepLength_.
737  nlStep_++;
738 
739  // Calculate the Jacobian for the current iterate.
741 
742 #ifdef Xyce_DEBUG_NONLINEAR
743  debugOutput1( *(lasSysPtr_->getJacobianMatrix()), *(lasSysPtr_->getRHSVector()));
744 #endif
745 
746  // Calculate a direction.
747  // Note: It is illegal to modify tmpVector between direction_ and
748  // computeStepLength_
749  direction_();
750 
751  // Perform update constraining...
754 
755  setX0_();
756 
757  // Calculate a step length (damping or backtracking for Newton, line search
758  // for others...) and take that step and recalulate the RHS.
760 
761 #ifdef Xyce_DEBUG_NONLINEAR
762  debugOutput3 ((**nextSolVectorPtrPtr_), *searchDirectionPtr_ );
763 #endif
764 
765  // Test for convergence based on the weighted norms of the RHS and DX.
766  // Also update convergence rate and counters. 1. Update the error
767  // weighting vector for use with the weighted norms.
768  updateWeights_();
769 
770  // Scale the search direction so that we can calculate some norms on it and
771  // use them in the convergence tests and output.
773 
774  // Weighted norm of the change (wtNormDX_) is used in converged_() and
775  // output by printStepInfo_().
776 
777 #ifdef Xyce_SPICE_NORMS
778  // TGKOLDA: Why is tmpVectorPtr_ used here and where was it set??
780  *tmpVectorPtr_,
781  &wtNormDX_);
782 #else
784 #endif
785 
786  // Relative RHS norm (normRHS_rel_) is output by printStepInfo_().
788 
789  // Max RHS norm (maxNormRHS_) is used in converged_() and output by
790  // printStepInfo_().
791  rhsVectorPtr_->infNorm(&maxNormRHS_);
792 
793  // Residual convergence rate (resConvRate) used by converged_(),
794  // also printed out below.
795  resConvRate_ = normRHS_ / normRHS_old;
796 
797  // Reset old norm for next nonlinear iteration.
798  normRHS_old = normRHS_;
799 
800 #ifdef Xyce_VERBOSE_NONLINEAR
801  printStepInfo_(Xyce::lout(), nlStep_);
802 #endif
803 
804  // Increment diagnostic step counters. These need to be updated after
805  // direction_, in which the current direction (NEWTON_DIR,
806  // MOD_NEWTON_DIR or GRADIENT_DIR) is set.
808  newtonStep_++;
809  else if (nlParams.getDirection() == MOD_NEWTON_DIR)
810  modNewtonStep_++;
811  else if (nlParams.getDirection() == GRADIENT_DIR)
812  descentStep_++;
813 
814  // Evaluate the need for a fresh Jacobian and/or preconditioner...
816  evalModNewton_();
817 
818  // Check our iteration status. Returns positive if done, negative if error,
819  // zero otherwise.
820  convergedStatus = converged_();
821 
822 #ifndef Xyce_SPICE_NORMS
823 #if 0
824  // Reset the tolerance which is used in converged_().
825  nlParams.setDeltaXTol(initialDeltaXTol);
826 #endif
827 #endif
828 
829  // Reset the linear solver tolerance.
830  lasSolverPtr_->setDefaultOption( "AZ_tol" );
831 
832  // Increment the number of calls to the nonlinear solver.
833  iNumCalls_++;
834 
835  // A positive value of convergedStatus indicates a successful nonlinear
836  // solutions. Otherwise, there was an error.
837  return convergedStatus;
838 }
839 
840 //-----------------------------------------------------------------------------
841 // Function : DampedNewton::takeOneSolveStep
842 //
843 // Purpose : Just like "::initializeSolve" is essentially the top
844 // of the "::solve" function, this function is similar
845 // to the bottom of the "::solve" function. This
846 // function also only takes a single NL step.
847 //
848 // The main differences between this function and solve are:
849 // 1) There is no while loop. - this is only a single
850 // step.
851 //
852 // 2) Most tasks performed prior to the while loop are not
853 // performed here, as they are initialization tasks.
854 // These include:
855 
856 // - step counters, such as "nlStep_" and "newtonStep_"
857 // are not reset to zero.
858 //
859 // - The function "resetCountersAndTimers" is not called.
860 // This may be fixed up later.
861 //
862 // - The resetting of the deltaXtol, which is performed
863 // on the "first call" to ::solve, is not done.
864 // I'll have to sort this out later as well.
865 //
866 // - the loadJacobianFlag is not updated - assumed to
867 // have not changed.
868 //
869 // Special Notes :
870 // Scope : public
871 // Creator : Eric R. Keiter, SNL, Computational Sciences
872 // Creation Date : 01/09/01
873 //-----------------------------------------------------------------------------
875 {
876  static const char *trace = "DampedNewton::takeOneSolveStep";
877 
878 #ifdef Xyce_DEBUG_NONLINEAR
879  setDebugFlags ();
880 #endif
881 
882  // Change the nlparams so that they are appropriate for the current
883  // time integration mode, if neccessary.
885 
887  lasSolverPtr_->setTolerance(1.0e-09);
888 
889  // Recall that prior to this solver step being called, a new solution has been
890  // predicted and resides in nextSolVectorPtr.
891  rhs_();
892 
893  // Used to calculate resConvRate which is used by converged_(). This changes
894  // each nonlinear iteration to be the norm of the RHS from the previous
895  // iteration.
896  double normRHS_old = normRHS_;
897 
898  // Nonlinear loop. Loop continues until convergedStatus is nonzero. A
899  // positive value indicates that the method has converged. A negative value
900  // that the method has failed. This is updated using the subroutine
901  // converged_() at the end of each nonlinear solver iteration.
902  int convergedStatus = 0;
903 
904  // Increment step counters. This needs to be updated *before* the rhs_
905  // call, which is inside of computeStepLength_.
906  nlStep_++;
907 
908  // Calculate the Jacobian for the current iterate.
909  if (loadJacobianFlag_)
910  jacobian_();
911 
912 #ifdef Xyce_DEBUG_NONLINEAR
913  debugOutput1( *(lasSysPtr_->getJacobianMatrix()), *(lasSysPtr_->getRHSVector()));
914 #endif
915 
916  // Calculate a direction.
917  // Note: It is illegal to modify tmpVector between direction_ and
918  // computeStepLength_
919  direction_();
920 
921  // Perform update constraining...
924 
925  setX0_();
926 
927  // Calculate a step length (damping or backtracking for Newton, line search
928  // for others...) and take that step and recalulate the RHS.
930 
931 #ifdef Xyce_DEBUG_NONLINEAR
932  debugOutput3 ((**nextSolVectorPtrPtr_), *searchDirectionPtr_ );
933 #endif
934 
935  // Test for convergence based on the weighted norms of the RHS and DX.
936  // Also update convergence rate and counters. 1. Update the error
937  // weighting vector for use with the weighted norms.
938  updateWeights_();
939 
940  // Scale the search direction so that we can calculate some norms on it and
941  // use them in the convergence tests and output.
943 
944  // Weighted norm of the change (wtNormDX_) is used in converged_() and
945  // output by printStepInfo_().
946 
947 #ifdef Xyce_SPICE_NORMS
948  // TGKOLDA: Why is tmpVectorPtr_ used here and where was it set??
950  *tmpVectorPtr_,
951  &wtNormDX_);
952 #else
954 #endif
955 
956  // Relative RHS norm (normRHS_rel_) is output by printStepInfo_().
958 
959  // Max RHS norm (maxNormRHS_) is used in converged_() and output by
960  // printStepInfo_().
961  rhsVectorPtr_->infNorm(&maxNormRHS_);
962 
963  // Residual convergence rate (resConvRate) used by converged_(),
964  // also printed out below.
965  resConvRate_ = normRHS_ / normRHS_old;
966 
967  // Reset old norm for next nonlinear iteration.
968  normRHS_old = normRHS_;
969 
970 #ifdef Xyce_VERBOSE_NONLINEAR
971  printStepInfo_(Xyce::lout(), nlStep_);
972 #endif
973 
974  // Increment diagnostic step counters. These need to be updated after
975  // direction_, in which the current direction (NEWTON_DIR,
976  // MOD_NEWTON_DIR or GRADIENT_DIR) is set.
978  newtonStep_++;
979  else if (nlParams.getDirection() == MOD_NEWTON_DIR)
980  modNewtonStep_++;
981  else if (nlParams.getDirection() == GRADIENT_DIR)
982  descentStep_++;
983 
984  // Evaluate the need for a fresh Jacobian and/or preconditioner...
986  evalModNewton_();
987 
988  // Check our iteration status. Returns positive if done, negative if error,
989  // zero otherwise.
990  convergedStatus = converged_();
991 
992 #ifdef Xyce_DEBUG_NONLINEAR
993  if (convergedStatus > 0)
994  {
995  // ERK Note: this one needs the nl solve step incremented by 1!
996  nlStep_++;
997  debugOutput1( *(lasSysPtr_->getJacobianMatrix()), *(lasSysPtr_->getRHSVector()));
998  nlStep_--; // restore original value.
999  }
1000 #endif
1001 
1002  // Reset the linear solver tolerance.
1003  lasSolverPtr_->setDefaultOption( "AZ_tol" );
1004 
1005  // Increment the number of calls to the nonlinear solver.
1006  iNumCalls_++;
1007 
1008  // A positive value of convergedStatus indicates a successful nonlinear
1009  // solutions. Otherwise, there was an error.
1010  return convergedStatus;
1011 }
1012 
1013 //-----------------------------------------------------------------------------
1014 // Function : DampedNewton::updateX_()
1015 //
1016 // Purpose : Update the value of nextSolVectorPtr_ using
1017 // steplength and the searchDirectionPtr_. On input, it is
1018 // assumed that tmpSolVectorPtr_ contains the intial guess for
1019 // the nonlinear iteration, i.e., setX0_ was called earlier.
1020 //
1021 // Special Notes : Also modifies errorEstVectorPtr_. Only this and
1022 // updateX_() should modify the solution directly.
1023 //
1024 // ERK: This is new version of this function, hopefully less
1025 // confusing. The roles of "tmp" and "next" have been
1026 // reversed, which removes the need for some of the switch
1027 // functions.
1028 //
1029 // Scope : private
1030 // Creator : Tamara G. Kolda, SNL, Compuational Sciences and
1031 // Mathematics Research Department
1032 //
1033 // Eric Keiter, SNL, Parallel Computational Sciences (9233)
1034 //
1035 // Creation Date : 01/24/02
1036 //-----------------------------------------------------------------------------
1038 {
1039  (*nextSolVectorPtrPtr_)->daxpy(**tmpSolVectorPtrPtr_,
1040  stepLength_,
1042 }
1043 
1044 //-----------------------------------------------------------------------------
1045 // Function : DampedNewton::rhs_()
1046 // Purpose : Updates the RHS based on nextSolVectorPtrPtr_ and
1047 // calculates normRHS_ based on nlParams.getNormLevel().
1048 // Special Notes : The rhsVectorPtr_ is really the NEGATIVE of F(x).
1049 // Scope : private
1050 // Creator : Tamara G. Kolda, SNL, Compuational Sciences and
1051 // Mathematics Research Department
1052 // Creation Date : 06/19/01
1053 //-----------------------------------------------------------------------------
1055 {
1056  bool status = NonLinearSolver::rhs_();
1057 
1058 #ifdef Xyce_DEBUG_NONLINEAR
1059  debugOutput3 ((**nextSolVectorPtrPtr_), *searchDirectionPtr_ );
1060 #endif
1061 
1063 
1064  return status;
1065 }
1066 
1067 //-----------------------------------------------------------------------------
1068 // Function : DampedNewton::newtonDirection_
1069 // Purpose : Computes the Newton Direction inexactly using
1070 // jacobianMatrixPtr_ and the rhsVectorPtr_.
1071 // On output, NewtonVectorPtr_ contains the Newton
1072 // Direction.
1073 // Special Notes : This is *not* meant to override newton_ in
1074 // NonLinearSolver.
1075 // Scope : private
1076 // Creator : Tamara G. Kolda, SNL, Compuational Sciences and
1077 // Mathematics Research Department
1078 // Creation Date : 06/19/01
1079 //-----------------------------------------------------------------------------
1081 {
1082  // Set the solver tolerance based on the residual norm.
1084 
1085  return NonLinearSolver::newton_();
1086 }
1087 
1088 //-----------------------------------------------------------------------------
1089 // Function : DampedNewton::direction_
1090 // Purpose : The function calculates the direction vector used in
1091 // the nonlinear solver (e.g., Newton direction).
1092 // Special Notes : If the search direction is a Newton direction, then a linear
1093 // system is solved inexactly. To ensure that the approximate
1094 // solution is a descent direction, the projection of the
1095 // solution along the gradient is computed to working
1096 // precision. On output, searchDirectionPtr_ points to the
1097 // search direction vector. The vector allocated to store the
1098 // NewtonVector is always used to store the
1099 // SearchDirectionVector, regardless of the actual algorithm
1100 // used to determine the search direction.
1101 // Scope : private
1102 // Creator : Scott A. Hutchinson, SNL, Computational Sciences Department
1103 // Creation Date : 02/21/01
1104 //-----------------------------------------------------------------------------
1106 {
1107  static const char *trace = "DampedNewton::direction_: ";
1108 
1109 #ifdef Xyce_DEBUG_NONLINEAR
1110 #if 0
1111  if (debugTimeFlag_ && nlParams.getDebugLevel() > 0 )
1112  {
1113  Xyce::dout() << "Search Method: " << nlParams.getSearchMethod() << " Strategy: " << nlParams.getNLStrategy() << std::endl;
1114  // Print status messages
1115  switch (nlParams.getNLStrategy())
1116  {
1117  case NEWTON:
1118  Xyce::dout() << "nonlinear strategy is Newton" << std::endl;
1119  break;
1120  case MOD_NEWTON:
1121  Xyce::dout() << "nonlinear strategy is Modified-Newton" << std::endl;
1122  break;
1123  case NEWTON_GRADIENT:
1124  Xyce::dout() << "\tnonlinear strategy is Inexact Newton" << std::endl
1125  << "\t\tconvergence rate: " << resConvRate_ << std::endl;
1126  break;
1127  case MOD_NEWTON_GRADIENT:
1128  Xyce::dout() << "\tnonlinear strategy is Inexact Modified-Newton\n\t\tconvergence rate: " << resConvRate_ << std::endl;
1129  break;
1130  case GRADIENT:
1131  Xyce::dout() << "nonlinear strategy is Gradient" << std::endl;
1132  break;
1133 
1134  default:
1135  Xyce::Report::DevelFatal0().in(trace) << "Invalid nonlinear solver strategy: " << static_cast<int> (nlParams.getNLStrategy());
1136  break;
1137  } // end switch on nlParams.getNLStrategy() - status messages
1138  }
1139 #endif
1140 #endif
1141 
1142  // --- Set the Direction Type ---
1143  switch (nlParams.getNLStrategy())
1144  {
1145  case NEWTON:
1147  break;
1148 
1149  case MOD_NEWTON:
1151  break;
1152 
1153  case NEWTON_GRADIENT:
1154  {
1155  const double minRedFac = 0.999;
1156  if ((nlStep_ > 15) || (resConvRate_ > minRedFac))
1158  else
1160 
1161  break;
1162  }
1163  case MOD_NEWTON_GRADIENT:
1164  {
1165  const double minRedFac = 0.999;
1166  if ((nlStep_ > 15) || (resConvRate_ > minRedFac))
1167  {
1169  }
1170  else
1171  {
1173  }
1174  break;
1175  }
1176  case GRADIENT:
1178  break;
1179 
1180  default:
1182  break;
1183 
1184  } // end switch on nlParams.getNLStrategy()
1185 
1186  // --- Calculate SearchDirection ---
1187  switch (nlParams.getDirection())
1188  {
1189  case NEWTON_DIR:
1190 
1191  // Compute the Newton direction
1192  newton_();
1193 
1194  // Copy the Newton direction into the search direction
1196 
1197  break;
1198 
1199  case MOD_NEWTON_DIR:
1200 
1201  // Compute the modified-Newton direction
1202  newton_();
1203 
1204  // Copy the modified-Newton direction into the search direction
1206 
1207  break;
1208 
1209  case GRADIENT_DIR:
1210 
1211  // Compute the gradient
1212  gradient_();
1213 
1214  // Copy Gradient into Search Direction, then Reverse and scale
1215  // search direction. (Can this be done in one operation??)
1217  searchDirectionPtr_->scale(-1.0 / normRHS_);
1218 
1219  break;
1220 
1221  default:
1222 
1223 #ifdef Xyce_DEBUG_NONLINEAR
1224  Xyce::Report::DevelFatal0().in("DampedNewton::direction_") << "Invalid search direction: " << static_cast<int>(nlParams.getDirection());
1225 #endif
1226 
1227  break;
1228 
1229  } // end switch on nlParams.getDirection()
1230 
1231 } // end direction_
1232 
1233 //-----------------------------------------------------------------------------
1234 // Function : DampedNewton::computeStepLength_
1235 // Purpose : This function calculates the step length for a given search
1236 // direction. If the search direction is Newton, it may use one
1237 // of a variety of backtracking methods including on suggested
1238 // by
1239 // Special Notes : !!!!!NOTE: In this method, we use the term FULLSTEP to mean
1240 // the maximum allowed Newton step (0..1). fullStep = 1.
1241 // corresponds to a Newton step, but when computeStepLength_
1242 // is coupled with a constraint back tracking, 0<fullStep<=1.
1243 // Scope : private
1244 // Creator : Scott A. Hutchinson, SNL, Computational Sciences Department
1245 // Creation Date : 02/21/01
1246 //-----------------------------------------------------------------------------
1248 
1249 {
1250  static const char *trace = "DampedNewton::computeStepLength_";
1251 
1252  searchStep_ = 0; // number of search steps
1253 
1254  // Choose step length approach based on the current solution method.
1255  switch (nlParams.getSearchMethod())
1256  {
1257  case DIVIDE:
1258 
1259  return divide_();
1260  break;
1261 
1262  case BACKTRACK:
1263 
1264  return backtrack_();
1265  break;
1266 
1267  case SIMPLE_BACKTRACK:
1268 
1269  return simpleBacktrack_();
1270  break;
1271 
1272  case BANK_ROSE:
1273 
1274  return bankRose_();
1275  break;
1276 
1277  case DESCENT:
1278 
1279  return descent_();
1280  break;
1281 
1282  default:
1283 
1284  return fullNewton_();
1285  break;
1286  }
1287 
1288  // If the code gets to this point in the function there has been an
1289  // error...
1290  return false;
1291 }
1292 
1293 //-----------------------------------------------------------------------------
1294 // Function : DampedNewton::divide_()
1295 // Purpose : Eric's "Divide and Conquer" method. Subsequently halves
1296 // the steplength until there is improvement in the RHS or
1297 // the steplength becomes smaller than N_UTL_MachineDependentParams::MachineEpsilon().
1298 // Special Notes : Does not verify that the search direction is a descent
1299 // direction.
1300 // Scope : private
1301 // Creator : Tamara G. Kolda, SNL, Compuational Sciences and
1302 // Mathematics Research Department
1303 // Creation Date : 06/29/01
1304 //-----------------------------------------------------------------------------
1306 {
1307  static const char *trace = "DampedNewton::divide__";
1308 
1309 #if defined(Xyce_DEBUG_NONLINEAR)
1310  char tmpStr[128];
1311 #endif
1312 
1313  // Want to get reduction in the norm, so we first copy the norm from
1314  // the previous iterate
1315  const double oldNormRHS = normRHS_;
1316 
1317  // The upper bound on the backtracking (damping) is determined by
1318  // the applied constraints.
1319  const double fullStep = constraintFactor_;
1320 
1321  // Starting value
1322  stepLength_ = fullStep;
1323 
1324  // Update the solution.
1325  updateX_();
1326 
1327  // Evaluate the new residual and take its norm:
1328  rhs_();
1329 
1330  if (Xyce::DEBUG_NONLINEAR && debugTimeFlag_ && nlParams.getDebugLevel() > 0 )
1331  {
1332  Xyce::dout() << "\n\tSearch Step: " << searchStep_ << ", Step Size: " << stepLength_ << std::endl
1333  << "\toldNormRHS: " << oldNormRHS << ", normRHS: " << normRHS_ << std::endl;
1334  }
1335 
1336  // Linesearch...
1337  bool searchDone = (normRHS_ < oldNormRHS);
1338  while (!searchDone)
1339  {
1340  // Cut the step length
1341  stepLength_ *= 0.5;
1342 
1343  // If the step length gets too small, try full step and finish.
1344  if (stepLength_ < N_UTL_MachineDependentParams::MachineEpsilon())
1345  {
1346 #ifdef Xyce_DEBUG_NONLINEAR
1347  if (debugTimeFlag_ && nlParams.getDebugLevel() > 0 )
1348  {
1349  Xyce::Report::UserWarning0() << "\tStep size too small: " << stepLength_ << "\n\tTaking a full step.\n";
1350  }
1351 #endif
1352 
1353  stepLength_ = fullStep;
1354  searchDone = true;
1355  }
1356 
1357  updateX_();
1358 
1359  // Evaluate the new RHS vector and take norm:
1360  rhs_();
1361 
1362  searchStep_++;
1363 
1364  // Check to see if we're done:
1365  searchDone = ((normRHS_ < oldNormRHS) || (searchDone) ||
1367 
1368  if (Xyce::DEBUG_NONLINEAR && debugTimeFlag_ && nlParams.getDebugLevel() > 0 )
1369  {
1370  Xyce::dout() << "\n\tSearch Step: " << searchStep_ << ", Step Size: " << stepLength_ << std::endl
1371  << "\toldNormRHS: " << oldNormRHS << ", normRHS: " << normRHS_ << std::endl;
1372  }
1373  }
1374 
1375  return (normRHS_ < oldNormRHS);
1376 }
1377 
1378 //-----------------------------------------------------------------------------
1379 // Function : DampedNewton::backtrack_
1380 // Purpose : Does a backtracking line search along the current search
1381 // direction in the hopes of satisfying the Wolfe conditions
1382 // for f(x) = 0.5 F(x)'F(x) = normrhs^2.
1383 // Special Notes :
1384 // Scope : private
1385 // Creator : Tamara G. Kolda, SNL, Compuational Sciences and
1386 // Mathematics Research Department
1387 // Creation Date : 07/02/01
1388 //-----------------------------------------------------------------------------
1390 {
1391  bool issuffdec;
1392 
1393  // Local references
1394  N_LAS_Vector& g(*gradVectorPtr_); // Gradient
1395  N_LAS_Vector& s(*searchDirectionPtr_); // Search Direction
1396 
1397  // Compute the gradient, i.e., fill in g
1398  gradient_();
1399 
1400  // Evaluate <s,g>
1401  double sdotg = s.dotProduct(g);
1402 
1403  double& normrhs(normRHS_); // Norm of RHS
1404  const double f = 0.5 * normrhs * normrhs; // f(x) = 0.5 F(x)'F(X)
1405 
1406  // If s is a descent direction, then do the line search
1407  if (sdotg < 0)
1408  {
1409  issuffdec = simpleBt_(sdotg, f);
1410  }
1411 
1412  // If s is not a descent direction or we cannot find descent along s, switch
1413  // to the gradient. We must recompute the maximum allowable step with the new
1414  // search direction.
1415  if (((sdotg >= 0) || (!issuffdec)) &&
1417  {
1418  Xyce::lout() << "Switching to Cauchy!" << std::endl;
1419  s = g;
1420  s.scale(-1.0 / normrhs);
1421  sdotg = s.dotProduct(g);
1422  double& maxstep(constraintFactor_);
1423  maxstep = constrain_();
1424  issuffdec = simpleBt_(sdotg, f);
1425  }
1426 
1427  if (!issuffdec)
1428  {
1429  double& step(stepLength_);
1430  step = 0.0;
1431  updateX_();
1432  rhs_();
1433 
1434  Xyce::lout() << "Linesearch failed!" << std::endl;
1435  }
1436 
1437  return issuffdec;
1438 }
1439 
1440 //-----------------------------------------------------------------------------
1441 // Function : DampedNewton::simpleBt_
1442 // Purpose : Does a backtracking line search along the current sear
1443 // Special Notes :
1444 // Scope : private
1445 // Creator : Tamara G. Kolda, SNL, Compuational Sciences and
1446 // Mathematics Research Department
1447 // Creation Date : 07/02/01
1448 //-----------------------------------------------------------------------------
1449 bool DampedNewton::simpleBt_(double gsinit, double finit)
1450 {
1451  // Local references
1452  double& normrhs(normRHS_); // Norm of RHS
1453  double& step(stepLength_); // Step length
1454  double& maxstep(constraintFactor_); // Maximum step length
1455  unsigned int& nstep(searchStep_); // Number of backtracking steps
1456 
1457  // Local values
1458  double f = 0.5 * normrhs * normrhs; // Function value
1459  const double alpha = 1.0e-6; // Suff dec cond parameter
1460 
1461  Xyce::dout().setf(std::ios::scientific);
1462  Xyce::dout().precision(10);
1463 
1464  Xyce::dout() << "\nIteration: " << 0
1465  << " Step: " << 0
1466  << " F(X): " << finit
1467  << " gsinit: " << gsinit
1468  << " alpha * gsinit: " << alpha * gsinit
1469  << std::endl;
1470 
1471  // Initialize search
1472  const int nstepmax = 20; // Max number of backtracking steps
1473  const double minstep = N_UTL_MachineDependentParams::MachineEpsilon(); // Min step length
1474  nstep = 0;
1475 
1476  // Loop until we obtain sufficient decrease or violate some condition.
1477  bool issuffdec = false;
1478  while (!issuffdec)
1479  {
1480  // Stop if we've taken too many steps
1481  if (nstep >= nstepmax)
1482  break;
1483 
1484  // If this is the first iteration, set the step to maxstep.
1485  // Otherwise, reduce the step by a factor of one-half.
1486  if (nstep == 0)
1487  step = maxstep;
1488  else
1489  step *= 0.5;
1490 
1491  // Stop if the step has gotten too small
1492  if (step < minstep)
1493  break;
1494 
1495  // Compute the new X based on the currect step
1496  updateX_();
1497 
1498  // Evaluate the RHS for the new X
1499  rhs_();
1500 
1501  // Compute the new function value
1502  f = 0.5 * normrhs * normrhs;
1503 
1504  // Check the sufficient decrease condition
1505  issuffdec = (f <= (finit + (alpha * step * gsinit)));
1506 
1507  // Increment the step counter
1508  nstep ++;
1509 
1510  Xyce::dout() << "Iteration: " << nstep
1511  << " Step: " << step
1512  << " F(X): " << f
1513  << " Test: " << (finit + (alpha * step * gsinit))
1514  << std::endl;
1515  }
1516 
1517  Xyce::dout().unsetf(std::ios::scientific);
1518 
1519  return issuffdec;
1520 }
1521 
1522 //-----------------------------------------------------------------------------
1523 // Function : DampedNewton::backtrack_()
1524 // Purpose : Backtracking method based on Dennis & Shnabel.
1525 // Special Notes :
1526 // Scope : private
1527 // Creator : Scott A. Hutchinson, SNL, Computational Sciences Department
1528 // Creation Date : 07/28/01
1529 //-----------------------------------------------------------------------------
1531 {
1532  static const char *trace = "DampedNewton::computeStepLength_";
1533 
1534 #if defined(Xyce_DEBUG_NONLINEAR)
1535  char tmpStr[128];
1536 #endif
1537 
1538  // Want to get reduction in the norm, so we first copy the norm from
1539  // the previous iterate
1540  const double oldNormRHS = normRHS_;
1541 
1542  // The upper bound on the backtracking (damping) is determined by
1543  // the applied constraints.
1544  const double fullStep = constraintFactor_;
1545 
1546  // Starting value
1547  stepLength_ = fullStep;
1548 
1549  // Update the solution.
1550  updateX_();
1551 
1552  // Evaluate the new residual and take its norm:
1553  rhs_();
1554 
1555  // Linesearch...
1556  double rho, delta, theta, lambda;
1557  const double t = 1.0e-04, thetaMin = 0.1, thetaMax = 0.5;
1558  const double minStep = pow(N_UTL_MachineDependentParams::MachineEpsilon(), 0.33);
1559 
1560  // If we're using the inexact-Newton forcing, we use the forcing value
1561  // (linear-solver convergence tolerance) to set the initial "lambda".
1562  if (nlParams.getForcingFlag())
1563  lambda = 1.0 - nlParams.getForcingTerm();
1564  else
1565  lambda = 1.0;
1566 
1567  rho = normRHS_ / oldNormRHS;
1568 
1569  while ((rho > 1.0 - t * lambda) &&
1571  {
1572  delta = rho * rho - 1.0 + 2.0 * lambda;
1573  if (delta <= 0.0)
1574  theta = thetaMax;
1575  else
1576  {
1577  theta = lambda / delta;
1578 
1579  if (theta > thetaMax)
1580  theta = thetaMax;
1581  else if (theta < thetaMin)
1582  theta = thetaMin;
1583  }
1584 
1585  lambda *= theta;
1586  stepLength_ *= theta;
1587 
1588  // Jump out if the step-size gets too small...
1589  if (stepLength_ < minStep)
1590  {
1591  stepLength_ = minStep;
1593  }
1594 
1595  // Update solution & evaluate the new RHS vector and take norm:
1596  updateX_();
1597  rhs_();
1598 
1599  rho = normRHS_ / oldNormRHS;
1600  searchStep_++;
1601 
1602  if (Xyce::DEBUG_NONLINEAR && debugTimeFlag_ && nlParams.getDebugLevel() > 0 )
1603  {
1604  Xyce::dout() << "\n\tSearch Step: " << searchStep_ << ", Step Size: " << stepLength_ << std::endl
1605  << "\toldNormRHS: " << oldNormRHS << ", normRHS: " << normRHS_ << std::endl;
1606  }
1607  }
1608 
1609  // If the linesearch fails, take a full step...
1610 // if (rho > 1.0 - t * lambda)
1611 // stepLength_ = fullStep;
1612 
1613  return (rho <= 1.0 - t * lambda);
1614 
1615 }
1616 
1617 //-----------------------------------------------------------------------------
1618 // Function : DampedNewton::fullNewton_()
1619 // Purpose :
1620 // Special Notes :
1621 // Scope : private
1622 // Creator : Tamara G. Kolda, SNL, Compuational Sciences and
1623 // Mathematics Research Department
1624 // Creation Date : 06/29/01
1625 //-----------------------------------------------------------------------------
1627 {
1628  // ***** Full Newton method *****
1629  stepLength_ = 1.0;
1630 
1631  updateX_();
1632 
1633  // Evaluate the new residual and take its norm:
1634  rhs_();
1635 
1636  return true;
1637 }
1638 
1639 //-----------------------------------------------------------------------------
1640 // Function : DampedNewton::spiceNewton_()
1641 //
1642 // Purpose : Mimics spice. Assumes you've loaded things so
1643 // that the system you are solving is:
1644 //
1645 // J x_{i+1} = -f + J x_i
1646 //
1647 // Special Notes : This is ONLY used for debugging!
1648 //
1649 // At the moment, this will not get the norm right.
1650 //
1651 // Scope : private
1652 // Creator : Eric R. Keite, SNL, Compuational Sciences
1653 // Creation Date : 12/16/04
1654 //-----------------------------------------------------------------------------
1656 {
1657  // ***** Spice Newton method *****
1658  stepLength_ = 1.0;
1659 
1660  // Create the real dx. Normally, the search vector contains dx, but not
1661  // in the spice-style iteration. The search vector is actually the new x,
1662  // so dx = search-nextSol.
1664  -1.0,
1666 
1667  // copy "dx" into x. "dx" is really the new solution.
1668  // - should really use the "=" operator for this.
1669  (*nextSolVectorPtrPtr_)->putScalar(0.0);
1670  (*nextSolVectorPtrPtr_)->addVec(1.0, *searchDirectionPtr_);
1671 
1672  // now copy the tmp vector back into searchDirection, to set up the real dx.
1673  (searchDirectionPtr_)->putScalar(0.0);
1674  (searchDirectionPtr_)->addVec(1.0, *tmpVectorPtr_);
1675 
1676  // Evaluate the new residual and take its norm:
1677  rhs_();
1678 
1679  return true;
1680 }
1681 
1682 //-----------------------------------------------------------------------------
1683 // Function : DampedNewton::evalModNewton_
1684 // Purpose : This method evaluates the modified Newton status to
1685 // determine if a fresh Jacobian and associated preconditioner
1686 // are required. If we take a step that dramatically degrades
1687 // the NLS convergence, we ensure that we're going to use a
1688 // "fresh" Jacobian and preconditioner for some number of steps
1689 // to get us back to a region of good convergence.
1690 // Special Notes :
1691 // Scope : private
1692 // Creator : Scott A. Hutchinson, SNL, Computational Sciences Department
1693 // Creation Date : 08/17/01
1694 //-----------------------------------------------------------------------------
1696 
1697 {
1698  static const char * trace = "DampedNewton::evalModNewton_";
1699 
1700  const double minConvFactor = 0.01;
1701 
1702  const double etaMax = 1.0;
1703  const double etaMin = 1.0e-12;
1704  const double etaInit = 1.0e-01;
1705  double eta;
1706  double nlResNorm = normRHS_;
1707 
1708  if (modNewtonStep_ <= 2)
1709  {
1710  eta = etaInit;
1711  nlResNormOld = nlResNorm;
1712  }
1713  else
1714  {
1715 // double denom = Xycemax(N_UTL_MachineDependentParams::DoubleMin(), nlResNormOld);
1716 // eta = fabs(nlResNorm - lasSolverPtr_->getLinearResidual()) / denom;
1717  eta = resConvRate_;
1718  }
1719 
1720 #ifdef Xyce_DEBUG_NONLINEAR
1721  if (debugTimeFlag_ && nlParams.getDebugLevel() > 0 )
1722  {
1723  N_UTL_Param linRes( "RESIDUAL", 0.0 );
1724  lasSolverPtr_->getInfo( linRes );
1725  Xyce::dout() << "\tnlResNorm: " << nlResNorm << " linRes: " << linRes.getImmutableValue<double>()
1726  << "nlResNormOld: " << nlResNormOld << "calculated eta: " << eta << std::endl;
1727  }
1728 #endif
1729 
1730  etaOld = eta;
1731  nlResNormOld = nlResNorm;
1732 
1733 #ifdef Xyce_DEBUG_NONLINEAR
1734  if (debugTimeFlag_ && nlParams.getDebugLevel() > 0 )
1735  {
1736  Xyce::dout() << "\t\teta:\t" << eta << "\n" << std::endl;
1737  }
1738 #endif
1739 
1740 
1741 // lasSolverPtr_->setUpdatedJacobian(true);
1742 // directSolverPtr_->setUpdatedJacobian(true);
1743  loadJacobianFlag_ = true;
1744 
1745  if (eta < minConvFactor)
1746  {
1747  loadJacobianFlag_ = false;
1748 // lasSolverPtr_->setUpdatedJacobian(false);
1749 // directSolverPtr_->setUpdatedJacobian(false);
1750  }
1751 
1752 #ifdef Xyce_VERBOSE_NONLINEAR
1753  if (loadJacobianFlag_)
1754  Xyce::lout() << " ***** Calculating new Jacobian and preconditioner\n" << std::endl;
1755  else
1756  Xyce::lout() << " ***** Using old Jacobian and preconditioner\n" << std::endl;
1757 #endif
1758 
1759 }
1760 
1761 //-----------------------------------------------------------------------------
1762 // Function : DampedNewton::converged_
1763 // Purpose : This method checks for convergence of the Newton solver.
1764 // Much of this is currently based on the methods in both
1765 // P.N. Brown, A.C. Hindmarsh and L.R. Petzold, "Consistent
1766 // initial condition calculations for differential-algebraic
1767 // systems", SIAM J. Sci. Comput., Vol. 19, No. 5,
1768 // pp. 1495-1512, Sep. 1998 and J.E. Dennis, Jr. and
1769 // R.B. Schnabel, "Numerical Methods for Unconstrained
1770 // Optimization and Nonlinear Equations", Prentice-Hall, 1983.
1771 //
1772 // Special Notes : Returns 0 if not converged, positive if converged, negative
1773 // if not converged but some other error.
1774 //
1775 // Eric Keiter, 9233, 03/27/03: I modified the code so that
1776 // it would use values in the ReturnCodes class as
1777 // the integer values for success/failure, rather than
1778 // hardwired numbers. The defaults of this class are the
1779 // same numbers as the old hardwired ones.
1780 //
1781 // The reason for doing this is to allow the user, or
1782 // calling code to set which circumstances should be
1783 // considered adequate convergence. If
1784 // running a two-level continuation loop, within the transient mode,
1785 // there needed to be a way to easily disable the
1786 // "nearConverged" senario.
1787 //
1788 // Scope : private
1789 // Creator : Scott A. Hutchinson, SNL, Computational Sciences Department
1790 // Creation Date : 01/13/01
1791 //-----------------------------------------------------------------------------
1793 
1794 {
1795  static const char *trace = "DampedNewton::converged_";
1796 
1797  const double rhsTol = N_UTL_MachineDependentParams::MachineEpsilon();
1798  const double convRateMax = 0.5 * N_UTL_MachineDependentParams::DoubleMax();
1799 // const double minResReduct = 1.0e-3;
1800  const double minResReduct = 9.0e-1;
1801  const double stagnationTol = 1.0e-3;
1802 // const double stagnationTol = pow(N_UTL_MachineDependentParams::MachineEpsilon(), 1.0/3.0);
1803 
1805 
1806  // This parameter "takes-out" any damping induced in the size of the norm by
1807  // a line-search or other globalization method.
1808  double updateSize = wtNormDX_ / stepLength_;
1809 
1810  // Devices need to satisfy their own convergence criteria
1812  bool allDevicesConverged_ = loaderPtr_->allDevsConverged();
1813  if (!allDevicesConverged_ &&
1815  {
1816  return 0;
1817  }
1818  }
1819 
1820  // This test is for 2-level solves ONLY.
1821  bool innerDevicesConverged_ = loaderPtr_->innerDevsConverged();
1822  if (!innerDevicesConverged_ )
1823  {
1824  return 0;
1825  }
1826 
1827  // Make sure the residual isn't too small (i.e., there's nothing more to do)
1828  if (maxNormRHS_ <= rhsTol && maxNormRHS_ <= nlParams.getRHSTol())
1829  return retCodes_.normTooSmall; // 1;
1830 
1831  // Check for "normal" convergence
1832  if (maxNormRHS_ <= nlParams.getRHSTol() &&
1833  updateSize <= nlParams.getDeltaXTol())
1834  return retCodes_.normalConvergence; // 2;
1835 
1836  // Transient "Near Converged"...
1837  if (nlStep_ > nlParams.getMaxNewtonStep() &&
1838  mode1 == TRANSIENT)
1839  {
1840  // Check for "near" convergence and let the time integrator handle the
1841  // error analysis and decide on whether or not to accept its step
1842  if ((normRHS_rel_ <= minResReduct) && (resConvRate_ <= 1.0))
1843  return retCodes_.nearConvergence; // 3;
1844  }
1845 
1846  // Check to see if update is really small
1847  if (updateSize <= nlParams.getSmallUpdateTol ())
1848  return retCodes_.smallUpdate; // 4;
1849 
1850  // Next check the number of steps
1852  return retCodes_.tooManySteps; // -1;
1853 
1854  // Make sure that we haven't had too big an update
1855  else if (resConvRate_ > convRateMax)
1856  return retCodes_.updateTooBig; // -2;
1857 
1858  // Check for a stall in the convergence rate - if it is near one for five
1859  // consecutive steps, we've stalled so jump out.
1860  if (mode1 == TRANSIENT && fabs(resConvRate_ - 1.0) <= stagnationTol)
1861  {
1863 
1864  count++;
1865 #ifdef Xyce_DEBUG_NONLINEAR
1866  if (debugTimeFlag_ && nlParams.getDebugLevel() > 0 )
1867  {
1868  Xyce::dout() << "\tcount:\t" << count << std::endl
1869  << "\tresConvRate_:\t" << resConvRate_ << "\n" << std::endl;
1870  }
1871 #endif
1872  }
1873  else
1874  count = 0;
1875 
1876  if (mode1 == TRANSIENT && count == 5)
1877  {
1878  count = 0;
1879  if ((normRHS_rel_ < minResReduct) && (tmpConvRate <= 1.0))
1880  {
1881 #ifdef Xyce_DEBUG_NONLINEAR
1882  if (debugTimeFlag_ && nlParams.getDebugLevel() > 0 )
1883  {
1884  Xyce::dout() << "\ttmpConvRate: " << tmpConvRate << std::endl
1885  << "\tnormRHS_rel_:\t" << normRHS_rel_ << std::endl
1886  << "\tReturning 3\n" << std::endl;
1887  }
1888 #endif
1889  return retCodes_.nearConvergence; // 3;
1890  }
1891  else
1892  {
1893 #ifdef Xyce_DEBUG_NONLINEAR
1894  if (debugTimeFlag_ && nlParams.getDebugLevel() > 0 )
1895  {
1896  Xyce::dout() << "\ttmpConvRate: " << tmpConvRate << std::endl
1897  << "\tnormRHS_rel_:\t" << normRHS_rel_ << std::endl
1898  << "\tReturning -3\n" << std::endl;
1899  }
1900 #endif
1901  return retCodes_.stalled; // -3;
1902  }
1903  }
1904 
1905  return 0;
1906 
1907 }
1908 
1909 //-----------------------------------------------------------------------------
1910 // Function : DampedNewton::constrain_
1911 // Purpose : This method performs constraint backtracking by determining
1912 // a damping factor:
1913 //
1914 // theta = min{1, theta-, theta+, theta_u}
1915 //
1916 // where theta-, theta+ and theta_u are determined according
1917 // the constraint algorithm of J.N. Shadid (SNL internal
1918 // communication).
1919 // Special Notes : See the notes in ConstraintBT.C for how these are
1920 // calculated.
1921 // Scope : private
1922 // Creator : Scott A. Hutchinson, SNL, Computational Sciences Department
1923 // Creation Date : 02/08/01
1924 //-----------------------------------------------------------------------------
1926 {
1927  double minValue;
1928 
1929  // First, update all the thetas...
1932 
1935 
1938 
1939  // Now, determine the minimum value
1940  minValue = Xycemin(1.0, nlConstraintPtr_->getThetaBoundNeg());
1941  minValue = Xycemin(minValue, nlConstraintPtr_->getThetaBoundPos());
1942  minValue = Xycemin(minValue, nlConstraintPtr_->getThetaChange());
1943 
1944  if (Xyce::DEBUG_NONLINEAR)
1945  Xyce::dout() << "DampedNewton::constrain_: minValue: " << minValue << std::endl;
1946 
1947 #ifdef Xyce_VERBOSE_NONLINEAR
1948  if (minValue < 1.0)
1949  Xyce::lout() << " ***** Constraining:\t" << minValue << std::endl;
1950 #endif
1951 
1952  return minValue;
1953 }
1954 
1955 //-----------------------------------------------------------------------------
1956 // Function : DampedNewton::setForcing_
1957 // Purpose : This method calculates the forcing term (i.e., linear
1958 // residual tolerance) for iterative solvers based on the
1959 // method of Walker and Pernice (RefXXX)
1960 // Special Notes :
1961 // Scope : private
1962 // Creator : Scott A. Hutchinson, SNL, Computational Sciences Department
1963 // Creation Date : 05/03/01
1964 //-----------------------------------------------------------------------------
1965 void DampedNewton::setForcing_(const double nlResNorm)
1966 {
1967  const double etaExp = 0.5 * (1.0 + sqrt(5.0));
1968  const double etaMax = 0.1;
1969  const double etaMin = 1.0e-12;
1970  const double etaInit = 1.0e-01;
1971  double eta;
1972  double etaSafe;
1973 
1974  if (newtonStep_ == 0 && modNewtonStep_ == 0)
1975  {
1976  eta = etaInit;
1977  nlResNormOld = nlResNorm;
1978  }
1979 
1980  else if (nlResNormOld > N_UTL_MachineDependentParams::DoubleMin())
1981  {
1982  N_UTL_Param linRes( "RESIDUAL", 0.0 );
1983  lasSolverPtr_->getInfo( linRes );
1984  eta = fabs(nlResNorm - linRes.getImmutableValue<double>()) / nlResNormOld;
1985  eta *= eta;
1986 
1987  // First safeguard...
1988  etaSafe = pow(etaOld, etaExp);
1989  if (etaSafe > 0.1)
1990  eta = Xycemax(eta, etaSafe);
1991 
1992  // Second safeguard...
1993  eta = Xycemin(etaMax, eta);
1994 
1995  // Not too small...
1996  eta = Xycemax(etaMin, eta);
1997 
1998 #ifdef Xyce_DEBUG_NONLINEAR
1999  Xyce::dout() << "\tnlResNorm: " << nlResNorm
2000  << " linRes: " << linRes.getImmutableValue<double>()
2001  << " nlResNormOld: " << nlResNormOld
2002  << " calculated eta: " << fabs(nlResNorm - linRes.getImmutableValue<double>())/nlResNormOld << std::endl;
2003 #endif
2004  }
2005 
2006  else
2007  eta = etaMax;
2008 
2009  etaOld = eta;
2010  nlResNormOld = nlResNorm;
2011 
2012 #ifdef Xyce_VERBOSE_NONLINEAR
2013  Xyce::lout() << "\t\teta:\t" << eta << "\n" << std::endl;
2014 #endif
2015 
2016  // Set the tolerance
2017  lasSolverPtr_->setTolerance(eta);
2018 
2019  // Set the parameter
2020  nlParams.setForcingTerm(eta);
2021 
2022 }
2023 
2024 //-----------------------------------------------------------------------------
2025 // Function : DampedNewton::bankRose_(nlInfo)
2026 // Purpose : Method of R.E. Bank and D.J. Rose in their paper "Global
2027 // Approximate Newton Methods", 1981, Numerische Mathematik,
2028 // Springer-Verlag. Modifications to their original method
2029 // have been made so that this implementation looks more like
2030 // that implemented in CAzM - W.M. Coughran, Jr., E. Grosse and
2031 // D.J. Rose, "CAzM: A Circuit Analyzer with Macromodeling",
2032 // IEEE Trans. on Electron Devices, Vol. ED-30, No. 9,
2033 // Sep. 1983.
2034 // Special Notes :
2035 // Scope : private
2036 // Creator : Tamara G. Kolda, SNL, Compuational Sciences and
2037 // Mathematics Research Department
2038 // Creation Date : 06/29/01
2039 //-----------------------------------------------------------------------------
2041 {
2042  Xyce::lout() << "THIS IS BROKEN. DO NOT USE THIS METHOD" << std::endl;
2043  return false;
2044 }
2045 
2046 //-----------------------------------------------------------------------------
2047 // Function : DampedNewton::descent_(nlInfo)
2048 // Purpose :
2049 // Special Notes :
2050 // Scope : private
2051 // Creator : Tamara G. Kolda, SNL, Compuational Sciences and
2052 // Mathematics Research Department
2053 // Creation Date : 06/29/01
2054 //-----------------------------------------------------------------------------
2056 {
2057  Xyce::lout() << "THIS IS BROKEN. DO NOT USE THIS METHOD" << std::endl;
2058  return false;
2059 }
2060 
2061 //-----------------------------------------------------------------------------
2062 // Dummy function since homotopy doesn't work with old solver
2063 //-----------------------------------------------------------------------------
2065 {
2066  return true;
2067 }
2068 
2069 //-----------------------------------------------------------------------------
2070 // Dummy function since homotopy doesn't work with old solver
2071 //-----------------------------------------------------------------------------
2073 {
2074  return 0;
2075 }
2076 
2077 //-----------------------------------------------------------------------------
2078 // Dummy function since homotopy doesn't work with old solver
2079 //-----------------------------------------------------------------------------
2081 {
2082  return 0;
2083 }
2084 
2085 } // namespace Nonlinear
2086 } // namespace Xyce