Xyce  6.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
N_NLS_NonLinearSolver.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_NonLinearSolver.C,v $
27 //
28 // Purpose : Body file which declares an interface common to all
29 // supported nonlinear solver algorithms. The Manager class
30 // uses this interface to call a concrete algorithm.
31 //
32 // Special Notes : This is the "Strategy" class in the Strategy design
33 // pattern.
34 //
35 // Creator : Scott A. Hutchinson, SNL, Parallel Computational Sciences
36 //
37 // Creation Date : 04/28/00
38 //
39 // Revision Information:
40 // ---------------------
41 //
42 // Revision Number: $Revision: 1.125 $
43 //
44 // Revision Date : $Date: 2014/02/24 23:49:25 $
45 //
46 // Current Owner : $Author: tvrusso $
47 //-------------------------------------------------------------------------
48 
49 #include <Xyce_config.h>
50 
51 
52 // ---------- Standard Includes ----------
53 
54 // ---------- Xyce Includes ----------
55 
56 #include <N_UTL_fwd.h>
57 #include <N_UTL_OptionBlock.h>
58 #include <N_NLS_Manager.h>
59 #include <N_NLS_NonLinearSolver.h>
60 #include <N_NLS_ConstraintBT.h>
61 #include <N_NLS_TwoLevelNewton.h>
62 
63 #include <N_LAS_Matrix.h>
64 #include <N_LAS_Vector.h>
65 
66 #include <N_LAS_Solver.h>
67 #include <N_LAS_Problem.h>
68 #include <N_LAS_SolverFactory.h>
69 #include <N_LAS_PrecondFactory.h>
70 
71 #include <N_LAS_System.h>
72 #include <N_LAS_Builder.h>
73 
74 #include <N_ERH_ErrorMgr.h>
75 
77 #include <N_ANP_AnalysisManager.h>
78 
79 #include <N_LOA_Loader.h>
80 #include <N_IO_OutputMgr.h>
81 
82 #include <N_IO_CmdParse.h>
83 #include <N_PDS_Manager.h>
84 #include <N_PDS_ParComm.h>
85 
86 // Harmonic Balance matrix free stuff
88 
89 //-----------------------------------------------------------------------------
90 // Function : N_NLS_NonLinearSolver::N_NLS_NonLinearSolver
91 // Purpose : Constructor
92 // Special Notes :
93 // Scope : public
94 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
95 // Creation Date : 5/01/00
96 //-----------------------------------------------------------------------------
98  : commandLine_(cp),
99  nextSolVectorPtrPtr_(0),
100  currSolVectorPtrPtr_(0),
101  tmpSolVectorPtrPtr_(0),
102  rhsVectorPtr_(0),
103  jacobianMatrixPtr_(0),
104  gradVectorPtr_(0),
105  NewtonVectorPtr_(0),
106  solWtVectorPtr_(0),
107  lasSysPtr_(0),
108  lasSolverPtr_(0),
109  petraOptionBlockPtr_(0),
110  tlnPtr_(0),
111  loaderPtr_(0),
112  nlpMgrPtr_(0),
113  outMgrPtr_(0),
114  pdsMgrPtr_(0),
115  anaIntPtr_(0),
116  netlistFileName_(""),
117  outputStepNumber_(0),
118  debugTimeFlag_(true),
119  contStep_(0),
120  matrixFreeFlag_(false)
121 {
123 
124  if (commandLine_.getArgumentValue("netlist") != "")
125  {
126  netlistFileName_ = commandLine_.getArgumentValue("netlist");
127  }
128 }
129 
130 
131 //-----------------------------------------------------------------------------
132 // Function : N_NLS_NonLinearSolver::~N_NLS_NonLinearSolver
133 // Purpose : Destructor
134 // Special Notes :
135 // Scope : public
136 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
137 // Creation Date : 5/01/00
138 //-----------------------------------------------------------------------------
140 {
141  if (NewtonVectorPtr_)
142  {
143  delete NewtonVectorPtr_;
144  NewtonVectorPtr_ = 0;
145  }
146 
147  if (gradVectorPtr_)
148  {
149  delete gradVectorPtr_;
150  gradVectorPtr_ = 0;
151  }
152 
153  if (solWtVectorPtr_)
154  {
155  delete solWtVectorPtr_;
156  solWtVectorPtr_ = 0;
157  }
158 
160  {
161  delete petraOptionBlockPtr_;
163  }
164 
165  if (lasSolverPtr_)
166  {
167  delete lasSolverPtr_;
168  lasSolverPtr_ = 0;
169  }
170 
171 }
172 
173 //-----------------------------------------------------------------------------
174 // Function : N_NLS_NonLinearSolver::setPetraOptions
175 //
176 // Purpose : Passes option block corresponding to "LINSOL" onto
177 // nonlinear solver.
178 //
179 // Special Notes: These options are saved to be passed onto the object from
180 // registerLinearSolver() in the initializeAll() function.
181 // This is only called if "NONLIN" is present in the
182 // circuit file.
183 //
184 // Return Type : boolean
185 //
186 // - Input Arguments -
187 //
188 // OB : Option block containing options corresponding to
189 // "LINSOL" in the netlist.
190 // Scope : public
191 // Creator : Robert Hoekstra, SNL, Parallel Computational Sciences
192 // Creation Date : 11/9/00
193 //-----------------------------------------------------------------------------
194 bool N_NLS_NonLinearSolver::setPetraOptions(const N_UTL_OptionBlock & OB)
195 {
196  petraOptionBlockPtr_ = new N_UTL_OptionBlock(OB);
197  return (petraOptionBlockPtr_ != 0);
198 }
199 
200 //-----------------------------------------------------------------------------
201 // Function : N_NLS_NonLinearSolver::setDCOPRestartOptions
202 // Purpose :
203 // Special Notes :
204 // Scope : public
205 // Creator : Eric Keiter
206 // Creation Date : 09/17/2007
207 //-----------------------------------------------------------------------------
208 bool N_NLS_NonLinearSolver::setDCOPRestartOptions(const N_UTL_OptionBlock& OB)
209 {
210  std::string msg = "DCOP restart options not supported for this solver. Use nox instead. ";
211  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
212  return true;
213 }
214 
215 //-----------------------------------------------------------------------------
216 // Function : N_NLS_NonLinearSolver::setICOptions
217 // Purpose :
218 // Special Notes :
219 // Scope : public
220 // Creator : Eric Keiter
221 // Creation Date : 09/17/2007
222 //-----------------------------------------------------------------------------
223 bool N_NLS_NonLinearSolver::setICOptions(const N_UTL_OptionBlock& OB)
224 {
225  std::string msg = ".IC options not supported for this nonlinear solver. Use nox instead. ";
226  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
227  return true;
228 }
229 
230 //-----------------------------------------------------------------------------
231 // Function : N_NLS_NonLinearSolver::setNodeSetOptions
232 // Purpose :
233 // Special Notes :
234 // Scope : public
235 // Creator : Eric Keiter
236 // Creation Date : 09/25/2007
237 //-----------------------------------------------------------------------------
238 bool N_NLS_NonLinearSolver::setNodeSetOptions(const N_UTL_OptionBlock& OB)
239 {
240  std::string msg = ".NODESET options not supported for this nonlinear solver. Use nox instead. ";
241  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
242  return true;
243 }
244 
245 //-----------------------------------------------------------------------------
246 // Function : N_NLS_NonLinearSolver::setLocaOptions
247 // Purpose :
248 // Special Notes :
249 // Scope : public
250 // Creator : Eric Keiter
251 // Creation Date :
252 //-----------------------------------------------------------------------------
253 bool N_NLS_NonLinearSolver::setLocaOptions (const N_UTL_OptionBlock& OB)
254 {
255  std::string msg = "N_NLS_NonLinearSolver::setLocaOptions - not implemented for this solver. Use nox instead. ";
256  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
257  return true;
258 }
259 
260 //---------------------------------------------------------------------------
261 // Function : N_NLS_NonLinearSolver::setTwoLevelLocaOptions
262 // Purpose :
263 // Special Notes :
264 // Scope : public
265 // Creator : Eric Keiter
266 // Creation Date :
267 //-----------------------------------------------------------------------------
268 bool N_NLS_NonLinearSolver::setTwoLevelLocaOptions (const N_UTL_OptionBlock& OB)
269 {
270  std::string msg = "N_NLS_NonLinearSolver::setTwoLevelLocaOptions - not implemented for this solver. Use nox instead.";
271  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
272  return true;
273 }
274 
275 //---------------------------------------------------------------------------
276 // Function : N_NLS_NonLinearSolver::setTwoLevelOptions
277 // Purpose :
278 // Special Notes :
279 // Scope : public
280 // Creator : Eric Keiter
281 // Creation Date :
282 //---------------------------------------------------------------------------
283 bool N_NLS_NonLinearSolver::setTwoLevelOptions (const N_UTL_OptionBlock& OB)
284 {
285  return true;
286 }
287 
288 //---------------------------------------------------------------------------
289 // Function : N_NLS_NonLinearSolver::setTwoLevelTranOptions
290 // Purpose :
291 // Special Notes :
292 // Scope : public
293 // Creator : Eric Keiter
294 // Creation Date :
295 //---------------------------------------------------------------------------
296 bool N_NLS_NonLinearSolver::setTwoLevelTranOptions (const N_UTL_OptionBlock& OB)
297 {
298  return true;
299 }
300 
301 //-----------------------------------------------------------------------------
302 // Function : N_NLS_NonLinearSolver::registerRHSVector
303 // Purpose :
304 // Special Notes :
305 // Return Type : boolean
306 // Scope : public
307 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
308 // Creation Date : 5/01/00
309 //-----------------------------------------------------------------------------
310 bool N_NLS_NonLinearSolver::registerRHSVector(N_LAS_Vector* tmp_RHSVecPtr)
311 {
312  rhsVectorPtr_ = tmp_RHSVecPtr;
313  return (rhsVectorPtr_ != 0);
314 }
315 
316 //-----------------------------------------------------------------------------
317 // Function : N_NLS_NonLinearSolver::registerLoader
318 // Purpose :
319 // Special Notes :
320 // Return Type : boolean
321 // Scope : public
322 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
323 // Creation Date : 5/01/00
324 //-----------------------------------------------------------------------------
325 bool N_NLS_NonLinearSolver::registerLoader(N_LOA_Loader* tmp_LoaderPtr)
326 {
327  loaderPtr_ = tmp_LoaderPtr;
328  return (loaderPtr_ != 0);
329 }
330 
331 //-----------------------------------------------------------------------------
332 // Function : N_NLS_NonLinearSolver::registerLinearSystem
333 // Purpose :
334 // Special Notes :
335 // Return Type : boolean
336 // Scope : public
337 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
338 // Creation Date : 6/09/00
339 //-----------------------------------------------------------------------------
340 
341 bool N_NLS_NonLinearSolver::registerLinearSystem(N_LAS_System* tmp_LasSysPtr)
342 {
343  lasSysPtr_ = tmp_LasSysPtr;
344  return (lasSysPtr_ != 0);
345 }
346 
347 //-----------------------------------------------------------------------------
348 // Function : N_NLS_NonLinearSolver::registerPrecondFactory
349 // Purpose :
350 // Special Notes :
351 // Return Type : boolean
352 // Scope : public
353 // Creator : Heidi Thornquist, SNL, Electrical & Microsystem Modeling
354 // Creation Date : 11/11/08
355 //-----------------------------------------------------------------------------
356 bool N_NLS_NonLinearSolver::registerPrecondFactory(const RefCountPtr<N_LAS_PrecondFactory>& tmp_LasPrecPtr)
357 {
358  lasPrecPtr_ = tmp_LasPrecPtr;
359  return (!Teuchos::is_null(lasPrecPtr_));
360 }
361 
362 //-----------------------------------------------------------------------------
363 // Function : N_NLS_NonLinearSolver::registerParallelMgr
364 // Purpose :
365 // Special Notes :
366 // Return Type : boolean
367 // Scope : public
368 // Creator : Eric Keiter, Sandia
369 // Creation Date : 6/8/2013
370 //-----------------------------------------------------------------------------
371 bool N_NLS_NonLinearSolver::registerParallelMgr(N_PDS_Manager * pdsMgrPtr)
372 {
373  pdsMgrPtr_ = pdsMgrPtr;
374  return (pdsMgrPtr_ != 0);
375 }
376 
377 //-----------------------------------------------------------------------------
378 // Function : N_NLS_NonLinearSolver::registerAnalysisInterface
379 // Purpose :
380 // Special Notes :
381 // Return Type : boolean
382 // Scope : public
383 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
384 // Creation Date : 5/03/00
385 //-----------------------------------------------------------------------------
387 {
388  anaIntPtr_ = tmp_anaIntPtr;
389  return (anaIntPtr_ != 0);
390 }
391 
392 //-----------------------------------------------------------------------------
393 // Function : N_NLS_NonLinearSolver::registerOutputMgr
394 // Purpose :
395 // Special Notes :
396 // Scope : public
397 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
398 // Creation Date : 9/23/03
399 //-----------------------------------------------------------------------------
400 bool N_NLS_NonLinearSolver::registerOutputMgr (N_IO_OutputMgr * outPtr)
401 {
402  outMgrPtr_ = outPtr;
403  return (outMgrPtr_ != 0);
404 }
405 
406 //-----------------------------------------------------------------------------
407 // Function : N_NLS_NonLinearSolver::registerTwoLevelSolver
408 // Purpose : This function is called in the event that the two-level
409 // Newton method has been invoked.
410 // Special Notes :
411 // Scope : public
412 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
413 // Creation Date : 10/24/02
414 //-----------------------------------------------------------------------------
416  (N_NLS_TwoLevelNewton * tmp_tlnPtr)
417 {
418  tlnPtr_ = tmp_tlnPtr;
419  return (tlnPtr_ != 0);
420 }
421 
422 //-----------------------------------------------------------------------------
423 // Function : N_NLS_NonLinearSolver::registerParamMgr
424 // Purpose : This function is called in the event that the two-level
425 // Newton method has been invoked.
426 // Special Notes :
427 // Scope : public
428 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
429 // Creation Date : 10/24/02
430 //-----------------------------------------------------------------------------
433 {
434  nlpMgrPtr_ = ptr;
435  return (nlpMgrPtr_ != 0);
436 }
437 
438 //-----------------------------------------------------------------------------
439 // Function : N_NLS_NonLinearSolver::registerTopology
440 // Purpose : This function is called to register a topology object
441 // Special Notes :
442 // Scope : public
443 // Creator : Richard Schiek, SNL, 1437
444 // Creation Date : 12/3/08
445 //-----------------------------------------------------------------------------
446 bool N_NLS_NonLinearSolver::registerTopology(N_TOP_Topology * ptr)
447 {
448  topologyRcp_ = rcp( ptr, false);
449  return (! Teuchos::is_null(topologyRcp_) );
450 }
451 
452 //-----------------------------------------------------------------------------
453 // Function : N_NLS_NonLinearSolver::initializeAll
454 //
455 // Purpose : Called after all register and set functions.
456 // Once the various registrations have taken place,
457 // this function sets the remaining pointers.
458 //
459 // Special Notes: This function obtains the solution, temporary solution and
460 // rhs vectors from the LAS system class.
461 //
462 // Return Type : boolean
463 // Scope : public
464 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
465 // Creation Date : 6/12/00
466 //-----------------------------------------------------------------------------
468 {
469  bool bsuccess = true;
470 
471 
472  // Check the registerLinearSystem has been successfully called.
473  if (lasSysPtr_ == 0)
474  return false;
475 
476  // get the temporaries:
477  tmpSolVectorPtrPtr_ = lasSysPtr_->getTmpSolVectorPtr();
478  bsuccess = bsuccess && (tmpSolVectorPtrPtr_ != 0);
479 
480  // get rhs vector:
481  rhsVectorPtr_ = lasSysPtr_->getRHSVector();
482  bsuccess = bsuccess && (rhsVectorPtr_ != 0);
483 
484  // get current solution vectors:
485  currSolVectorPtrPtr_ = lasSysPtr_->getCurrSolVectorPtr();
486  bsuccess = bsuccess && (currSolVectorPtrPtr_ != 0);
487 
488  // get next solution vectors:
489  nextSolVectorPtrPtr_ = lasSysPtr_->getNextSolVectorPtr();
490  bsuccess = bsuccess && (nextSolVectorPtrPtr_ != 0);
491 
492  // get jacobian:
493  jacobianMatrixPtr_ = lasSysPtr_->getJacobianMatrix();
494  bsuccess = bsuccess && (jacobianMatrixPtr_ != 0);
495 
496  N_LAS_Builder & bld_ = lasSysPtr_->builder();
497 
498  // create gradient vector:
499  //gradVectorPtr_ = lasSysPtr_->builder().createVector();
500  gradVectorPtr_ = bld_.createVector();
501  bsuccess = bsuccess && (gradVectorPtr_ != 0);
502 
503  // create Newton update vector:
504  //NewtonVectorPtr_ = lasSysPtr_->builder().createVector();
505  NewtonVectorPtr_ = bld_.createVector();
506  bsuccess = bsuccess && (NewtonVectorPtr_ != 0);
507 
508  // create solution weighting vector:
509  //solWtVectorPtr_ = lasSysPtr_->builder().createVector();
510  solWtVectorPtr_ = bld_.createVector();
511  bsuccess = bsuccess && (solWtVectorPtr_ != 0);
512 
513  if( !petraOptionBlockPtr_ )
514  {
515  petraOptionBlockPtr_ = new N_UTL_OptionBlock();
516  petraOptionBlockPtr_->getParams().push_back( N_UTL_Param( "TYPE", "DEFAULT" ) );
517  }
518 
519  Teuchos::RefCountPtr<N_LAS_Vector> NewtonVectorRCPtr = Teuchos::rcp(NewtonVectorPtr_, false);
520  Teuchos::RefCountPtr<N_LAS_Vector> rhsVectorRCPtr = Teuchos::rcp(rhsVectorPtr_, false);
521 
522  if (!matrixFreeFlag_)
523  {
524  Teuchos::RefCountPtr<N_LAS_Matrix> jacobianMatrixRCPtr = Teuchos::rcp(jacobianMatrixPtr_,false);
525  // Normal full matrix linear solver options
526  lasProblemRCPtr_ = rcp( new N_LAS_Problem( jacobianMatrixRCPtr,
527  NewtonVectorRCPtr,
528  rhsVectorRCPtr) );
529  }
530  else
531  {
532  // Matrix free harmonic balance linear solver options
533  // Create MatrixFreeLinearProblem
534  Teuchos::RefCountPtr<N_NLS_NonLinearSolver> NonlinearSolverRCPtr = Teuchos::rcp(this, false);
535  Teuchos::RefCountPtr<N_NLS_MatrixFreeEpetraOperator>
536  matFreeOp = matrixFreeEpetraOperator(
537  NonlinearSolverRCPtr,
538  NewtonVectorRCPtr,
539  rhsVectorRCPtr,
540  //lasSysPtr_->builder().getSolutionMap()
541  bld_.getSolutionMap()
542  );
543 
544  Teuchos::RefCountPtr<Epetra_Operator> epetraOperator = Teuchos::rcp_dynamic_cast<Epetra_Operator>(matFreeOp,true);
545  // Create N_LAS_Problem
546  lasProblemRCPtr_ = rcp( new N_LAS_Problem(
547  epetraOperator,
548  NewtonVectorRCPtr,
549  rhsVectorRCPtr
550  )
551  );
552 
553  }
554 
555  lasSolverPtr_ = N_LAS_SolverFactory::create( *petraOptionBlockPtr_,
557 
558  // If a preconditioner factory has been provided by the analysis package,
559  // use it to generate a preconditioner for the linear solver.
560  if (!Teuchos::is_null(lasPrecPtr_)) {
561  Teuchos::RefCountPtr<N_LAS_Preconditioner> precond = lasPrecPtr_->create( Teuchos::rcp( lasSysPtr_, false ) );
562  lasSolverPtr_->setPreconditioner( precond );
563  }
564 
565 #ifdef Xyce_DEBUG_NONLINEAR
566  Xyce::dout() << "size of solution vector: " << lasSysPtr_->getGlobalSolutionSize() << std::endl
567  << "size of state vector: " << lasSysPtr_->getGlobalStateSize() << std::endl
568  << "End of N_NLS_NonLinearSolver::initializeAll\n";
569 #endif
570 
571  return bsuccess;
572 }
573 
574 #ifdef Xyce_DEBUG_NONLINEAR
575 //-----------------------------------------------------------------------------
576 // Function : N_NLS_NonLinearSolver::debugOutput1_
577 // Purpose : Write Jacobian to the file matrix.(n).txt, and the residual
578 // to the file rhs.(n).txt
579 //
580 // Special Notes : These are objects that will be availabled prior to the
581 // linear solve performed at each Newton step.
582 //
583 // Scope : public
584 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
585 // Creation Date : 10/30/00
586 //-----------------------------------------------------------------------------
587 void N_NLS_NonLinearSolver::debugOutput1
588  (N_LAS_Matrix & jacobian, N_LAS_Vector & rhs)
589 {
590  int debugLevel = getDebugLevel();
591  int newtStep = getNumIterations();
592  int screenOutput = getScreenOutputFlag ();
593  int contStep = getContinuationStep();
594  int paramNumber = getParameterNumber ();
595 
596  if (!debugTimeFlag_ || debugLevel < 1) return;
597 
598  char filename1[256]; for (int ich = 0; ich < 256; ++ich) filename1[ich] = 0;
599  char filename2[256]; for (int ich = 0; ich < 256; ++ich) filename2[ich] = 0;
600 
601  if (debugLevel >= 3)
602  {
603  sprintf(filename1, "matrix_%03d_%03d_%03d_%03d.txt",outputStepNumber_,paramNumber,contStep,newtStep);
604  }
605  else if (debugLevel == 2)
606  {
607  sprintf(filename1, "matrix_%03d_%03d.txt",outputStepNumber_,newtStep);
608  }
609  else
610  {
611  sprintf(filename1, "matrix_%03d.txt", newtStep);
612  }
613 
614  jacobian.writeToFile(filename1,false, getMMFormat () );
615 
616  if (screenOutput == 1)
617  {
618  Xyce::dout() << "\n\t***** Jacobian matrix:" << std::endl;
619  jacobian.printPetraObject(Xyce::dout());
620  }
621 
622  if (debugLevel >= 3)
623  {
624  sprintf(filename2, "rhs_%03d_%03d_%03d_%03d.txt",outputStepNumber_,paramNumber,contStep,newtStep);
625  }
626  else if (debugLevel == 2)
627  {
628  sprintf(filename2, "rhs_%03d_%03d.txt",outputStepNumber_,newtStep);
629  }
630  else
631  {
632  sprintf(filename2, "rhs_%03d.txt", newtStep);
633  }
634 
635  if (screenOutput == 1)
636  {
637  Xyce::dout() << "\n\t***** RHS vector:" << std::endl;
638 
639  rhs.printPetraObject(Xyce::dout());
640  }
641 
642  rhs.writeToFile(filename2);
643 
644 #ifdef Xyce_DEBUG_VOLTLIM
645  debugOutputJDX_VOLTLIM ();
646 #endif
647 
648  debugOutputDAE ();
649 
650 }
651 #endif // Xyce_DEBUG_NONLINEAR
652 
653 #ifdef Xyce_DEBUG_VOLTLIM
654 //-----------------------------------------------------------------------------
655 // Function : N_NLS_NonLinearSolver::debugOutputJDX_VOLTLIM_
656 // Purpose : Write JDX vector to output files.
657 // Special Notes : This requires a matvec.
658 // Scope : public
659 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
660 // Creation Date : 01/05/05
661 //-----------------------------------------------------------------------------
662 void N_NLS_NonLinearSolver::debugOutputJDX_VOLTLIM()
663 {
664  int debugLevel = getDebugLevel();
665  int newtStep = getNumIterations();
666  int contStep = getContinuationStep();
667  int paramNumber = getParameterNumber ();
668 
669  char filename1[256]; for (int ich = 0; ich < 256; ++ich) filename1[ich] = 0;
670  char filename2[256]; for (int ich = 0; ich < 256; ++ich) filename2[ich] = 0;
671  char filename3[256]; for (int ich = 0; ich < 256; ++ich) filename3[ich] = 0;
672 
673  if (debugLevel >= 3)
674  {
675  sprintf(filename1, "jdxVL_%03d_%03d_%03d_%03d.txt",outputStepNumber_,paramNumber,contStep,newtStep);
676  sprintf(filename2, "fdxVL_%03d_%03d_%03d_%03d.txt",outputStepNumber_,paramNumber,contStep,newtStep);
677  sprintf(filename3, "qdxVL_%03d_%03d_%03d_%03d.txt",outputStepNumber_,paramNumber,contStep,newtStep);
678  }
679  else if (debugLevel == 2)
680  {
681  sprintf(filename1, "jdxVL_%03d_%03d.txt",outputStepNumber_,newtStep);
682  sprintf(filename2, "fdxVL_%03d_%03d.txt",outputStepNumber_,newtStep);
683  sprintf(filename3, "qdxVL_%03d_%03d.txt",outputStepNumber_,newtStep);
684  }
685  else
686  {
687  sprintf(filename1, "jdxVL_%03d.txt", newtStep);
688  sprintf(filename2, "fdxVL_%03d.txt", newtStep);
689  sprintf(filename3, "qdxVL_%03d.txt", newtStep);
690  }
691 
692  bool Transpose = false; // if set to true, the matvec does the transpose.
693 
694  jdxVLVectorPtr_->putScalar(0.0);
695  fdxVLVectorPtr_->putScalar(0.0);
696  qdxVLVectorPtr_->putScalar(0.0);
697 
698  jacTestMatrixPtr_->matvec( Transpose , *dxVoltlimVectorPtr_, *jdxVLVectorPtr_);
699  dFdxTestMatrixPtr_->matvec( Transpose , *dxVoltlimVectorPtr_, *fdxVLVectorPtr_);
700  dQdxTestMatrixPtr_->matvec( Transpose , *dxVoltlimVectorPtr_, *qdxVLVectorPtr_);
701 
702  jdxVLVectorPtr_->writeToFile(filename1);
703  fdxVLVectorPtr_->writeToFile(filename2);
704  qdxVLVectorPtr_->writeToFile(filename3);
705 
706  if (debugLevel >= 3)
707  {
708  sprintf(filename1, "jtest_%03d_%03d_%03d_%03d.txt",outputStepNumber_,paramNumber,contStep,newtStep);
709  sprintf(filename2, "ftest_%03d_%03d_%03d_%03d.txt",outputStepNumber_,paramNumber,contStep,newtStep);
710  sprintf(filename3, "qtest_%03d_%03d_%03d_%03d.txt",outputStepNumber_,paramNumber,contStep,newtStep);
711  }
712  else if (debugLevel == 2)
713  {
714  sprintf(filename1, "jtest_%03d_%03d.txt",outputStepNumber_,newtStep);
715  sprintf(filename2, "ftest_%03d_%03d.txt",outputStepNumber_,newtStep);
716  sprintf(filename3, "qtest_%03d_%03d.txt",outputStepNumber_,newtStep);
717  }
718  else
719  {
720  sprintf(filename1, "jtest_%03d.txt", newtStep);
721  sprintf(filename2, "ftest_%03d.txt", newtStep);
722  sprintf(filename3, "qtest_%03d.txt", newtStep);
723  }
724 
725  jacTestMatrixPtr_->writeToFile(filename1);
726  dFdxTestMatrixPtr_->writeToFile(filename2);
727  dQdxTestMatrixPtr_->writeToFile(filename3);
728 
729 }
730 #endif
731 
732 #ifdef Xyce_DEBUG_NONLINEAR
733 //-----------------------------------------------------------------------------
734 // Function : N_NLS_NonLinearSolver::debugOutputDAE
735 // Purpose : Write DAE vectors and matrices to output files.
736 // Special Notes :
737 // Scope : public
738 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
739 // Creation Date : 01/27/03
740 //-----------------------------------------------------------------------------
741 void N_NLS_NonLinearSolver::debugOutputDAE()
742 {
743  int debugLevel = getDebugLevel();
744  int newtStep = getNumIterations();
745  int contStep = getContinuationStep();
746  int paramNumber = getParameterNumber ();
747 
748  char filename1[256]; for (int ich = 0; ich < 256; ++ich) filename1[ich] = 0;
749  char filename2[256]; for (int ich = 0; ich < 256; ++ich) filename2[ich] = 0;
750 
751  char filename4[256]; for (int ich = 0; ich < 256; ++ich) filename4[ich] = 0;
752  char filename6[256]; for (int ich = 0; ich < 256; ++ich) filename6[ich] = 0;
753  char filename7[256]; for (int ich = 0; ich < 256; ++ich) filename7[ich] = 0;
754  char filename8[256]; for (int ich = 0; ich < 256; ++ich) filename8[ich] = 0;
755  char filename9[256]; for (int ich = 0; ich < 256; ++ich) filename9[ich] = 0;
756 
757  N_LAS_Matrix *dQdx = lasSysPtr_->getDAEdQdxMatrix ();
758  N_LAS_Matrix *dFdx = lasSysPtr_->getDAEdFdxMatrix ();
759 
760  N_LAS_Vector *daeQ = lasSysPtr_->getDAEQVector();
761  N_LAS_Vector *daeF = lasSysPtr_->getDAEFVector();
762 
763  N_LAS_Vector *daeFlim = lasSysPtr_->getdFdxdVpVector ();
764  N_LAS_Vector *daeQlim = lasSysPtr_->getdQdxdVpVector ();
765 
766  //Xyce::dout() << "In debugOutputDAE" << std::endl;
767 
768  if (debugLevel >= 3)
769  {
770  sprintf(filename1, "dQdx_%03d_%03d_%03d_%03d.txt" ,outputStepNumber_,paramNumber,contStep,newtStep);
771  sprintf(filename2, "dFdx_%03d_%03d_%03d_%03d.txt" ,outputStepNumber_,paramNumber,contStep,newtStep);
772 
773  sprintf(filename4, "daeQ_%03d_%03d_%03d_%03d.txt" ,outputStepNumber_,paramNumber,contStep,newtStep);
774  sprintf(filename6, "daeF_%03d_%03d_%03d_%03d.txt" ,outputStepNumber_,paramNumber,contStep,newtStep);
775 
776  sprintf(filename8, "daeQlim_%03d_%03d_%03d_%03d.txt" ,outputStepNumber_,paramNumber,contStep,newtStep);
777  sprintf(filename9, "daeFlim_%03d_%03d_%03d_%03d.txt" ,outputStepNumber_,paramNumber,contStep,newtStep);
778  }
779  else if (debugLevel >= 2)
780  {
781  sprintf(filename1, "dQdx_%03d_%03d.txt" ,outputStepNumber_,newtStep);
782  sprintf(filename2, "dFdx_%03d_%03d.txt" ,outputStepNumber_,newtStep);
783 
784  sprintf(filename4, "daeQ_%03d_%03d.txt" ,outputStepNumber_,newtStep);
785  sprintf(filename6, "daeF_%03d_%03d.txt" ,outputStepNumber_,newtStep);
786 
787  sprintf(filename8, "daeQlim_%03d_%03d.txt" ,outputStepNumber_,newtStep);
788  sprintf(filename9, "daeFlim_%03d_%03d.txt" ,outputStepNumber_,newtStep);
789  }
790  else
791  {
792  sprintf(filename1, "dQdx_%03d.txt" , newtStep);
793  sprintf(filename2, "dFdx_%03d.txt" , newtStep);
794 
795  sprintf(filename4, "daeQ_%03d.txt" , newtStep);
796  sprintf(filename6, "daeF_%03d.txt" , newtStep);
797 
798  sprintf(filename8, "daeQlim_%03d.txt" , newtStep);
799  sprintf(filename9, "daeFlim_%03d.txt" , newtStep);
800  }
801 
802  // write the matrices:
803  dQdx->writeToFile (filename1,false, getMMFormat () );
804  dFdx->writeToFile (filename2,false, getMMFormat () );
805 
806  // write the vectors:
807  daeQ->writeToFile(filename4);
808  daeF->writeToFile(filename6);
809  daeQlim->writeToFile(filename8);
810  daeFlim->writeToFile(filename9);
811 }
812 #endif // Xyce_DEBUG_NONLINEAR
813 
814 #ifdef Xyce_DEBUG_NONLINEAR
815 //-----------------------------------------------------------------------------
816 // Function : N_NLS_NonLinearSolver::debugOutput3_
817 // Purpose : Write out the update vector and the new solution.
818 //
819 // Special Notes : These are objects that will be availabled *after* the
820 // linear solve performed at each Newton step. That
821 // differentiates this function from debugOutput1_.
822 //
823 // Scope : public
824 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
825 // Creation Date : 10/30/00
826 //-----------------------------------------------------------------------------
827 void N_NLS_NonLinearSolver::debugOutput3
828  (N_LAS_Vector & dxVector, N_LAS_Vector & xVector)
829 {
830  int debugLevel = getDebugLevel();
831  int nlStep = getNumIterations();
832  int contStep = getContinuationStep();
833  int paramNumber = getParameterNumber ();
834 
835  if (!debugTimeFlag_ || debugLevel < 1) return;
836 
837  char filename[256]; for (int ich = 0; ich < 256; ++ich) filename[ich] = 0;
838 
839  if (debugLevel >= 3)
840  {
841  sprintf(filename, "update_%03d_%03d_%03d_%03d.txt",outputStepNumber_,paramNumber,contStep,nlStep);
842  }
843  else if (debugLevel == 2)
844  {
845  sprintf(filename, "update_%03d_%03d.txt",outputStepNumber_,nlStep);
846  }
847  else
848  {
849  sprintf(filename, "update_%03d.txt", nlStep);
850  }
851  dxVector.writeToFile(filename);
852 
853 #if 0
854  if (nlParams.getScreenOutputFlag () )
855  {
856  Xyce::dout() << "\n\t***** Update vector:" << std::endl << std::endl;
857  //searchDirectionPtr_->printPetraObject();
858  dxVector.printPetraObject();
859  }
860 #endif
861 
862 
863  if (debugLevel >= 3)
864  {
865  sprintf(filename, "solution_%03d_%03d_%03d_%03d.txt",outputStepNumber_,paramNumber,contStep,nlStep);
866  }
867  if (debugLevel == 2)
868  {
869  sprintf(filename,"solution_%03d_%03d.txt",outputStepNumber_,nlStep);
870  }
871  else
872  {
873  sprintf(filename, "solution_%03d.txt", nlStep);
874  }
875  xVector.writeToFile(filename);
876 
877 #if 0
878  if (nlParams.getScreenOutputFlag () )
879  {
880  Xyce::dout() << "\n\t***** Solution vector:" << std::endl;
881  //x->printPetraObject();
882  xVector.printPetraObject();
883  }
884 #endif
885 
886 }
887 #endif
888 
889 //-----------------------------------------------------------------------------
890 // Function : N_NLS_NonLinearSolver::resetCountersAndTimers_
891 //
892 // Purpose : Resets all the counters and timers in this object.
893 //
894 // Scope : protected
895 // Creator : Tamara G. Kolda, SNL, CSMR (8950)
896 // Eric Keiter, SNL, Parallel Computational Sciences. (9233)
897 // Creation Date : 01/24/02
898 //-----------------------------------------------------------------------------
900 {
901  numJacobianLoads_ = 0;
903  numLinearSolves_ = 0;
905  numResidualLoads_ = 0;
907  totalLinearSolveTime_ = 0.0;
910 }
911 
912 //-----------------------------------------------------------------------------
913 // Function : N_NLS_NonLinearSolver::setX0_()
914 //
915 // Purpose : This should be called at the beginning of each nonlinear
916 // iteration. Copies information from nextSolVector (and
917 // related vectors that are important but hidden from the
918 // nonlinear solver) into tmpSolVector.
919 //
920 // Return Type : boolean
921 // Scope : protected
922 // Creator : Tamara G. Kolda, SNL, CSMR (8950)
923 // Eric Keiter, SNL, Parallel Computational Sciences. (9233)
924 // Creation Date : 01/24/02
925 //-----------------------------------------------------------------------------
927 {
929  return true;
930 }
931 
932 //-----------------------------------------------------------------------------
933 // Function : N_NLS_NonLinearSolver::rhs_()
934 //
935 // Purpose : Calculates the RHS corresponding to the current solution
936 // vector. More specifically, it fills in the content of
937 // RHSVectorPtr_ based on the contents of nextSolVectorPtr_.
938 //
939 // Special Notes : The rhsVectorPtr_ is really the NEGATIVE of F(x).
940 //
941 // Scope : private
942 // Creator : Tamara G. Kolda, SNL, Compuational Sciences and
943 // Mathematics Research Department
944 // Creation Date : 06/19/01
945 //-----------------------------------------------------------------------------
946 
948 {
949  loaderPtr_->loadRHS();
951  totalResidualLoadTime_ += loaderPtr_->getResidualTime();
952 
953  return true;
954 }
955 
956 //-----------------------------------------------------------------------------
957 // Function : N_NLS_NonLinearSolver::jacobian_()
958 //
959 // Purpose : Calculates the Jacobian corresponding to the current
960 // solution vector. More specifically, it fills in the
961 // content of jacobianMatrixPtr_ based on the contents of
962 // nextSolVectorPtr_.
963 //
964 // Return Type : boolean
965 // Scope : private
966 // Creator : Tamara G. Kolda, SNL, Compuational Sciences and
967 // Mathematics Research Department
968 // Creation Date : 06/19/01
969 //-----------------------------------------------------------------------------
971 {
972  loaderPtr_->loadJacobian();
974  totalJacobianLoadTime_ += loaderPtr_->getJacobianTime();
975  return true;
976 }
977 
978 //-----------------------------------------------------------------------------
979 // Function : N_NLS_NonLinearSolver::applyJacobian()
980 //
981 // Purpose : Applies the Jacobian corresponding to the current
982 // solution vector.
983 //
984 // Return Type : boolean
985 // Scope : private
986 // Creator : Todd Coffey, Ting Mei
987 // Creation Date : 07/29/08
988 //-----------------------------------------------------------------------------
989 bool N_NLS_NonLinearSolver::applyJacobian(const N_LAS_Vector& input, N_LAS_Vector& result)
990 {
991  loaderPtr_->applyJacobian(input,result);
993  totalJacobianLoadTime_ += loaderPtr_->getJacobianTime();
994  return true;
995 }
996 
997 //-----------------------------------------------------------------------------
998 // Function : N_NLS_NonLinearSolver::newton_
999 //
1000 // Purpose : Calculates the Newton direction corresponding to the
1001 // current RHS and Jacobian matrix.
1002 //
1003 // Return Type : boolean
1004 //
1005 // Scope : private
1006 // Creator : Tamara G. Kolda, SNL, Compuational Sciences and
1007 // Mathematics Research Department
1008 // Creation Date : 06/19/01
1009 //-----------------------------------------------------------------------------
1011 {
1012  int solutionStatus = lasSolverPtr_->solve( false );
1013 
1014  totalLinearSolveTime_ += lasSolverPtr_->solutionTime();
1015  ++numLinearSolves_;
1016 
1017  if( lasSolverPtr_->isIterative() )
1018  {
1019  N_UTL_Param param( "Iterations", 0 );
1020  lasSolverPtr_->getInfo( param );
1021  totalNumLinearIters_ += param.getImmutableValue<int>();
1022 
1023  if( solutionStatus ) ++numFailedLinearSolves_;
1024  }
1025  else
1026  {
1027  N_UTL_Param param( "Refactored", 0 );
1028  lasSolverPtr_->getInfo( param );
1029  if( param.getImmutableValue<int>() ) ++numJacobianFactorizations_;
1030  if( solutionStatus ) ++numFailedLinearSolves_;
1031  }
1032 
1033  if( solutionStatus ) return false;
1034 
1035  return true;
1036 }
1037 
1038 //-----------------------------------------------------------------------------
1039 // Function : N_NLS_NonLinearSolver::gradient_()
1040 // Purpose : Calculates the Gradient direction corresponding to the
1041 // current RHS and Jacobian matrix.
1042 // Computes gradient using jacobianMatrixPtr_ and
1043 // the rhsVectorPtr_. On output, gradVectorPtr_ contains
1044 // the gradient of 0.5 * ||F(x)||^2.
1045 //
1046 // Special Notes : The rhsVectorPtr_ is really the NEGATIVE of F(x).
1047 //
1048 // Return Type : boolean
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  // Compute gradient = Jacobian' * RHS
1057  bool transpose = true;
1058  jacobianMatrixPtr_->matvec(transpose, *rhsVectorPtr_, *gradVectorPtr_);
1059 
1060  // We have to scale by -1 because the rhsVectorPtr_ is really the
1061  // NEGATIVE of F(X). This gives us gradient = Jacobian' * F(X).
1062  gradVectorPtr_->scale(-1.0);
1063 
1064  return true;
1065 }
1066 
1067 //-----------------------------------------------------------------------------
1068 // Function : N_NLS_NonLinearSolver::getCouplingMode ()
1069 // Purpose :
1070 // Special Notes :
1071 // Scope : public
1072 // Creator : Eric R. Keiter, SNL, Compuational Sciences and
1073 // Creation Date : 12/05/02
1074 //-----------------------------------------------------------------------------
1076 {
1077  return FULL_PROBLEM;
1078 }
1079 
1080 //-----------------------------------------------------------------------------
1081 // Function : N_NLS_NonLinearSolver::setMatrixFreeFlag
1082 // Purpose :
1083 // Special Notes :
1084 // Scope : public
1085 // Creator : Todd Coffey, Ting Mei
1086 // Creation Date : 7/29/08
1087 //-----------------------------------------------------------------------------
1089 {
1090  matrixFreeFlag_ = matrixFreeFlag;
1091 }
1092 
1093 //-----------------------------------------------------------------------------
1094 // Function : N_NLS_NonLinearSolver::getMatrixFreeFlag
1095 // Purpose :
1096 // Special Notes :
1097 // Scope : public
1098 // Creator : Todd Coffey, Ting Mei
1099 // Creation Date : 7/29/08
1100 //-----------------------------------------------------------------------------
1102 {
1103  return matrixFreeFlag_;
1104 }
1105 
1106 #ifdef Xyce_DEBUG_NONLINEAR
1107 //-----------------------------------------------------------------------------
1108 // Function : N_NLS_NonLinearSolver::setDebugFlags
1109 // Purpose :
1110 // Special Notes :
1111 // Scope : private
1112 // Creator : Eric R. Keiter, SNL, Compuational Sciences
1113 // Creation Date : 09/23/01
1114 //-----------------------------------------------------------------------------
1115 void N_NLS_NonLinearSolver::setDebugFlags()
1116 {
1118 
1119  double currTime = anaIntPtr_->getTime();
1120 
1121  debugTimeFlag_ =
1122  (currTime >= getDebugMinTime() &&
1123  currTime <= getDebugMaxTime() ) &&
1124  (outputStepNumber_ >= getDebugMinTimeStep() &&
1125  outputStepNumber_ <= getDebugMaxTimeStep());
1126 
1127  bool steadyStateFlag = anaIntPtr_->getAnalysisMgr ()->getSteadyStateFlag ();
1128 
1129  if (tlnPtr_ != 0)
1131  else
1132  contStep_ = 0;
1133 
1134  if (steadyStateFlag == true)
1135  {
1136  outputStepNumber_ -= 1;
1137  }
1138 }
1139 #endif
1140