Xyce  6.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
N_ANP_AnalysisBase.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_ANP_AnalysisBase.C,v $
27 // Purpose : Base class analysis functions.
28 // Special Notes :
29 // Creator : Richard Schiek, SNL, Electrical and Microsystem Modeling
30 // Creation Date : 01/24/08
31 //
32 // Revision Information:
33 // ---------------------
34 // Revision Number: $Revision: 1.35 $
35 // Revision Date : $Date: 2014/02/24 23:49:12 $
36 // Current Owner : $Author: tvrusso $
37 //-----------------------------------------------------------------------------
38 #include <Xyce_config.h>
39 
40 #include <N_ANP_AnalysisBase.h>
41 
42 #include <N_ANP_AnalysisManager.h>
43 #include <N_ERH_Message.h>
44 
45 namespace Xyce {
46 namespace Analysis {
47 
48 //-----------------------------------------------------------------------------
49 // Function : AnalysisBase::AnalysisBase
50 // Purpose : Constructor
51 // Special Notes :
52 // Scope : public
53 // Creator : Todd S. Coffey, SNL.
54 // Creation Date : 01/29/08
55 //-----------------------------------------------------------------------------
57  : tiaParams(anaManagerPtr->tiaParams),
58  beginningIntegration(true),
59  integrationMethod_(TIAMethod_NONE),
60  stepNumber(0),
61  tranStepNumber(0),
62  totalNumberSuccessfulStepsTaken_(0),
63  totalNumberSuccessStepsThisParameter_(0),
64  totalNumberFailedStepsAttempted_(0),
65  totalNumberJacobiansEvaluated_(0),
66  totalNumberIterationMatrixFactorizations_(0),
67  totalNumberLinearSolves_(0),
68  totalNumberFailedLinearSolves_(0),
69  totalNumberLinearIters_(0),
70  totalNumberResidualEvaluations_(0),
71  totalNonlinearConvergenceFailures_(0),
72  totalLinearSolutionTime_(0.0),
73  totalResidualLoadTime_(0.0),
74  totalJacobianLoadTime_(0.0),
75  doubleDCOPFlag_(false),
76  doubleDCOPStep_(0),
77  sensFlag_(false),
78  inputOPFlag_(false),
79  commandLine_(anaManagerPtr->getCommandLine())
80 {
81  anaManagerRCPtr_ = rcp(anaManagerPtr, false );
82  assemblerRCPtr_ = anaManagerPtr->assemblerPtr;
83  lasSystemRCPtr_ = anaManagerPtr->lasSysPtr;
84  loaderRCPtr_ = anaManagerPtr->loaderPtr;
85  nlsMgrRCPtr_ = anaManagerPtr->nlsMgrPtr;
87  secRCPtr_ = anaManagerPtr->secPtr_;
88  wimRCPtr_ = anaManagerPtr->wimPtr;
89 }
90 
91 //-----------------------------------------------------------------------------
92 // Function : AnalysisBase::~AnalysisBase()
93 // Purpose : destructor
94 // Special Notes :
95 // Scope :
96 // Creator : Todd S. Coffey, SNL.
97 // Creation Date : 01/29/08
98 //-----------------------------------------------------------------------------
100 {
101 }
102 
103 
104 //-----------------------------------------------------------------------------
105 // Function : AnalysisBase::printStepHeader()
106 // Purpose : Prints out step information.
107 // Special Notes :
108 // Scope : public
109 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
110 // Creation Date : 6/26/00
111 //-----------------------------------------------------------------------------
112 void AnalysisBase::printStepHeader(std::ostream &os)
113 {}
114 
115 //-----------------------------------------------------------------------------
116 // Function : AnalysisBase::printProgress()
117 // Purpose : Outputs run completion percentage and estimated
118 // time-to-completion.
119 // Special Notes :
120 // Scope : public
121 // Creator : Scott A. Hutchinson, SNL, Computational Sciences
122 // Creation Date : 06/07/2002
123 //-----------------------------------------------------------------------------
124 void AnalysisBase::printProgress(std::ostream &os)
125 {}
126 
127 //-----------------------------------------------------------------------------
128 // Function : AnalysisBase::resetForStepAnalysis()
129 // Purpose : When doing a .STEP sweep, some data must be reset to its
130 // initial state.
131 // Special Notes :
132 // Scope : public
133 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
134 // Creation Date : 8/26/04
135 //-----------------------------------------------------------------------------
137 {
139  stepNumber = 0;
140  beginningIntegration = true;
141 
142  return true;
143 }
144 
145 
146 //-----------------------------------------------------------------------------
147 // Function : AnalysisBase::setupSweepLoop
148 // Purpose : Processes sweep parameters.
149 // Special Notes : Used for DC and STEP analysis classes.
150 // Scope : public
151 // Creator : Eric R. Keiter, SNL.
152 // Creation Date : 08/21/04
153 //-----------------------------------------------------------------------------
154 int AnalysisBase::setupSweepLoop_( std::vector<SweepParam> & sweepParamVec )
155 {
156  double pinterval = 1.0;
157  double pcount = 0.0, pstart, pstop, pstep;
158 
159  // loop over the param containers, and check that all the params exist.
160  // (the device package will complain if it can't find the param)
161  double tmp = 0.0;
162  for (int iparam=0;iparam<(int)sweepParamVec.size();++iparam)
163  {
164  tmp = loaderRCPtr_->getParamAndReduce(sweepParamVec[iparam].name);
165  }
166 
167 
168 #ifdef Xyce_DEBUG_ANALYSIS
169  if (anaManagerRCPtr_->tiaParams.debugLevel > 0)
170  {
171  Xyce::dout() << std::endl << std::endl
172  << Xyce::subsection_divider << std::endl
173  << "AnalysisManager::setupSweepLoop" << std::endl;
174  }
175 #endif
176 
177  // loop over the param containers:
178  for (int iparam=0;iparam<(int)sweepParamVec.size();++iparam)
179  {
180 #ifdef Xyce_DEBUG_ANALYSIS
181  if (anaManagerRCPtr_->tiaParams.debugLevel > 0)
182  {
183  Xyce::dout() << "name = " << sweepParamVec[iparam].name << std::endl;
184  }
185 #endif
186  // set interval:
187  sweepParamVec[iparam].interval = static_cast<int> (pinterval);
188 
189  // This stuff should probably be moved up into the SweepParam class.
190 
191  // obtain next pinterval:
192  if (sweepParamVec[iparam].type=="LIN")
193  {
194  pstart = sweepParamVec[iparam].startVal;
195  pstop = sweepParamVec[iparam].stopVal;
196  pstep = sweepParamVec[iparam].stepVal;
197  // ----------
198  // pcount = floor(((pstop - pstart)/pstep) + 1.0);
199  // The computation of "pcount" above is notoriously prone to roundoff
200  // error, especially if the "floor" function is an in-line function
201  // and is subject to high levels of optimization on x86 processors.
202  // The symptom is that the last step of a DC sweep (or other sweep)
203  // gets lost for very specific combinations of start/stop/step.
204  // The next few lines are an attempt to mitigate this roundoff issue,
205  // which was present in Xyce for years, and was inherited from SPICE3F5,
206  // from which the above expression was taken verbatim.
207 
208  // Compute the number of steps of size pstep between pstart and pstop
209  pcount = floor(((pstop - pstart)/pstep));
210  // Here we're checking that adding one more step doesn't pass pstop
211  // by more than machine precision. If we're within machine precision
212  // of pstop by taking one more step, that must mean we undercounted
213  // due to roundoff in the division --- if we hadn't undercounted, we'd
214  // exceed pstop by (nearly) a full pstep.
215  if ( fabs(pstop-(pstart+(pcount+1.0)*pstep)) < 2.0*N_UTL_MachineDependentParams::MachinePrecision())
216  {
217  pcount += 1.0;
218  }
219 
220  // Pcount is now the exact number of steps of size pstep between pstart
221  // and pstop, with roundoff handled mostly cleanly.
222 
223  // finally, because our actual loop does a loop from zero to maxStep-1,
224  // we have to pad maxStep (as is done in the original pcount expression
225  // above) to get the full range.
226  pcount += 1.0;
227 
228  // done this way, we should no longer miss final steps of DC sweeps.
229  // Fixed 31 Jul 2012. This was bug 695 in Bugzilla, and had plagued
230  // us since Xyce was first ported to Linux with GCC.
231  // ----------
232 
233  sweepParamVec[iparam].maxStep = static_cast<int>(pcount);
234 #ifdef Xyce_DEBUG_ANALYSIS
235  if (anaManagerRCPtr_->tiaParams.debugLevel > 0)
236  {
237  Xyce::dout() << "pstart = " << pstart << std::endl;
238  Xyce::dout() << "pstop = " << pstop << std::endl;
239  Xyce::dout() << "pstep = " << pstep << std::endl;
240  Xyce::dout() << "pstop-pstart/pstep = " << ((pstop - pstart)/pstep) << std::endl;
241  Xyce::dout() << "floor ()= " << floor(((pstop - pstart)/pstep)+1.0) << std::endl;
242  Xyce::dout() << "pcount = " << pcount << std::endl;
243  Xyce::dout() << "maxStep = " << sweepParamVec[iparam].maxStep << std::endl;
244  }
245 #endif
246  }
247  else if(sweepParamVec[iparam].type=="DEC")
248  {
249  double numSteps = static_cast<double>(sweepParamVec[iparam].numSteps);
250  // stepMult could also be calculated as pow(10,(1/numSteps))
251  double stepMult = exp(log(10.0)/numSteps);
252  sweepParamVec[iparam].stepMult = stepMult;
253 
254  pstart = sweepParamVec[iparam].startVal;
255  pstop = sweepParamVec[iparam].stopVal;
256  pcount = floor(fabs(log10(pstart) - log10(pstop)) * numSteps + 1);
257  sweepParamVec[iparam].maxStep = static_cast<int>(pcount);
258  }
259  else if(sweepParamVec[iparam].type=="OCT")
260  {
261  double numSteps = static_cast<double>(sweepParamVec[iparam].numSteps);
262  // stepMult could also be calculated as pow(2,1/(numSteps))
263  double stepMult = exp(log(2.0)/numSteps);
264 
265  // changed to remove dependence on "log2" function, which apparently
266  // doesn't exist in the math libraries of FreeBSD or the mingw
267  // cross-compilation suite. Log_2(x)=log_e(x)/log_e(2.0)
268  double ln2=log(2.0);
269 
270  sweepParamVec[iparam].stepMult = stepMult;
271  pstart = sweepParamVec[iparam].startVal;
272  pstop = sweepParamVec[iparam].stopVal;
273  pcount = floor(fabs(log(pstart) - log(pstop))/ln2 * numSteps + 1);
274  sweepParamVec[iparam].maxStep = static_cast<int>(pcount);
275  }
276  else if(sweepParamVec[iparam].type=="LIST")
277  {
278  pcount = sweepParamVec[iparam].valList.size();
279  sweepParamVec[iparam].maxStep = sweepParamVec[iparam].valList.size();
280  }
281  else
282  {
283  Report::UserError0() << " Unsupported STEP type";
284  }
285  pinterval *= pcount;
286 
287 #ifdef Xyce_DEBUG_ANALYSIS
288  if (anaManagerRCPtr_->tiaParams.debugLevel > 0)
289  {
290  Xyce::dout() << "parameter = " << sweepParamVec[iparam].name << std::endl;
291  Xyce::dout() << "pcount = " << pcount << std::endl;
292  Xyce::dout() << "pinterval = " << pinterval << std::endl;
293  }
294 #endif
295  }
296 
297  // At this point, pinterval equals the total number of steps
298  // for the step loop.
299  return static_cast<int> (pinterval);
300 }
301 
302 //-----------------------------------------------------------------------------
303 // Function : AnalysisBase::updateSweepParams_
304 // Purpose : Update parameters either for DC or STEP sweeps
305 // Special Notes :
306 // Scope : public
307 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
308 // Creation Date : 8/26/04
309 //-----------------------------------------------------------------------------
311  (int loopIter, std::vector<SweepParam> & sweepParamVec)
312 {
313  std::vector<SweepParam>::iterator iterParam;
314  std::vector<SweepParam>::iterator firstParam = sweepParamVec.begin();
315  std::vector<SweepParam>::iterator lastParam = sweepParamVec.end ();
316  bool resetFlag=false;
317 
318  // set parameter(s)
319  for (iterParam=firstParam; iterParam != lastParam;++iterParam)
320  {
321  iterParam->updateCurrentVal (loopIter);
322  resetFlag = resetFlag || iterParam->getSweepResetFlag();
323  loaderRCPtr_->setParam (iterParam->name, iterParam->currentVal);
324  }
325 
326  // Tell the manager if any of our sweeps are being reset in this loop
327  // iteration.
328  anaManagerRCPtr_->setSweepSourceResetFlag(resetFlag);
329 
330  return true;
331 }
332 
333 //-----------------------------------------------------------------------------
334 // Function : AnalysisBase::resetAll
335 // Purpose :
336 // Special Notes :
337 // Scope : public
338 // Creator : Eric R. Keiter
339 // Creation Date : 12/16/10
340 //-----------------------------------------------------------------------------
342 {
343  stepNumber = 0;
344  tranStepNumber = 0;
346 
359 }
360 
361 //-----------------------------------------------------------------------------
362 // Function : AnalysisBase::saveLoopInfo
363 // Purpose :
364 // Special Notes :
365 // Scope : public
366 // Creator : Eric R. Keiter
367 // Creation Date : 12/16/10
368 //-----------------------------------------------------------------------------
370 {
371  int pos = saveTimeI.size();
372  saveTimeI.resize(pos+1);
373  saveTimeD.resize(pos+1);
374 
379  saveTimeI[pos].push_back(totalNumberLinearSolves_);
381  saveTimeI[pos].push_back(totalNumberLinearIters_);
384  saveTimeD[pos].push_back(totalResidualLoadTime_);
385  saveTimeD[pos].push_back(totalJacobianLoadTime_);
386  saveTimeD[pos].push_back(totalLinearSolutionTime_);
387 
388  if (pos == 0)
389  {
390  int i,j;
391 
392  pos++;
393  saveTimeI.resize(pos+1);
394  saveTimeD.resize(pos+1);
395 
396  j = saveTimeI[0].size();
397  saveTimeI[1].resize(j);
398  for (i=0 ; i<j ; ++i)
399  {
400  saveTimeI[1][i] = saveTimeI[0][i];
401  saveTimeI[0][i] = 0;
402  }
403  j = saveTimeD[0].size();
404  saveTimeD[1].resize(j);
405  for (i=0 ; i<j ; ++i)
406  {
407  saveTimeD[1][i] = saveTimeD[0][i];
408  saveTimeD[0][i] = 0;
409  }
410  }
411 
412  return pos;
413 }
414 
415 //-----------------------------------------------------------------------------
416 // Function : AnalysisBase::printLoopInfo
417 // Purpose : Prints out time loop information.
418 // Special Notes : Prints stats from save point start to save point finish.
419 // Special case 0,0 is entire run to this point
420 // Scope : public
421 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
422 // Creation Date : 6/26/00
423 //-----------------------------------------------------------------------------
424 bool AnalysisBase::printLoopInfo(int start, int finish)
425 {
426  bool bsuccess = true;
427  int t1, t2, indI = 0, indD = 0;
428  if (start == 0 && finish == 0)
429  {
430  t2 = saveLoopInfo();
431  t1 = 0;
432  }
433  else
434  {
435  t1 = start;
436  t2 = finish;
437  }
438 
439  lout() << "\tNumber Successful Steps Taken:\t\t" << saveTimeI[t2][indI]-saveTimeI[t1][indI] << std::endl;
440  indI++;
441 
442  lout() << "\tNumber Failed Steps Attempted:\t\t" << saveTimeI[t2][indI]-saveTimeI[t1][indI] << std::endl;
443  indI++;
444 
445  lout() << "\tNumber Jacobians Evaluated:\t\t" << saveTimeI[t2][indI]-saveTimeI[t1][indI] << std::endl;
446  indI++;
447 
448  lout() << "\tNumber Iteration Matrix Factorizations:\t" << saveTimeI[t2][indI]-saveTimeI[t1][indI] << std::endl;
449  indI++;
450 
451  lout() << "\tNumber Linear Solves:\t\t\t" << saveTimeI[t2][indI]-saveTimeI[t1][indI] << std::endl;
452  indI++;
453 
454  lout() << "\tNumber Failed Linear Solves:\t\t" << saveTimeI[t2][indI]-saveTimeI[t1][indI] << std::endl;
455  indI++;
456 
457  if (saveTimeI[t2][indI] > saveTimeI[t1][indI])
458  {
459  lout() << "\tNumber Linear Solver Iterations:\t" << saveTimeI[t2][indI]-saveTimeI[t1][indI] << std::endl;
460  }
461  indI++; // still needs to be incremented, even if info is not printed
462 
463  lout() << "\tNumber Residual Evaluations:\t\t" << saveTimeI[t2][indI]-saveTimeI[t1][indI] << std::endl;
464  indI++;
465 
466  lout() << "\tNumber Nonlinear Convergence Failures:\t" << saveTimeI[t2][indI]-saveTimeI[t1][indI] << std::endl;
467  indI++;
468 
469  lout() << std::endl;
470 
471  lout() << "\tTotal Residual Load Time:\t\t" << saveTimeD[t2][indD]-saveTimeD[t1][indD] << " seconds" << std::endl;
472  indD++;
473 
474  lout() << "\tTotal Jacobian Load Time:\t\t" << saveTimeD[t2][indD]-saveTimeD[t1][indD] << " seconds" << std::endl;
475  indD++;
476 
477  lout() << "\tTotal Linear Solution Time:\t\t" << saveTimeD[t2][indD]-saveTimeD[t1][indD] << " seconds" << std::endl;
478  indD++;
479 
480  return bsuccess;
481 }
482 
484 {
485  if (secRCPtr_->newtonConvergenceStatus <= 0)
486  {
488  }
489 
490  totalNumberJacobiansEvaluated_ += nlsMgrRCPtr_->getNumJacobianLoads();
491  totalNumberLinearSolves_ += nlsMgrRCPtr_->getNumLinearSolves();
492  totalNumberFailedLinearSolves_ += nlsMgrRCPtr_->getNumFailedLinearSolves();
493  totalNumberLinearIters_ += nlsMgrRCPtr_->getTotalNumLinearIters();
494  totalNumberResidualEvaluations_ += nlsMgrRCPtr_->getNumResidualLoads();
495  totalNumberIterationMatrixFactorizations_ += nlsMgrRCPtr_->getNumJacobianFactorizations();
496  totalLinearSolutionTime_ += nlsMgrRCPtr_->getTotalLinearSolveTime();
497  totalResidualLoadTime_ += nlsMgrRCPtr_->getTotalResidualLoadTime();
498  totalJacobianLoadTime_ += nlsMgrRCPtr_->getTotalJacobianLoadTime();
499 }
500 
501 } // namespace Analysis
502 } // namespace Xyce
503