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