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