Xyce  6.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
N_TIA_OneStep.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_TIA_OneStep.C,v $
27 //
28 // Purpose : This file contains the functions which define the
29 // Trap method, order 1-2, class.
30 //
31 // Special Notes :
32 //
33 // Creator : Ting Mei
34 //
35 // Creation Date : 2/16/04
36 //
37 // Revision Information:
38 // ---------------------
39 //
40 // Revision Number: $Revision: 1.69.2.6 $
41 //
42 // Revision Date : $Date: 2014/09/03 11:14:50 $
43 //
44 // Current Owner : $Author: erkeite $
45 //-----------------------------------------------------------------------------
46 
47 #include <Xyce_config.h>
48 
49 // ---------- Standard Includes ----------
50 #include <N_UTL_Misc.h>
51 
52 #ifdef HAVE_IOSTREAM
53 #include <iostream>
54 #else
55 #include <iostream.h>
56 #endif
57 
58 // ---------- Xyce Includes ----------
59 #include <N_TIA_OneStep.h>
60 
61 #include <N_ANP_OutputMgrAdapter.h>
62 #include <N_ERH_ErrorMgr.h>
63 #include <N_LAS_BlockVector.h>
64 #include <N_LAS_Matrix.h>
65 #include <N_LAS_System.h>
66 #include <N_LAS_Vector.h>
67 #include <N_PDS_Comm.h>
68 #include <N_PDS_Manager.h>
69 #include <N_TIA_DataStore.h>
70 #include <N_TIA_StepErrorControl.h>
71 #include <N_TIA_TIAParams.h>
73 
74 using std::min;
75 using std::max;
76 using std::abs;
77 
78 //-----------------------------------------------------------------------------
79 // Function : N_TIA_OneStep::N_TIA_OneStep
80 // Purpose : constructor
81 // Special Notes :
82 // Scope : public
83 // Creator : Ting Mei, SNL
84 // Creation Date : 11/16/07
85 //-----------------------------------------------------------------------------
88  N_TIA_StepErrorControl & secTmp,
89  N_TIA_DataStore & dsTmp)
90 : N_TIA_TimeIntegrationMethod(tiaP,secTmp,dsTmp),
91  timeStepForHistory2_(0.0)
92 {
93  leadingCoeff = 1;
94  sec.maxOrder_=(Xycemin(2,tiaParams.maxOrder));
95  sec.minOrder_=(Xycemax(1,tiaParams.minOrder));
96 
97  if (sec.minOrder_ > sec.maxOrder_)
98  {
99  sec.minOrder_ = sec.maxOrder_;
100  }
101  // sec.maxOrder_ = 2;
102  timept_ = -1.0;
103  return ;
104 }
105 
106 //-----------------------------------------------------------------------------
107 // Function : N_TIA_OneStep::obtainPredictor
108 // Purpose : Calculate predictor
109 // Special Notes : stored in ds.xn0Ptr,qn0Ptr,qpn0Ptr
110 // Scope : public
111 // Creator : Ting Mei, SNL
112 // Creation Date : 11/16/07
113 //-----------------------------------------------------------------------------
115 {
116  // evaluate predictor
117  *ds.xn0Ptr = *(ds.xHistory[0]);
118  *ds.qn0Ptr = *(ds.qHistory[0]);
119  *ds.sn0Ptr = *(ds.sHistory[0]);
120  *ds.ston0Ptr = *(ds.stoHistory[0]);
122 
123  ds.spn0Ptr->putScalar(0.0);
124  ds.stopn0Ptr->putScalar(0.0); // need to verify that these are needed
125  ds.stoQpn0Ptr->putScalar(0.0); // need to verify that these are needed
126  for (int i=1;i<=sec.currentOrder_;++i)
127  {
128  ds.xn0Ptr->linearCombo(sec.beta_[i],*(ds.xHistory[i]),1.0,*ds.xn0Ptr);
129  }
130 
131 #ifdef Xyce_DEBUG_TIME
132  Xyce::dout().width(21); Xyce::dout().precision(13); Xyce::dout().setf(std::ios::scientific);
133 
134  if (tiaParams.debugLevel > 1)
135  {
136  Xyce::dout() << std::endl;
137  Xyce::dout() << Xyce::section_divider << std::endl;
138  Xyce::dout() <<
139  " N_TIA_OneStep::obtainPredictor" << std::endl;
140  Xyce::dout() << "\n currentOrder = " << sec.currentOrder_ << std::endl;
141  Xyce::dout() << "\n sec.nscsco_: " << sec.nscsco_ << std::endl;
142  for (int i=0; i<=sec.currentOrder_ ; ++i)
143  Xyce::dout() << "\n sec.beta_[" << i << "] = " << sec.beta_[i] << "\n" << std::endl;
144  for (int i=0; i<=sec.currentOrder_ ; ++i)
145  {
146  Xyce::dout() << "\n xHistory["<< i << "]: \n" << std::endl;
147  (ds.xHistory[i])->printPetraObject(Xyce::dout());
148  Xyce::dout() << std::endl;
149  }
150  for (int i=0; i<=sec.currentOrder_ ; ++i)
151  {
152  Xyce::dout() << "\n qHistory["<< i << "]: \n" << std::endl;
153  (ds.qHistory[i])->printPetraObject(Xyce::dout());
154  Xyce::dout() << std::endl;
155  }
156  for (int i=0; i<=sec.currentOrder_ ; ++i)
157  {
158  Xyce::dout() << "\n sHistory["<< i << "]: \n" << std::endl;
159  (ds.sHistory[i])->printPetraObject(Xyce::dout());
160  Xyce::dout() << std::endl;
161  }
162  for (int i=0; i<=sec.currentOrder_ ; ++i)
163  {
164  Xyce::dout() << "\n stoHistory["<< i << "]: \n" << std::endl;
165  (ds.stoHistory[i])->printPetraObject(Xyce::dout());
166  Xyce::dout() << std::endl;
167  }
168  for (int i=0; i<=sec.currentOrder_ ; ++i)
169  {
170  Xyce::dout() << "\n stoLeadCurrQHistory["<< i << "]: \n" << std::endl;
171  (ds.stoLeadCurrQHistory[i])->printPetraObject(Xyce::dout());
172  Xyce::dout() << std::endl;
173  }
174  Xyce::dout() << "\n xn0: \n" << std::endl;
175  ds.xn0Ptr->printPetraObject(Xyce::dout());
176  Xyce::dout() << std::endl;
177  Xyce::dout() << "\n qn0: \n" << std::endl;
178  ds.qn0Ptr->printPetraObject(Xyce::dout());
179  Xyce::dout() << std::endl;
180  Xyce::dout() << "\n qpn0: \n" << std::endl;
181  ds.qpn0Ptr->printPetraObject(Xyce::dout());
182  Xyce::dout() << std::endl;
183  Xyce::dout() << "\n sn0: \n" << std::endl;
184  ds.sn0Ptr->printPetraObject(Xyce::dout());
185  Xyce::dout() << std::endl;
186  Xyce::dout() << "\n spn0: \n" << std::endl;
187  ds.spn0Ptr->printPetraObject(Xyce::dout());
188  Xyce::dout() << std::endl;
189  Xyce::dout() << "\n stoQn0Ptr: \n" << std::endl;
190  ds.stoQn0Ptr->printPetraObject(Xyce::dout());
191  Xyce::dout() << std::endl;
192  Xyce::dout() << "\n stoQpn0Ptr: \n" << std::endl;
193  ds.stoQpn0Ptr->printPetraObject(Xyce::dout());
194  Xyce::dout() << std::endl;
195  Xyce::dout() << Xyce::section_divider << std::endl;
196  }
197 #endif // Xyce_DEBUG_TIME
198  // copy the prediction into the next solution:
199  *(ds.nextSolutionPtr) = *(ds.xn0Ptr);
200 
202 
203  return;
204 }
205 
206 //-----------------------------------------------------------------------------
207 // Function : N_TIA_OneStep::obtainSensitivityPredictors
208 // Purpose : Calculate predictor
209 // Special Notes :
210 // Scope : public
211 // Creator : Eric Keiter, SNL
212 // Creation Date :
213 //-----------------------------------------------------------------------------
215 {
216  // commenting out b/c for sensitivities these aren't needed
217 #if 1
218  int numParams = ds.sensRHSPtrVector.size();
219  for (int ip=0; ip<numParams;++ip)
220  {
221  std::vector<N_LAS_Vector*> & dfdpHistory = ds.dfdpHistory[ip];
222  std::vector<N_LAS_Vector*> & dqdpHistory = ds.dqdpHistory[ip];
223  std::vector<N_LAS_Vector*> & dbdpHistory = ds.dbdpHistory[ip];
224  std::vector<N_LAS_Vector*> & dXdpHistory = ds.dXdpHistory[ip];
225 
226  N_LAS_Vector* dfdpn0Ptr = ds.dfdpn0PtrVector[ip];
227  N_LAS_Vector* dqdpn0Ptr = ds.dqdpn0PtrVector[ip];
228  N_LAS_Vector* dbdpn0Ptr = ds.dbdpn0PtrVector[ip];
229  N_LAS_Vector* dXdpn0Ptr = ds.dXdpn0PtrVector[ip];
230 
231  N_LAS_Vector* dfdppn0Ptr = ds.dfdppn0PtrVector[ip];
232  N_LAS_Vector* dqdppn0Ptr = ds.dqdppn0PtrVector[ip];
233  N_LAS_Vector* dbdppn0Ptr = ds.dbdppn0PtrVector[ip];
234  N_LAS_Vector* dXdppn0Ptr = ds.dXdppn0PtrVector[ip];
235 
236  // evaluate predictor
237  *dfdpn0Ptr = *(dfdpHistory[0]);
238  *dqdpn0Ptr = *(dqdpHistory[0]);
239  *dbdpn0Ptr = *(dbdpHistory[0]);
240  *dXdpn0Ptr = *(dXdpHistory[0]);
241 
242  // not used, only for BDF15 integrator.
243  dfdppn0Ptr->putScalar(0.0);
244  dqdppn0Ptr->putScalar(0.0);
245  dbdppn0Ptr->putScalar(0.0);
246  dXdppn0Ptr->putScalar(0.0);
247 
248 #ifdef DEBUG_SENS
249  Xyce::dout() << "obtainSensitivityPredictors: dqdpn0Ptr, param = " << ip << "\n";
250  dqdpn0Ptr->printPetraObject(Xyce::dout());
251 #endif
252 
253  // this probably isn't needed, but following the pattern.
254  for (int i=1;i<=sec.currentOrder_;++i)
255  {
256  dXdpn0Ptr->linearCombo(sec.beta_[i],*(dXdpHistory[i]),1.0,*dXdpn0Ptr);
257  }
258  }
259 #endif
260 
261  return;
262 }
263 
264 //-----------------------------------------------------------------------------
265 // Function : N_TIA_OneStep::obtainResidual
266 // Purpose : Calculate Residual
267 // Special Notes :
268 // Scope : public
269 // Creator : Ting Mei, SNL
270 // Creation Date : 11/16/07
271 //-----------------------------------------------------------------------------
273 {
274  // output: ds.RHSVectorPtr
275  // Note: ds.nextSolutionPtr is used to get Q,F,B in N_ANP_AnalysisManager::loadRHS.
276  ds.RHSVectorPtr->linearCombo(1.0,*ds.daeQVectorPtr,-1.0,*ds.qn0Ptr);
277 
278 #ifdef Xyce_DEBUG_TIME
279  if (tiaParams.debugLevel > 1)
280  {
281  Xyce::dout() << std::endl;
282  Xyce::dout() << Xyce::section_divider << std::endl;
283  Xyce::dout() <<
284  " N_TIA_OneStep::obtainResidual" << std::endl;
285  Xyce::dout() << "\n t = " << sec.nextTime << "\n" << std::endl;
286  Xyce::dout() << "\n solution: \n" << std::endl;
287  ds.nextSolutionPtr->printPetraObject(Xyce::dout());
288  Xyce::dout() << "\n daeQVector: \n" << std::endl;
289  ds.daeQVectorPtr->printPetraObject(Xyce::dout());
290  Xyce::dout() << "\n qn0: \n" << std::endl;
291  ds.qn0Ptr->printPetraObject(Xyce::dout());
292  Xyce::dout() << "\n qpn0: \n" << std::endl;
293  ds.qpn0Ptr->printPetraObject(Xyce::dout());
294  Xyce::dout() << "\n sec.alphas_/hn: " << sec.alphas_/sec.currentTimeStep << "\n" << std::endl;
295  Xyce::dout() << "\n daeFVector: \n" << std::endl;
296  ds.daeFVectorPtr->printPetraObject(Xyce::dout());
297 
298  Xyce::dout() << "\n dQdt-vector: \n" << std::endl;
299  ds.RHSVectorPtr->printPetraObject(Xyce::dout());
300  Xyce::dout() << std::endl;
301  }
302 #endif
303 
304  if (sec.currentOrder_ == 2)
305  {
306  ds.RHSVectorPtr->linearCombo(1.0/sec.currentTimeStep,*ds.RHSVectorPtr,+1.0/2.0,*ds.daeFVectorPtr);
307 #ifdef SEPARATE_F_AND_B
308  ds.RHSVectorPtr->linearCombo(1.0,*ds.RHSVectorPtr,-1.0/2.0,*ds.daeBVectorPtr);
309 #endif
310  ds.RHSVectorPtr->linearCombo(1.0,*ds.RHSVectorPtr,+1.0/2.0,*ds.qHistory[2]);
311  }
312  else
313  {
315 #ifdef SEPARATE_F_AND_B
316  ds.RHSVectorPtr->linearCombo(1.0,*ds.RHSVectorPtr,-1.0,*ds.daeBVectorPtr);
317 #endif
318  }
319  // since the nonlinear solver is expecting a -f, scale by -1.0:
320  ds.RHSVectorPtr->scale(-1.0);
321 
322  // if voltage limiting is on, add it in:
323  if (ds.limiterFlag)
324  {
326 
327  (ds.RHSVectorPtr)->daxpy(
328  *(ds.RHSVectorPtr), +1.0, *(ds.dQdxdVpVectorPtr));
329 
330  double fscalar(1.0);
331 
332  if (sec.currentOrder_ == 2)
333  fscalar =1.0/2.0;
334 
335  (ds.RHSVectorPtr)->daxpy(
336  *(ds.RHSVectorPtr), fscalar, *(ds.dFdxdVpVectorPtr));
337  }
338 
339 #ifdef Xyce_DEBUG_TIME
340  if (tiaParams.debugLevel > 1)
341  {
342  Xyce::dout() << "\n Residual-vector: \n" << std::endl;
343  Xyce::dout() << "-(qpn0-(sec.alpha_s/h)*(Q-qn0)+F-B) \n" << std::endl;
344  ds.RHSVectorPtr->printPetraObject(Xyce::dout());
345  Xyce::dout() << Xyce::section_divider << std::endl;
346  Xyce::dout() << std::endl;
347  }
348 #endif
349 }
350 
351 //-----------------------------------------------------------------------------
352 // Function : N_TIA_OneStep::obtainSensitivityResiduals
353 // Purpose : Calculate Sensitivity Residual
354 // Special Notes :
355 // Scope : public
356 // Creator : Eric Keiter, SNL
357 // Creation Date : 09/02/2014
358 //-----------------------------------------------------------------------------
360 {
361  int numParams = ds.sensRHSPtrVector.size();
362  for (int ip=0; ip<numParams;++ip)
363  {
364  N_LAS_Vector & RHSVec = *(ds.sensRHSPtrVector[ip]);
365  N_LAS_Vector & dfdpVec = *(ds.nextDfdpPtrVector[ip]);
366  N_LAS_Vector & dqdpVec = *(ds.nextDqdpPtrVector[ip]);
367  N_LAS_Vector & dbdpVec = *(ds.nextDbdpPtrVector[ip]);
368 
369  N_LAS_Vector & currDXdpVec = *(ds.currDXdpPtrVector[ip]);
370 
371  N_LAS_Vector & dqdpn0Vec = *(ds.dqdpn0PtrVector[ip]);
372  std::vector<N_LAS_Vector*> & dqdpHistory = ds.dqdpHistory[ip];
373 
374  N_LAS_Vector & currDQdxDXdpVec = *(ds.currDQdxDXdpPtrVector[ip]);
375 
376  N_LAS_Vector & nextDQdxDXdpDerivVec = *(ds.nextDQdxDXdpDerivPtrVector[ip]);
377  N_LAS_Vector & currDQdxDXdpDerivVec = *(ds.currDQdxDXdpDerivPtrVector[ip]);
378 
379  N_LAS_Vector & nextDqdpDerivVec = *(ds.nextDqdpDerivPtrVector[ip]);
380  N_LAS_Vector & currDqdpDerivVec = *(ds.currDqdpDerivPtrVector[ip]);
381 
382 
383  if (sec.currentOrder_ == 2)
384  {
385  nextDqdpDerivVec.linearCombo(1.0,dqdpVec,-1.0, dqdpn0Vec);
386  nextDqdpDerivVec.scale(2.0/sec.currentTimeStep);
387  nextDqdpDerivVec.linearCombo(-1.0,currDqdpDerivVec,1.0,nextDqdpDerivVec);
388 
389  RHSVec.linearCombo(0.5,nextDqdpDerivVec,+0.5,dfdpVec);
390 #ifdef SEPARATE_F_AND_B
391  RHSVec.linearCombo(1.0,RHSVec,-0.5,dbdpVec);
392 #endif
393  }
394  else
395  {
396  nextDqdpDerivVec.linearCombo(1.0,dqdpVec,-1.0, dqdpn0Vec);
397  nextDqdpDerivVec.scale(1.0/sec.currentTimeStep);
398 
399  RHSVec.linearCombo(1.0,nextDqdpDerivVec,+1.0,dfdpVec);
400 #ifdef SEPARATE_F_AND_B
401  RHSVec.linearCombo(1.0,RHSVec,-1.0,dbdpVec);
402 #endif
403  }
404 
405  // since the nonlinear solver is expecting a -dFdp, scale by -1.0:
406  RHSVec.scale(-1.0);
407 
408  // correction terms
409  double qscalar(-sec.alphas_/sec.currentTimeStep);
410  RHSVec.linearCombo(1.0,RHSVec, qscalar, currDQdxDXdpVec);
411 
412  // second order "correction" term
413  if (sec.currentOrder_ == 2)
414  {
415  RHSVec.linearCombo(1.0,RHSVec, 0.5, currDQdxDXdpDerivVec);
416  }
417 
418 #ifdef DEBUG_SENS
419  Xyce::dout() << "obtainSensitivityResiduals: RHS Vector, ip = " << ip << ":\n";
420  RHSVec.printPetraObject(Xyce::dout());
421 #endif
422  }
423 }
424 
425 //-----------------------------------------------------------------------------
426 // Function : N_TIA_OneStep::loadFinalSensitivityDerivatives
427 // Purpose : obtain final derivatives w.r.t. time for transient
428 // sensitivity analysis
429 // Special Notes :
430 // Scope : public
431 // Creator : Eric Keiter, SNL
432 // Creation Date : 09/02/2014
433 //-----------------------------------------------------------------------------
435 {
436  int numParams = ds.sensRHSPtrVector.size();
437  for (int ip=0; ip<numParams;++ip)
438  {
439  N_LAS_Vector & nextDQdxDXdpVec = *(ds.nextDQdxDXdpPtrVector[ip]);
440  N_LAS_Vector & currDQdxDXdpVec = *(ds.currDQdxDXdpPtrVector[ip]);
441 
442  N_LAS_Vector & nextDQdxDXdpDerivVec = *(ds.nextDQdxDXdpDerivPtrVector[ip]);
443  N_LAS_Vector & currDQdxDXdpDerivVec = *(ds.currDQdxDXdpDerivPtrVector[ip]);
444 
445  double dt = sec.currentTimeStep;
446  if (sec.currentOrder_ == 2)
447  {
448  nextDQdxDXdpDerivVec.putScalar(0.0);
449  nextDQdxDXdpDerivVec.linearCombo(2.0,nextDQdxDXdpVec,-2.0,currDQdxDXdpVec);
450  nextDQdxDXdpDerivVec.scale(1.0/dt);
451  nextDQdxDXdpDerivVec.linearCombo(-1.0,currDQdxDXdpDerivVec,1.0,nextDQdxDXdpDerivVec);
452  }
453  else
454  {
455  nextDQdxDXdpDerivVec.putScalar(0.0);
456  nextDQdxDXdpDerivVec.linearCombo(1.0,nextDQdxDXdpVec,-1.0,currDQdxDXdpVec);
457  nextDQdxDXdpDerivVec.scale(1.0/dt);
458  }
459  }
460 }
461 
462 //-----------------------------------------------------------------------------
463 // Function : N_TIA_OneStep::obtainJacobian
464 // Purpose : Calculate Jacobian
465 // Special Notes :
466 // Scope : public
467 // Creator : Ting Mei, SNL
468 // Creation Date : 11/16/07
469 //-----------------------------------------------------------------------------
471 {
472 #ifdef Xyce_DEBUG_TIME
473  if (tiaParams.debugLevel > 1)
474  {
475  Xyce::dout() << std::endl;
476  Xyce::dout() << Xyce::section_divider << std::endl;
477  Xyce::dout() <<
478  " N_TIA_OneStep::obtainJacobian" << std::endl;
479  }
480 #endif
481  // output: ds.JMatrixPtr
482 
483  // This function returns the following matrix:
484  // $-(sec.alphas_/hn)dQdx(x)+dFdx$
485 
486  // Note: ds.nextSolutionPtr is used to get dQdx,dFdx in N_ANP_AnalysisManager::loadJacobian.
487 
488  N_LAS_Matrix & dQdx = *(ds.dQdxMatrixPtr);
489  N_LAS_Matrix & dFdx = *(ds.dFdxMatrixPtr);
490  N_LAS_Matrix & Jac = *(ds.JMatrixPtr);
491 
492  double qscalar(-sec.alphas_/sec.currentTimeStep);
493  double fscalar(1.0);
494  if (sec.currentOrder_ == 2)
495  fscalar =1.0/2.0;
496 
497  Jac.linearCombo( qscalar, dQdx, fscalar, dFdx );
498 
499 #ifdef Xyce_DEBUG_TIME
500  if (tiaParams.debugLevel > 1)
501  {
502  Xyce::dout() << "fscalar = " << fscalar << " qscalar = " << qscalar << std::endl;
503 
504  Xyce::dout() << "\n dFdx:" <<std::endl;
505  dFdx.printPetraObject(Xyce::dout());
506  Xyce::dout() << "\n dQdx:" <<std::endl;
507  dQdx.printPetraObject(Xyce::dout());
508  Xyce::dout() << "\n Total Jacobian:" <<std::endl;
509  Jac.printPetraObject(Xyce::dout());
510  // for (int i=0;i<3;++i)
511  // {
512  // printf("[ %25.20g\t%25.20g\t%25.20g ]\n",Jac[i][0],Jac[i][1],Jac[i][2]);
513  // }
514 
515  Xyce::dout() << Xyce::section_divider << std::endl;
516  Xyce::dout() << std::endl;
517  }
518 #endif
519 
520 }
521 
522 //-----------------------------------------------------------------------------
523 // Function : N_TIA_OneStep::interpolateSolution
524 // Purpose : Interpolate solution approximation at prescribed time point.
525 // Scope : public
526 // Creator : Ting Mei, SNL
527 // Creation Date : 11/16/07
528 //-----------------------------------------------------------------------------
529 bool N_TIA_OneStep::interpolateSolution(double timepoint,
530  N_LAS_Vector * tmpSolVectorPtr, std::vector<N_LAS_Vector*> & historyVec)
531 {
532  // this is a very course approximation to determine if we are too
533  // close to the actual time step to do an interpolation.
534  // it could use more work.
535  double dtr = timepoint - sec.currentTime; // the delta time requested.
536  *tmpSolVectorPtr = *(historyVec[0]);
537 
538  if( -dtr < 100 * N_UTL_MachineDependentParams::MachinePrecision() )
539  {
540  return false;
541  }
542 
543  if( sec.usedOrder_ <= 2)
544  {
545  // do first order interpolation
546  // X_interp = X + delta_t_requested * delta_X/delta_t[last step]
547  dtr = dtr / sec.lastTimeStep;
548  tmpSolVectorPtr->linearCombo(1.0,*tmpSolVectorPtr,dtr,*(historyVec[1]));
549  }
550 
551  return true;
552 }
553 
554 //-----------------------------------------------------------------------------
555 // Function : N_TIA_OneStep::interpolateMPDESolution
556 // Purpose : Interpolate solution approximation at prescribed time points.
557 // Special Notes : This routine computes the solution at the output
558 // : timepoints by intepolation of the history using the order
559 // : used for the most recent completed step, orderUsed.
560 // : The output is put into provided N_LAS_Vector pointer.
561 // : The interpolation is as follows:
562 // : tmpSolVectorPtr->block(i) is interpolated at timepoint(i)
563 // : Therefore, if you want them all interpolated at the same time,
564 // : then use timepoint(i) = timepoint(0) forall i
565 // : or use interpolateSolution.
566 // Scope : public
567 // Creator : Ting Mei, Eric Keiter, SNL
568 // Creation Date : 11/28/06
569 //-----------------------------------------------------------------------------
570 bool N_TIA_OneStep::interpolateMPDESolution(std::vector<double>& timepoint,
571  N_LAS_Vector * tmpSolVectorPtr)
572 {
573  N_LAS_BlockVector * blockTempSolVectorPtr =
574  dynamic_cast<N_LAS_BlockVector*>(tmpSolVectorPtr);
575  if (blockTempSolVectorPtr == NULL)
576  {
577  std::string msg = "N_TIA_OneStep::interpolateMPDESolution: ";
578  msg += "N_LAS_Vector tmpSolVectorPtr is not of type N_LAS_BlockVector";
579  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
580  return(false);
581  }
582 
583  double tfuzz; // fuzz factor to check for valid output time
584  double tp; // approximately t_{n-1}
585  int numblocks = timepoint.size();
586  int blockCount = blockTempSolVectorPtr->blockCount();
587  if (numblocks > blockCount)
588  {
589  std::string msg = "N_TIA_OneStep::interpolateMPDESolution: ";
590  msg += "Number of time points requested is greater than number of fast time points in MPDE block vector";
591  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
592  return(false);
593  }
594  double delt;
595  double c = 1.0;
596  double gam;
597  int kord; // order of interpolation
598  double tn = sec.currentTime;
599  double hh = sec.currentTimeStep;
600  double hused = sec.usedStep_;
601  int kused = sec.usedOrder_;
602  double uround = 0.0; // unit round-off (set to zero for now)
603 
604  tfuzz = 100 * uround * (tn + hh);
605  tp = tn - hused - tfuzz;
606  for (int i=0; i<numblocks ; ++i)
607  {
608  if ( (timepoint[i] - tp)*hh < 0.0 )
609  return false;
610  }
611 
612  *tmpSolVectorPtr = *(ds.xHistory[0]);
613 
614  N_LAS_Vector * solVectorPtr;
615  N_LAS_Vector * xHistoryVectorPtr;
616  // Loop over blocks first so that maximal order can be maintained
617  for (int i=0; i < numblocks ; ++i)
618  {
619  if ((kused == 0) || (timepoint[i] == tn)) { kord = 1; }
620  else { kord = kused; }
621  solVectorPtr = &(blockTempSolVectorPtr->block(i));
622  c = 1.0;
623  delt = timepoint[i] - tn;
624  gam = delt/sec.psi_[0];
625  for (int j=1 ; j <= kord ; ++j)
626  {
627  c = c*gam;
628  gam = (delt + sec.psi_[j-1])/sec.psi_[j];
629  N_LAS_BlockVector * blockXHistoryVectorPtr =
630  dynamic_cast<N_LAS_BlockVector*>(ds.xHistory[j]);
631  if (blockXHistoryVectorPtr == NULL)
632  {
633  Xyce::Report::DevelFatal0().in("N_TIA_OneStep::interpolateMPDESolution") << "N_LAS_Vector ds.xHistory[j] is not of type N_LAS_BlockVector\n j = " << j;
634  return(false);
635  }
636  xHistoryVectorPtr = &(blockXHistoryVectorPtr->block(i));
637  solVectorPtr->linearCombo(1.0,*solVectorPtr,c,*xHistoryVectorPtr);
638  }
639  }
640  return true;
641 }
642 
643 //-----------------------------------------------------------------------------
644 // Function : N_TIA_OneStep::printMPDEOutputSolution()
645 // Purpose : Print transient output from MPDE simulation
646 // Special Notes : This routine uses interpolateMPDESolution.
647 // Scope : public
648 // Creator : Ting Mei, SNL
649 // Creation Date : 11/28/06
650 //-----------------------------------------------------------------------------
652  N_ANP_OutputMgrAdapter & outputManagerAdapter,
653  const double time,
654  N_LAS_Vector * solnVecPtr,
655  const std::vector<double> & fastTimes )
656 {
657 #ifdef Xyce_DEBUG_TIME
658  if (tiaParams.debugLevel > 0)
659  {
660  Xyce::dout() << std::endl;
661  Xyce::dout() << Xyce::section_divider << std::endl;
662  Xyce::dout() <<
663  " N_TIA_OneStep::printMPDEOutputSolution" << std::endl;
664  }
665 #endif // Xyce_DEBUG_TIME
666  double timestep = sec.lastAttemptedTimeStep;
667  double lasttime = sec.currentTime - timestep;
668  double tn = sec.currentTime;
669  // Set these values up to read output time intervals. FIXME
670  double beg_of_output_time_interval = lasttime;
671  double end_of_output_time_interval = tn;
672  double start_time = max(lasttime,beg_of_output_time_interval);
673  double stop_time = min(tn,end_of_output_time_interval);
674 #ifdef Xyce_DEBUG_TIME
675  if (tiaParams.debugLevel > 0)
676  {
677  Xyce::dout() << "timestep = " << timestep << std::endl;
678  Xyce::dout() << "lasttime = " << lasttime << std::endl;
679  Xyce::dout() << "tn = " << tn << std::endl;
680  Xyce::dout() << "beg_of_output_time_interval = " << beg_of_output_time_interval << std::endl;
681  Xyce::dout() << "end_of_output_time_interval = " << end_of_output_time_interval << std::endl;
682  Xyce::dout() << "start_time = " << start_time << std::endl;
683  Xyce::dout() << "stop_time = " << stop_time << std::endl;
684  }
685 #endif // Xyce_DEBUG_TIME
686 
687 
688  N_LAS_BlockVector * blockTmpSolVectorPtr =
689  dynamic_cast<N_LAS_BlockVector*>(ds.tmpSolVectorPtr);
690  if (blockTmpSolVectorPtr == NULL)
691  {
692  std::string msg = "N_TIA_OneStep::printMPDEOutputSolution: ";
693  msg += "N_LAS_Vector ds.tmpSolVectorPtr is not of type N_LAS_BlockVector";
694  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
695  return(false);
696  }
697  int blockCount = blockTmpSolVectorPtr->blockCount();
698 
699  // Create list of timepoints to interpolate (along characteristic curve)
700  double T2 = fastTimes.back();
701  //double charcross = start_time - floor(start_time/T2)*T2; // (start_time mod T2)
702  double charcross = fmod(start_time,T2); // (start_time mod T2)
703  int s_ind_0 = -1;
704  // find s_ind_0 = first fast time point >= charcross.
705  // This could probably be made faster: FIXME
706  for (int i=0 ; i<=blockCount ; ++i)
707  {
708  if (fastTimes[i] >= charcross)
709  {
710  s_ind_0 = i;
711  break;
712  }
713  }
714  if (s_ind_0 == -1)
715  {
716  std::string msg = "N_TIA_OneStep::printMPDEOutputSolution: ";
717  msg += "Cannot find where characteristic curve crosses fast time slice at start_time";
718  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
719  return(false);
720  }
721  std::vector<double> h2(blockCount,0);
722  for (int j=0 ; j < blockCount ; ++j)
723  {
724  h2[j] = fastTimes[j+1] - fastTimes[j];
725  }
726  std::vector<double> ti;
727  //double first_interp = floor(start_time/T2)*T2 + fastTimes[s_ind_0];
728  double first_interp = start_time - charcross + fastTimes[s_ind_0];
729 #ifdef Xyce_DEBUG_TIME
730  if (tiaParams.debugLevel > 0)
731  {
732  Xyce::dout() << "first_interp = " << first_interp << std::endl;
733  }
734 #endif // Xyce_DEBUG_TIME
735  if (s_ind_0 == blockCount) { s_ind_0 = 0; };
736  // Don't interpolate the first point
737  double eps = fabs(start_time)*1.0e-6;
738  if ( fabs(first_interp-timept_) <= eps )
739  {
740  first_interp += h2[s_ind_0];
741  s_ind_0++;
742  if (s_ind_0 == blockCount) { s_ind_0 = 0; };
743 #ifdef Xyce_DEBUG_TIME
744  if (tiaParams.debugLevel > 0)
745  {
746  Xyce::dout() << "Moving first_interp forward to avoid duplicate outputs: " << first_interp << std::endl;
747  }
748 #endif // Xyce_DEBUG_TIME
749  }
750  int sn = s_ind_0;
751  double t = first_interp;
752  while (t <= stop_time)
753  {
754  ti.push_back(t);
755  t += h2[sn];
756  sn++;
757  if (sn >= blockCount) { sn = 0; }
758  }
759 #ifdef Xyce_DEBUG_TIME
760  if (tiaParams.debugLevel > 0)
761  {
762  Xyce::dout() << "T2 = " << T2 << std::endl;
763  Xyce::dout() << "charcross = " << charcross << std::endl;
764  Xyce::dout() << "s_ind_0 = " << s_ind_0 << std::endl;
765  Xyce::dout() << "Expecting to interpolate the following points:" << std::endl;
766  unsigned int numinterp = ti.size();
767  for (unsigned int i=0 ; i < numinterp ; ++i)
768  {
769  Xyce::dout() << ti[i] << std::endl;
770  }
771  Xyce::dout() << "Total of " << numinterp << " points" << std::endl;
772  }
773 #endif // Xyce_DEBUG_TIME
774  timept_ = start_time; // used later for interpolating stop_time
775  unsigned int tinum = ti.size();
776  int total_interp = 0;
777  std::vector<double> timepoint_vec(blockCount,stop_time);
778  int num_interp_this_cycle = 0;
779  int s_ind = s_ind_0;
780  for (unsigned int i=0; i < tinum ; ++i)
781  {
782  timepoint_vec[s_ind] = ti[i];
783  num_interp_this_cycle++;
784  s_ind++;
785  if (s_ind >= blockCount) { s_ind = 0; };
786  // If we're back to s_ind_0 or out of ti points, then interpolate:
787  if ((s_ind == s_ind_0) || (i == tinum-1))
788  {
790  // Now we print these points out
791  int s = s_ind_0;
792  for (int j=0 ; j < num_interp_this_cycle ; ++j)
793  {
794  timept_ = timepoint_vec[s];
795  outputManagerAdapter.tranOutput(
796  timept_, blockTmpSolVectorPtr->block(s),
800  total_interp++;
801  s++;
802  if (s >= blockCount) { s = 0; }
803 #ifdef Xyce_DEBUG_TIME
804  Xyce::dout() << "Interpolated to t = " << timept_ << std::endl;
805 #endif // Xyce_DEBUG_TIME
806  }
807  num_interp_this_cycle = 0;
808  }
809  }
810 #ifdef Xyce_DEBUG_TIME
811  Xyce::dout() << "Total of " << total_interp << " points" << std::endl;
812 #endif // Xyce_DEBUG_TIME
813 
814  // Now we interpolate stop_time unless its too close to the last timept interpolated.
815  eps = fabs(stop_time)*1.0e-8;
816  // fudge factor for printing, this should come from elsewhere FIXME
817  if (fabs(timept_ - stop_time) >= eps)
818  {
819 #ifdef Xyce_DEBUG_TIME
820  if (tiaParams.debugLevel > 0)
821  {
822  Xyce::dout() << "Previous timept = " << timept_ << std::endl;
823  Xyce::dout() << "Expecting to interpolate the following point: " << stop_time << std::endl;
824  }
825 #endif // Xyce_DEBUG_TIME
826  N_LAS_Vector* tmpSolnVecPtr = solnVecPtr;
827  N_LAS_Vector* tmpVecPtr = ds.tmpXn0APtr;
828  if (stop_time < tn)
829  {
831  tmpSolnVecPtr = ds.tmpXn0BPtr;
832  }
833  N_LAS_BlockVector * blockTmpSolnVecPtr =
834  dynamic_cast<N_LAS_BlockVector*>(tmpSolnVecPtr);
835  N_LAS_BlockVector * blockTmpVecPtr =
836  dynamic_cast<N_LAS_BlockVector*>(tmpVecPtr);
837  if (blockTmpSolnVecPtr == NULL)
838  {
839  std::string msg = "N_TIA_OneStep::printMPDEOutputSolution: ";
840  msg += "N_LAS_Vector tmpSolnVecPtr is not of type N_LAS_BlockVector";
841  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
842  return(false);
843  }
844  if (blockTmpVecPtr == NULL)
845  {
846  std::string msg = "N_TIA_OneStep::printMPDEOutputSolution: ";
847  msg += "N_LAS_Vector tmpVecPtr is not of type N_LAS_BlockVector";
848  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
849  return(false);
850  }
851  // Interpolate where the caracteristic crosses stop_time (instead of start_time).
852  //charcross = stop_time - floor(stop_time/T2)*T2;
853  charcross = fmod(stop_time,T2);
854  int s_ind_1 = -1;
855  // Find index just before charcross:
856  if( charcross < fastTimes[0] )
857  {
858  // by periodicity, the one before in this case is the one at the end
859  s_ind_1 = blockCount-1;
860  }
861  else
862  {
863  for (int i=blockCount-1 ; i>=0 ; --i)
864  {
865  if (fastTimes[i] <= charcross)
866  {
867  s_ind_1 = i;
868  break;
869  }
870  }
871  }
872  if (s_ind_1 == -1)
873  {
874  std::string msg = "N_TIA_OneStep::printMPDEOutputSolution: ";
875  msg += "Cannot find where characteristic curve crosses fast time slice at stop_time";
876  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
877  return(false);
878  }
879  int sm = s_ind_1;
880  int sp = s_ind_1+1;
881  double coeff_sm = fastTimes[sp]-charcross;
882  double coeff_sp = charcross-fastTimes[sm];
883  if (sp == blockCount) { sp = 0; }
884  double dt = h2[s_ind_1];
885  timept_ = stop_time;
886 #ifdef Xyce_DEBUG_TIME
887  Xyce::dout() << "charcross = " << charcross << std::endl;
888  Xyce::dout() << "s_ind_1 = " << s_ind_1 << std::endl;
889  Xyce::dout() << "sp = " << sp << std::endl;
890  Xyce::dout() << "sm = " << sm << std::endl;
891  Xyce::dout() << "dt = " << dt << std::endl;
892  Xyce::dout() << "timept = " << timept_ << std::endl;
893  Xyce::dout() << "coeff_sm = " << coeff_sm << std::endl;
894  Xyce::dout() << "coeff_sp = " << coeff_sp << std::endl;
895 #endif // Xyce_DEBUG_TIME
896  blockTmpVecPtr->block(0).linearCombo(
897  coeff_sm/dt, blockTmpSolnVecPtr->block(sm),
898  coeff_sp/dt, blockTmpSolnVecPtr->block(sp)
899  );
900  outputManagerAdapter.tranOutput(timept_, blockTmpVecPtr->block(0),
904 #ifdef Xyce_DEBUG_TIME
905  Xyce::dout() << "Interpolated to t = " << timept_ << std::endl;
906 #endif // Xyce_DEBUG_TIME
907  }
908 #ifdef Xyce_DEBUG_TIME
909  else // no extra interpolation
910  {
911  if (tiaParams.debugLevel > 0)
912  {
913  Xyce::dout() << "No further interpolation required." << std::endl;
914  }
915  }
916 #endif // Xyce_DEBUG_TIME
917 
918 #ifdef Xyce_DEBUG_TIME
919  Xyce::dout() << Xyce::section_divider << std::endl;
920 #endif // Xyce_DEBUG_TIME
921  return true;
922 }
923 
924 //-----------------------------------------------------------------------------
925 // Function : N_TIA_OneStep::printWaMPDEOutputSolution()
926 // Purpose : Print transient output from WaMPDE simulation
927 // Special Notes : This routine uses interpolateSolution.
928 // Scope : public
929 // Creator : Ting Mei, SNL
930 // Creation Date : 12/15/06
931 //-----------------------------------------------------------------------------
933  N_ANP_OutputMgrAdapter & outputManagerAdapter,
934  const double time,
935  N_LAS_Vector * solnVecPtr,
936  const std::vector<double> & fastTimes,
937  const int phiGID )
938 {
939 #ifdef Xyce_DEBUG_TIME
940  if (tiaParams.debugLevel > 0)
941  {
942  Xyce::dout() << std::endl;
943  Xyce::dout() << Xyce::section_divider << std::endl;
944  Xyce::dout() <<
945  " N_TIA_OneStep::printWaMPDEOutputSolution" << std::endl;
946  }
947 #endif // Xyce_DEBUG_TIME
948  double timestep = sec.lastAttemptedTimeStep;
949  double lasttime = sec.currentTime - timestep;
950  double tn = sec.currentTime;
951  // Set these values up to read output time intervals. FIXME
952  double beg_of_output_time_interval = lasttime;
953  double end_of_output_time_interval = tn;
954  double start_time = max(lasttime,beg_of_output_time_interval);
955  double stop_time = min(tn,end_of_output_time_interval);
956 #ifdef Xyce_DEBUG_TIME
957  if (tiaParams.debugLevel > 0)
958  {
959  Xyce::dout() << "start_time = " << start_time << std::endl;
960  Xyce::dout() << "stop_time = " << stop_time << std::endl;
961  }
962 #endif // Xyce_DEBUG_TIME
963 
964  // 12/15/06 tscoffe:
965  // First break the interval up as printOutputSolution does:
966  // hh = timestep/(sec.usedOrder_) and interpolate in the intervals:
967  // [tn+(i-1)*hh,tn+i*hh], i=1..usedOrder
968  // Assume phi(t_1) is linear in these intervals and approximate with:
969  // phi(t) = (1/hh)*(phi(tn)(tn+hh-t)+phi(tn+hh)(t-tn))
970  // Then the t_1 values we want to interpolate are:
971  // n2 = number of fast time points.
972  // T2 = fast time period
973  // h2 = T2/n2 = average spacing on fast time scale
974  // t1_vals = [tn:h2:tn+hh]
975  // And the t_2 values we want to interpolate are:
976  // t2_vals = phi(t1_vals) mod T2
977  // Then take the N_LAS blocks and do 2D linear interpolation on the intervals:
978  // (t1,s1), (t1,s2), (t2,s1), (t2,s2)
979  // x(t) = x(t,s) approx =
980  // (1/(t2-t1))(1/(s2-s1))[ x(t1,s1)(t2-t)(s2-s)
981  // +x(t1,s2)(t2-t)(s-s1)
982  // +x(t2,s1)(t-t1)(s2-s)
983  // +x(t2,s2)(t-t1)(s-s1) ]
984  // where t = t1_vals and s = t2_vals
985 
986  N_LAS_BlockVector * blockTmpSolVectorPtr =
987  dynamic_cast<N_LAS_BlockVector*>(ds.tmpSolVectorPtr);
988  N_LAS_BlockVector * blockTmpXn0APtr =
989  dynamic_cast<N_LAS_BlockVector*>(ds.tmpXn0APtr);
990  N_LAS_BlockVector * blockTmpXn0BPtr =
991  dynamic_cast<N_LAS_BlockVector*>(ds.tmpXn0BPtr);
992  if (blockTmpSolVectorPtr == NULL)
993  {
994  std::string msg = "N_TIA_OneStep::printWaMPDEOutputSolution: ";
995  msg += "N_LAS_Vector ds.tmpSolVectorPtr is not of type N_LAS_BlockVector";
996  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
997  return(false);
998  }
999  if (blockTmpXn0APtr == NULL)
1000  {
1001  std::string msg = "N_TIA_OneStep::printWaMPDEOutputSolution: ";
1002  msg += "N_LAS_Vector ds.tmpXn0APtr is not of type N_LAS_BlockVector";
1003  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
1004  return(false);
1005  }
1006  if (blockTmpXn0BPtr == NULL)
1007  {
1008  std::string msg = "N_TIA_OneStep::printWaMPDEOutputSolution: ";
1009  msg += "N_LAS_Vector ds.tmpXn0BPtr is not of type N_LAS_BlockVector";
1010  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
1011  return(false);
1012  }
1013  int phiLID = blockTmpSolVectorPtr->pmap()->globalToLocalIndex(phiGID);
1014  double hh = timestep/(sec.usedOrder_);
1015  double timeA = -1.0;
1016  double timeB = -1.0;
1017 #ifdef Xyce_DEBUG_TIME
1018  if (tiaParams.debugLevel > 0)
1019  {
1020  Xyce::dout() << " sec.usedOrder_ = " << sec.usedOrder_ << std::endl;
1021  Xyce::dout() << " sec.currentTime_ = " << sec.currentTime << std::endl;
1022  Xyce::dout() << " lasttime = " << lasttime << std::endl;
1023  }
1024 #endif // Xyce_DEBUG_TIME
1025 
1026  for (int i=0 ; i < sec.usedOrder_ ; ++i)
1027  {
1028  if (i == 0)
1029  {
1030  bool junk;
1031  timeA = lasttime + hh*i;
1032  junk = interpolateSolution(timeA,ds.tmpXn0APtr, ds.xHistory);
1033  if (!junk)
1034  {
1035  std::string msg = "interpolateSolution returned false!";
1036  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0,msg);
1037  }
1038  }
1039  else
1040  {
1041  // We don't need to interpolate this again.
1042  *ds.tmpXn0APtr = *ds.tmpXn0BPtr;
1043  timeA = timeB;
1044  }
1045  timeB = lasttime + hh*(i+1);
1047 #ifdef Xyce_DEBUG_TIME
1048  if (tiaParams.debugLevel > 0)
1049  {
1050  Xyce::dout() << "Interpolating in [ " << timeA << ", " << timeB << " ]" << std::endl;
1051  Xyce::dout() << "timeA = " << timeA << std::endl;
1052  Xyce::dout() << "timeB = " << timeB << std::endl;
1053  }
1054 #endif // Xyce_DEBUG_TIME
1055 
1056  // Now we can interpolate [tmpXn0APtr,tmpXn0BPtr] in [timeA,timeB].
1057  std::vector<double> t1vals;
1058  double T2 = fastTimes.back();
1059  int blockCount = blockTmpSolVectorPtr->blockCount();
1060  double h2 = T2/blockCount; // Average mesh width in fast time-scale
1061  double tval = timeA+h2;
1062  while (tval <= timeB)
1063  {
1064  t1vals.push_back(tval);
1065  tval += h2;
1066  }
1067  // fudge factor for printing, this should come from elsewhere FIXME
1068  double eps = fabs(timeB)*1.0e-8;
1069  if ( (t1vals.size() == 0) || (fabs(t1vals.back() - timeB) >= eps) )
1070  {
1071  t1vals.push_back(timeB);
1072  }
1073  std::vector<double> t2vals, phiAB(2);
1074  std::vector<double> tmpPhiAB(2, 0.0);
1075  if (phiLID >= 0)
1076  {
1077  tmpPhiAB[0] = (*ds.tmpXn0APtr)[phiLID]; // Get from MPDE Manager
1078  tmpPhiAB[1] = (*ds.tmpXn0BPtr)[phiLID]; // Get from MPDE Manager
1079  }
1080  blockTmpSolVectorPtr->pmap()->pdsComm()->sumAll( &tmpPhiAB[0], &phiAB[0], 2 );
1081 
1082  double phiA = phiAB[0], phiB = phiAB[1];
1083  for (unsigned int j=0 ; j<t1vals.size() ; ++j)
1084  {
1085  double phi = (1/(timeB-timeA))*(phiA*(timeB-t1vals[j])+phiB*(t1vals[j]-timeA));
1086  t2vals.push_back(fmod(phi,T2)); // phi(t1vals[j]) mod T2
1087  }
1088 #ifdef Xyce_DEBUG_TIME
1089  if (tiaParams.debugLevel > 0)
1090  {
1091  Xyce::dout() << "t1vals = " << std::endl;
1092  for (unsigned int j=0 ; j < t1vals.size() ; ++j)
1093  {
1094  Xyce::dout() << t1vals[j] << std::endl;
1095  }
1096  Xyce::dout() << "phi(" << timeA << ") = " << phiA << std::endl;
1097  Xyce::dout() << "phi(" << timeB << ") = " << phiB << std::endl;
1098  Xyce::dout() << "t2vals = " << std::endl;
1099  for (unsigned int j=0 ; j< t2vals.size() ; ++j)
1100  {
1101  Xyce::dout() << t2vals[j] << std::endl;
1102  }
1103  }
1104 #endif // Xyce_DEBUG_TIME
1105  // Now we can do our block 2D interpolations
1106  // loop through t1vals and move the fast time blocks as we go
1107  double t1 = timeA; // slow time indices
1108  double t2 = timeB;
1109  double t = t1vals[0]; // current time indices
1110  double s = t2vals[0];
1111  int b1,b2; // fast time block indices corresponding to s1,s2
1112  // Find the block surrounding s:
1113  b1 = -2;
1114  for (int j=0 ; j < blockCount ; ++j)
1115  {
1116  if ((fastTimes[j] <= s) && (s < fastTimes[j+1]))
1117  {
1118  b1 = j;
1119  }
1120  }
1121  b2 = b1+1;
1122  if (b2 == blockCount)
1123  {
1124  b2 = 0;
1125  }
1126  double s1 = fastTimes[b1];
1127  double s2 = fastTimes[b1+1]; // Note: fastTimes[blockCount] = T2
1128  if ((s < s1) || (s > s2))
1129  {
1130  std::string msg = "N_TIA_OneStep::printWaMPDEOutputSolution: ";
1131  msg += " Interpolator cannot find a fast time block containing the first point ";
1132  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
1133  }
1134 #ifdef Xyce_DEBUG_TIME
1135  Xyce::dout() << "Found s = " << s << " in block " << b1;
1136  Xyce::dout() << " with boundary = [" << s1 << "," << s2 << "]" << std::endl;
1137 #endif // Xyce_DEBUG_TIME
1138  for (unsigned int j=0 ; j < t1vals.size() ; ++j)
1139  {
1140  t = t1vals[j];
1141  s = t2vals[j];
1142  if (t > t2) break; // This should never occur
1143  // If s has moved outside our block, then increment block.
1144  if ( (s < s1) || (s > s2) )
1145  {
1146 #ifdef Xyce_DEBUG_TIME
1147  Xyce::dout() << "Incrementing fast time block for next interpolation." << std::endl;
1148 #endif // Xyce_DEBUG_TIME
1149  b1++;
1150  if (b1 == blockCount)
1151  {
1152  b1 = 0;
1153  }
1154  b2 = b1+1;
1155  if (b2 == blockCount)
1156  {
1157  b2 = 0;
1158  }
1159  s1 = fastTimes[b1];
1160  s2 = fastTimes[b1+1];
1161  }
1162  // If s isn't in the next block, then search for it.
1163  if ((s < s1) || (s > s2))
1164  {
1165 #ifdef Xyce_DEBUG_TIME
1166  Xyce::dout() << "Searching for fast time block for next interpolation." << std::endl;
1167 #endif // Xyce_DEBUG_TIME
1168  b1 = -2;
1169  for (int j2=0 ; j2 < blockCount ; ++j2)
1170  {
1171  if ((fastTimes[j2] <= s) && (s < fastTimes[j2+1]))
1172  {
1173  b1 = j2;
1174  }
1175  }
1176  b2 = b1+1;
1177  if (b2 == blockCount)
1178  {
1179  b2 = 0;
1180  }
1181  s1 = fastTimes[b1];
1182  s2 = fastTimes[b1+1];
1183  }
1184  // If a block surrounding s can't be found, then quit.
1185  if ((s < s1) || (s > s2))
1186  {
1187  std::string msg = "N_TIA_OneStep::printWaMPDEOutputSolution: ";
1188  msg += " Interpolator moved fast time block but new point is not in this block ";
1189  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
1190  }
1191  if (s > T2) break; // Just double checking...
1192  //blockTmpXn0APtr->block(b1) // x(t1,s1)
1193  //blockTmpXn0APtr->block(b2) // x(t1,s2)
1194  //blockTmpXn0BPtr->block(b1) // x(t2,s1)
1195  //blockTmpXn0BPtr->block(b2) // x(t2,s2)
1196  // (1/(t2-t1))(1/(s2-s1))[ x(t1,s1)(t2-t)(s2-s)
1197  // +x(t1,s2)(t2-t)(s-s1)
1198  // +x(t2,s1)(t-t1)(s2-s)
1199  // +x(t2,s2)(t-t1)(s-s1) ]
1200 #ifdef Xyce_DEBUG_TIME
1201  if (tiaParams.debugLevel > 0)
1202  {
1203  Xyce::dout() << "Interpolating in the block:" << std::endl;
1204  Xyce::dout() << "(t1,t2) = (" << t1 << "," << t2 << ")" << std::endl;
1205  Xyce::dout() << "(s1,s2) = (" << s1 << "," << s2 << ")" << std::endl;
1206  }
1207 #endif // Xyce_DEBUG_TIME
1208  double denom = (t2-t1)*(s2-s1);
1209  double coeff0 = (t2-t)*(s2-s)/denom;
1210  double coeff1 = (t2-t)*(s-s1)/denom;
1211  double coeff2 = (t-t1)*(s2-s)/denom;
1212  double coeff3 = (t-t1)*(s-s1)/denom;
1213  (blockTmpSolVectorPtr->block(b1)).linearCombo(
1214  coeff0, blockTmpXn0APtr->block(b1),
1215  coeff1, blockTmpXn0APtr->block(b2),
1216  coeff2, blockTmpXn0BPtr->block(b1) );
1217  (blockTmpSolVectorPtr->block(b1)).update(
1218  coeff3, blockTmpXn0BPtr->block(b2), 1.0 );
1219 
1220  // erkeite 2/24/07. This is needed because currently the interpolation goes back to t=-1.0.
1221  if (t >= 0.0)
1222  {
1223  outputManagerAdapter.tranOutput(t, blockTmpSolVectorPtr->block(b1),
1227 #ifdef Xyce_DEBUG_TIME
1228  Xyce::dout() << "Interpolated to (t,phi(t)) = (" << t << "," << s << ")" << std::endl;
1229 #endif // Xyce_DEBUG_TIME
1230  }
1231  }
1232  }
1233 
1234 #ifdef Xyce_DEBUG_TIME
1235  Xyce::dout() << Xyce::section_divider << std::endl;
1236 #endif // Xyce_DEBUG_TIME
1237  return true;
1238 }
1239 
1240 //-----------------------------------------------------------------------------
1241 // Function : N_TIA_OneStep::printOutputSolution()
1242 // Purpose : Print output that is dumbed down in order.
1243 // Special Notes : This routine picks smaller time steps to approximate first
1244 // : order integration from the perspective of the output.
1245 // Scope : public
1246 // Creator : Ting Mei, SNL
1247 // Creation Date : 11/16/07
1248 //-----------------------------------------------------------------------------
1250  N_ANP_OutputMgrAdapter & outputManagerAdapter,
1251  const double time,
1252  N_LAS_Vector * solnVecPtr,
1253  const bool doNotInterpolate,
1254  const std::vector<double> &outputInterpolationTimes,
1255  bool skipPrintLineOutput)
1256 {
1257 #ifdef Xyce_DEBUG_TIME
1258  if (tiaParams.debugLevel > 0)
1259  {
1260  Xyce::dout() << std::endl;
1261  Xyce::dout() << Xyce::section_divider << std::endl;
1262  Xyce::dout() <<
1263  " N_TIA_OneStep::printOutputSolution" << std::endl;
1264  Xyce::dout() << "usedOrder_ = " << sec.usedOrder_ << std::endl;
1265  }
1266 #endif // Xyce_DEBUG_TIME
1267  double timestep = sec.lastAttemptedTimeStep;
1268  double lasttime = sec.currentTime - timestep;
1269  bool dointerp = true;
1270  double hh = timestep/(sec.usedOrder_);
1271 
1272  if (hh <= 10*sec.minTimeStep)
1273  {
1274  dointerp = false;
1275  }
1276 
1277  if (!(tiaParams.interpOutputFlag))
1278  {
1279  dointerp = false;
1280  }
1281 
1282  if (doNotInterpolate)
1283  {
1284  dointerp = false;
1285  }
1286 
1287  if (dointerp && !outputInterpolationTimes.empty())
1288  {
1289  for (unsigned int i=0;i<outputInterpolationTimes.size();++i)
1290  {
1291  interpolateSolution(outputInterpolationTimes[i], ds.tmpSolVectorPtr, ds.xHistory); // interpolate solution
1292  interpolateSolution(outputInterpolationTimes[i], ds.tmpStaVectorPtr, ds.sHistory); // interpolate state
1293  interpolateSolution(outputInterpolationTimes[i], ds.tmpStoVectorPtr, ds.stoHistory);// interpolate store
1294 
1295  outputManagerAdapter.tranOutput(outputInterpolationTimes[i], *ds.tmpSolVectorPtr,
1299  }
1300  }
1301 
1302  // Either way, do an output on the actual computed time step, but only
1303  // if we weren't given a list of specific times *or* we were told not to
1304  // interpoloate.
1305  if (outputInterpolationTimes.empty() || doNotInterpolate)
1306  {
1307  outputManagerAdapter.tranOutput(time, *ds.currSolutionPtr,
1311  skipPrintLineOutput);
1312  }
1313 
1314 #ifdef Xyce_DEBUG_TIME
1315  Xyce::dout() << Xyce::section_divider << std::endl;
1316 #endif // Xyce_DEBUG_TIME
1317  return true;
1318 }
1319 
1320 //-----------------------------------------------------------------------------
1321 // Function : N_TIA_OneStep::saveOutputSolution
1322 // Purpose : This is similar to printOutputSolution, but is in support of
1323 // the .SAVE capability, rather than .PRINT.
1324 // Special Notes :
1325 // Scope : public
1326 // Creator : Eric Keiter, SNL
1327 // Creation Date : 10/21/07
1328 //-----------------------------------------------------------------------------
1330  N_ANP_OutputMgrAdapter & outputManagerAdapter,
1331  N_LAS_Vector * solnVecPtr,
1332  const double saveTime,
1333  const bool doNotInterpolate)
1334 {
1335 #ifdef Xyce_DEBUG_TIME
1336  if (tiaParams.debugLevel > 0)
1337  {
1338  Xyce::dout() << std::endl;
1339  Xyce::dout() << Xyce::section_divider << std::endl;
1340  Xyce::dout() <<
1341  " N_TIA_OneStep::saveOutputSolution" << std::endl;
1342  }
1343 #endif // Xyce_DEBUG_TIME
1344 
1345  double timestep = sec.lastAttemptedTimeStep;
1346  double lasttime = sec.currentTime - timestep;
1347  bool dointerp = true;
1348  double hh = timestep/(sec.usedOrder_);
1349 
1350  outputManagerAdapter.outputDCOP( *(ds.currSolutionPtr) );
1351 
1352 #ifdef Xyce_DEBUG_TIME
1353  Xyce::dout() << Xyce::section_divider << std::endl;
1354 #endif // Xyce_DEBUG_TIME
1355  return true;
1356 }
1357 
1358 //-----------------------------------------------------------------------------
1359 // Function : N_TIA_OneStep::updateHistory
1360 // Purpose : Update history array after a successful step
1361 // Special Notes :
1362 // Scope : public
1363 // Creator : Ting Mei, SNL
1364 // Creation Date : 11/16/07
1365 //-----------------------------------------------------------------------------
1367 {
1368 #ifdef Xyce_DEBUG_TIME
1369  if (tiaParams.debugLevel > 1)
1370  {
1371  Xyce::dout() << std::endl;
1372  Xyce::dout() << Xyce::section_divider << std::endl;
1373  Xyce::dout() <<
1374  " N_TIA_OneStep::updateHistory" << std::endl;
1375  Xyce::dout() << "\n Before updates \n" << std::endl;
1376  for (int i=0; i<=sec.maxOrder_ ; ++i)
1377  {
1378  Xyce::dout() << "\n xHistory["<< i << "]: \n" << std::endl;
1379  (ds.xHistory[i])->printPetraObject(Xyce::dout());
1380  Xyce::dout() << std::endl;
1381  }
1382  for (int i=0; i<=sec.maxOrder_ ; ++i)
1383  {
1384  Xyce::dout() << "\n qHistory["<< i << "]: \n" << std::endl;
1385  (ds.qHistory[i])->printPetraObject(Xyce::dout());
1386  Xyce::dout() << std::endl;
1387  }
1388  for (int i=0; i<=sec.maxOrder_ ; ++i)
1389  {
1390  Xyce::dout() << "\n sHistory["<< i << "]: \n" << std::endl;
1391  (ds.sHistory[i])->printPetraObject(Xyce::dout());
1392  Xyce::dout() << std::endl;
1393  }
1394  Xyce::dout() << Xyce::section_divider << std::endl;
1395  }
1396 #endif // Xyce_DEBUG_TIME
1397 
1398 
1399  if (sec.currentOrder_ == 2)
1400  {
1401  *(ds.xHistory[2]) = *(ds.xHistory[1]);
1402 #ifdef SEPARATE_F_AND_B
1403  (ds.qHistory[2])->linearCombo(1.0, *(ds.daeFVectorPtr), -1.0, *(ds.daeBVectorPtr));
1404 #else
1405  *(ds.qHistory[2]) = *(ds.daeFVectorPtr);
1406 #endif
1407  *(ds.stoHistory[2]) = *(ds.stoHistory[1]);
1409  }
1410 
1411  (ds.xHistory[1])->linearCombo(1.0, *ds.nextSolutionPtr, -1.0,*(ds.xHistory[0]));
1412  (ds.qHistory[1])->linearCombo(1.0, *ds.daeQVectorPtr, -1.0,*(ds.qHistory[0]));
1413 
1414  (ds.stoHistory[1])->linearCombo(1.0, *ds.nextStorePtr, -1.0,*(ds.stoHistory[0]));
1415  (ds.stoLeadCurrQHistory[1])->linearCombo(1.0, *ds.nextStoreLeadCurrQPtr, -1.0,*(ds.stoLeadCurrQHistory[0]));
1416 
1417  *(ds.xHistory[0]) = *ds.nextSolutionPtr;
1418  *(ds.qHistory[0]) = *ds.daeQVectorPtr;
1419  *(ds.sHistory[0]) = *ds.nextStatePtr;
1420  *(ds.stoHistory[0]) = *ds.nextStorePtr;
1422 
1423 #ifdef Xyce_DEBUG_TIME
1424  if (tiaParams.debugLevel > 1)
1425  {
1426  Xyce::dout() << "\n After updates \n" << std::endl;
1427  Xyce::dout() << "\n newtonCorrectionPtr: " << std::endl;
1428  ds.newtonCorrectionPtr->printPetraObject(Xyce::dout());
1429  Xyce::dout() << "\n qnewtonCorrectionPtr: " << std::endl;
1430  ds.qNewtonCorrectionPtr->printPetraObject(Xyce::dout());
1431  for (int i=0; i<=sec.maxOrder_ ; ++i)
1432  {
1433  Xyce::dout() << "\n xHistory["<< i << "]: \n" << std::endl;
1434  (ds.xHistory[i])->printPetraObject(Xyce::dout());
1435  Xyce::dout() << std::endl;
1436  }
1437  for (int i=0; i<=sec.maxOrder_ ; ++i)
1438  {
1439  Xyce::dout() << "\n qHistory["<< i << "]: \n" << std::endl;
1440  (ds.qHistory[i])->printPetraObject(Xyce::dout());
1441  Xyce::dout() << std::endl;
1442  }
1443  Xyce::dout() << "\n sNewtonCorrectionPtr: " << std::endl;
1444  ds.sNewtonCorrectionPtr->printPetraObject(Xyce::dout());
1445  Xyce::dout() << std::endl;
1446  Xyce::dout() << "\n nextStatePtr: " << std::endl;
1447  ds.nextStatePtr->printPetraObject(Xyce::dout());
1448  Xyce::dout() << std::endl;
1449  for (int i=0; i<=sec.maxOrder_ ; ++i)
1450  {
1451  Xyce::dout() << "\n sHistory["<< i << "]: \n" << std::endl;
1452  (ds.sHistory[i])->printPetraObject(Xyce::dout());
1453  Xyce::dout() << std::endl;
1454  }
1455  Xyce::dout() << Xyce::section_divider << std::endl;
1456  }
1457 #endif // Xyce_DEBUG_TIME
1458 
1460 }
1461 
1462 //-----------------------------------------------------------------------------
1463 // Function : N_TIA_OneStep::updateSensitivityHistory
1464 // Purpose : Update sensitivity history array after a successful step
1465 // Special Notes :
1466 // Scope : public
1467 // Creator : Eric Keiter, SNL
1468 // Creation Date :
1469 //-----------------------------------------------------------------------------
1471 {
1472  int numParams = ds.sensRHSPtrVector.size();
1473  for (int ip=0; ip<numParams;++ip)
1474  {
1475  std::vector<N_LAS_Vector*> & dqdpHistory = ds.dqdpHistory[ip];
1476  std::vector<N_LAS_Vector*> & dXdpHistory = ds.dXdpHistory[ip];
1477 
1478  N_LAS_Vector & nextDXdpVec = *(ds.nextDXdpPtrVector[ip]);
1479  N_LAS_Vector & nextDfdpVec = *(ds.nextDfdpPtrVector[ip]);
1480  N_LAS_Vector & nextDqdpVec = *(ds.nextDqdpPtrVector[ip]);
1481  N_LAS_Vector & nextDbdpVec = *(ds.nextDbdpPtrVector[ip]);
1482 
1483  if (sec.currentOrder_ == 2)
1484  {
1485  *(dXdpHistory[2]) = *(dXdpHistory[1]);
1486 #ifdef SEPARATE_F_AND_B
1487  (dqdpHistory[2])->linearCombo(1.0, nextDfdpVec, -1.0, nextDbdpVec);
1488 #else
1489  *(dqdpHistory[2]) = nextDfdpVec;
1490 #endif
1491  }
1492 
1493  (dXdpHistory[1])->linearCombo(1.0, nextDXdpVec, -1.0,*(dXdpHistory[0]));
1494  (dqdpHistory[1])->linearCombo(1.0, nextDqdpVec, -1.0,*(dqdpHistory[0]));
1495 
1496  *(dXdpHistory[0]) = nextDXdpVec;
1497  *(dqdpHistory[0]) = nextDqdpVec;
1498 
1499 #ifdef DEBUG_SENS
1500  Xyce::dout() << "dxdp history ip = " << ip << ":\n";
1501  for (int i=0; i<=sec.maxOrder_ ; ++i)
1502  {
1503  dXdpHistory[i]->printPetraObject(Xyce::dout());
1504  }
1505 #endif
1506  }
1507 }
1508 
1509 //-----------------------------------------------------------------------------
1510 // Function : N_TIA_OneStep::restoreHistory
1511 // Purpose : Restore history array after a failed step
1512 // Special Notes :
1513 // Scope : public
1514 // Creator : Ting Mei, SNL
1515 // Creation Date : 10/21/07
1516 //-----------------------------------------------------------------------------
1518 {
1519  for (int i=1;i<=sec.currentOrder_;++i)
1520  {
1521  sec.psi_[i-1] = sec.psi_[i];
1522  }
1523 #ifdef Xyce_DEBUG_TIME
1524  if (tiaParams.debugLevel > 1)
1525  {
1526  Xyce::dout() << std::endl;
1527  Xyce::dout() << Xyce::section_divider << std::endl;
1528  Xyce::dout() <<
1529  " N_TIA_OneStep::restoreHistory" << std::endl;
1530  for (int i=1;i<=sec.currentOrder_;++i)
1531  Xyce::dout() << "\n sec.psi_[i] = " << sec.psi_[i] << std::endl;
1532  for (int i=0; i<=sec.maxOrder_ ; ++i)
1533  {
1534  Xyce::dout() << "\n xHistory["<< i << "]: \n" << std::endl;
1535  (ds.xHistory[i])->printPetraObject(Xyce::dout());
1536  Xyce::dout() << std::endl;
1537  }
1538  for (int i=0; i<=sec.maxOrder_ ; ++i)
1539  {
1540  Xyce::dout() << "\n qHistory["<< i << "]: \n" << std::endl;
1541  (ds.qHistory[i])->printPetraObject(Xyce::dout());
1542  Xyce::dout() << std::endl;
1543  }
1544  for (int i=0; i<=sec.maxOrder_ ; ++i)
1545  {
1546  Xyce::dout() << "\n sHistory["<< i << "]: \n" << std::endl;
1547  (ds.sHistory[i])->printPetraObject(Xyce::dout());
1548  Xyce::dout() << std::endl;
1549  }
1550  Xyce::dout() << Xyce::section_divider << std::endl;
1551  }
1552 #endif // Xyce_DEBUG_TIME
1553 }
1554 
1555 //-----------------------------------------------------------------------------
1556 // Function : N_TIA_OneStep::updateCoeffs
1557 // Purpose : Update method coefficients
1558 // Special Notes :
1559 // Scope : public
1560 // Creator : Ting Mei, SNL
1561 // Creation Date : 10/21/07
1562 //-----------------------------------------------------------------------------
1564 {
1565  // synchronize with Step Error Control
1566  // sec.psi_[0] = sec.currentTimeStep;
1567 #ifdef Xyce_DEBUG_TIME
1568  if (tiaParams.debugLevel > 0)
1569  {
1570  Xyce::dout() << std::endl;
1571  Xyce::dout() << Xyce::section_divider << std::endl;
1572  Xyce::dout() << " N_TIA_OneStep::updateCoeffs" << std::endl;
1573  Xyce::dout() << " currentTimeStep = " << sec.currentTimeStep << std::endl;
1574  Xyce::dout() << " numberOfSteps_ = " << sec.numberOfSteps_ << std::endl;
1575  Xyce::dout() << " currentOrder_ = " << sec.currentOrder_ << std::endl;
1576  Xyce::dout() << " nscsco_ = " << sec.nscsco_ << std::endl;
1577  Xyce::dout() << " psi_[0] = " << sec.psi_[0] << std::endl;
1578  }
1579 #endif
1580 
1581  double temp1 = sec.currentTimeStep;
1582 
1583  if (sec.currentOrder_ == 2)
1584  {
1585  sec.psi_[2] = sec.psi_[1];
1586  }
1587  sec.psi_[1] = sec.psi_[0];
1588  sec.psi_[0] = temp1;
1589 
1590  sec.beta_[0] = 1.0;
1591  sec.alpha_[0] = 1.0;
1592  // sec.ck_ = 1.0;
1593  sec.alphas_ = -1.0;
1594 
1595  if (sec.currentOrder_ == 2)
1596  {
1597  double temp2 = sec.psi_[1];
1598  sec.beta_[1] = temp1/temp2 + (temp1/temp2)*(temp1/temp2)/2;
1599  sec.beta_[2] = -1.0*temp1*temp1/temp2/sec.psi_[2]/2;
1601  }
1602  else
1603  {
1604  sec.beta_[1] = temp1/sec.psi_[1];
1606  }
1607 
1608 #ifdef Xyce_DEBUG_TIME
1609  if (tiaParams.debugLevel > 0)
1610  {
1611  Xyce::dout() << " nscsco_ = " << sec.nscsco_ << std::endl;
1612  Xyce::dout() << " beta_[0] = " << sec.beta_[0] << std::endl;
1613  Xyce::dout() << " beta_[1] = " << sec.beta_[1] << std::endl;
1614  Xyce::dout() << " beta_[2] = " << sec.beta_[2] << std::endl;
1615  Xyce::dout() << " beta_[3] = " << sec.beta_[3] << std::endl;
1616  Xyce::dout() << " beta_[4] = " << sec.beta_[4] << std::endl;
1617  Xyce::dout() << " alpha_[0] = " << sec.alpha_[0] << std::endl;
1618  Xyce::dout() << " alpha_[1] = " << sec.alpha_[1] << std::endl;
1619  Xyce::dout() << " alpha_[2] = " << sec.alpha_[2] << std::endl;
1620  Xyce::dout() << " alpha_[3] = " << sec.alpha_[3] << std::endl;
1621  Xyce::dout() << " alpha_[4] = " << sec.alpha_[4] << std::endl;
1622  Xyce::dout() << " alphas_ = " << sec.alphas_ << std::endl;
1623  Xyce::dout() << " alpha0_ = " << sec.alpha0_ << std::endl;
1624  Xyce::dout() << " gamma_[0] = " << sec.gamma_[0] << std::endl;
1625  Xyce::dout() << " gamma_[1] = " << sec.gamma_[1] << std::endl;
1626  Xyce::dout() << " gamma_[2] = " << sec.gamma_[2] << std::endl;
1627  Xyce::dout() << " gamma_[3] = " << sec.gamma_[3] << std::endl;
1628  Xyce::dout() << " gamma_[4] = " << sec.gamma_[4] << std::endl;
1629  Xyce::dout() << " psi_[0] = " << sec.psi_[0] << std::endl;
1630  Xyce::dout() << " psi_[1] = " << sec.psi_[1] << std::endl;
1631  Xyce::dout() << " psi_[2] = " << sec.psi_[2] << std::endl;
1632  Xyce::dout() << " psi_[3] = " << sec.psi_[3] << std::endl;
1633  Xyce::dout() << " psi_[4] = " << sec.psi_[4] << std::endl;
1634  Xyce::dout() << " sigma_[0] = " << sec.sigma_[0] << std::endl;
1635  Xyce::dout() << " sigma_[1] = " << sec.sigma_[1] << std::endl;
1636  Xyce::dout() << " sigma_[2] = " << sec.sigma_[2] << std::endl;
1637  Xyce::dout() << " sigma_[3] = " << sec.sigma_[3] << std::endl;
1638  Xyce::dout() << " sigma_[4] = " << sec.sigma_[4] << std::endl;
1639  Xyce::dout() << " ck_ = " << sec.ck_ << std::endl;
1640  Xyce::dout() << Xyce::section_divider << std::endl;
1641  }
1642 #endif // Xyce_DEBUG_TIME
1643 }
1644 
1645 //-----------------------------------------------------------------------------
1646 // Function : N_TIA_OneStep::initialize
1647 // Purpose : Initialize method with initial solution & step-size
1648 // Special Notes :
1649 // Scope : public
1650 // Creator : Ting Mei, SNL
1651 // Creation Date : 10/21/07
1652 //-----------------------------------------------------------------------------
1654 {
1655 
1656  // we assume the solution vector is available here
1657  // Note that I'm using currSolutionPtr instead of
1658  // nextSolutionPtr because this is the first step.
1659 
1660  // Update next stop time from StepErrorControl:
1661  // ERK. Commenting this out, as it is already called from N_ANP_AnalysisManager,
1662  // right before this initialize call. It should not be called 2x, as
1663  // it is history dependent (unfortunately), so calling it 2x in a row changes
1664  // the stop time to a different number.
1665  // sec.updateStopTime();
1666 
1667  // Choose initial step-size
1668  double time_to_stop = sec.stopTime - sec.currentTime;
1669  double currentTimeStep;
1670 
1671  tiaParams.TimeStepLimitedbyBP = false;
1672 
1674  {
1675  currentTimeStep = 0.1 * time_to_stop;
1676  currentTimeStep = Xycemin(sec.startingTimeStep, currentTimeStep);
1677  sec.currentTimeStep = currentTimeStep;
1678  }
1679  else
1680  {
1681  // compute an initial step-size based on rate of change in the
1682  // solution initially
1683 #ifdef Xyce_INCOMPLETE_2LEVEL_NORMS
1684  double dnorm_q = 0.0;
1685  (ds.qHistory[1])->wRMSNorm(*ds.qErrWtVecPtr, &dnorm_q);
1686 #else
1687  double dnorm_q = ds.delta_x_errorNorm_q1();
1688 #endif
1689  if (dnorm_q > 0.0) // time-dependent DAE
1690  {
1691  if (sec.currentTime == sec.initialTime)
1692  currentTimeStep = Xycemin(sec.h0_max_factor_*abs(time_to_stop),sqrt(2.0)/(sec.h0_safety_*dnorm_q));
1693  else
1694  currentTimeStep = 0.1* Xycemin(sec.savedTimeStep, abs(time_to_stop));
1695  }
1696  else // non-time-dependent DAE
1697  {
1698  if (sec.currentTime == sec.initialTime)
1699  currentTimeStep = sec.h0_max_factor_*abs(time_to_stop);
1700  else
1701  currentTimeStep = 0.1* Xycemin(sec.savedTimeStep, abs(time_to_stop));
1702  }
1703  // choose min of user specified value and our value:
1704  if (sec.startingTimeStep > 0.0 && (sec.currentTime == sec.initialTime))
1705  currentTimeStep = Xycemin(sec.startingTimeStep, currentTimeStep);
1706  // check for maximum step-size:
1707  double rh = abs(currentTimeStep)*sec.h_max_inv_;
1708  if (rh>1.0) currentTimeStep = currentTimeStep/rh;
1709 
1710 
1711  // Apply this new stepsize only if it is smaller than the one preceding
1712  // the breakpoint, but only do this if this is a non-DCOP breakpoint.
1713  if (sec.currentTime != sec.initialTime) // if not DCOP:
1714  {
1715  sec.currentTimeStep = Xycemin(sec.currentTimeStep, currentTimeStep);
1716  }
1717  else // else if DCOP:
1718  {
1719  sec.currentTimeStep = currentTimeStep;
1720  }
1721  }
1722 
1723  sec.currentTimeStepRatio = 1.0;
1725 
1729 
1731  sec.stepAttemptStatus = true;
1732 
1733 #ifdef Xyce_VERBOSE_TIME
1734  if (tiaParams.errorAnalysisOption == 1)
1735  {
1736  Xyce::dout() << "ERROROPTION=1: DeltaT Grow = 2" << "\n" << std::endl;
1737  Xyce::dout() << "ERROROPTION=1: DeltaT Cut = 0.125" << "\n" << std::endl;
1738  Xyce::dout() << "ERROROPTION=1: NL MIN = " << tiaParams.NLmin << "\n" << std::endl;
1739  Xyce::dout() << "ERROROPTION=1: NL MAX = " << tiaParams.NLmax << "\n" << std::endl;
1740  Xyce::dout() << "ERROROPTION=1: DELMAX = " << sec.maxTimeStep << "\n" << std::endl;
1741  }
1742 #endif //Xyce_VERBOSE_TIME
1743 
1744  // sec.tolAimFac_ = 0.5;
1745 
1747 
1748  // x history
1749  *(ds.xHistory[0]) = *(ds.currSolutionPtr);
1750  (ds.xHistory[1])->putScalar(0.0); // no need to multiply by dt here
1751 
1752  // q history
1753  *(ds.qHistory[0]) = *(ds.daeQVectorPtr);
1754 #ifdef SEPARATE_F_AND_B
1755  (ds.qHistory[1])->linearCombo(1.0, *(ds.daeFVectorPtr), -1.0, *(ds.daeBVectorPtr));
1756 #else
1757  *(ds.qHistory[1]) = *(ds.daeFVectorPtr);
1758 #endif
1759  (ds.qHistory[1])->scale(-sec.currentTimeStep);
1760 
1761  // state history
1762  *(ds.sHistory[0]) = *(ds.currStatePtr);
1763  (ds.sHistory[1])->putScalar(0.0);
1764 
1765  // store history
1766  *(ds.stoHistory[0]) = *(ds.currStorePtr);
1767  (ds.stoHistory[1])->putScalar(0.0);
1768 
1769  // lead current Q compontent history
1771  (ds.stoLeadCurrQHistory[1])->putScalar(0.0);
1772 
1773  // Coefficient initialization
1774  sec.numberOfSteps_ = 0; // number of total time integration steps taken
1775  sec.currentOrder_ = 1;
1776  sec.usedOrder_ = 1;
1777  sec.psi_[0] = sec.currentTimeStep;
1778  sec.cj_ = 1/sec.psi_[0];
1779  sec.nscsco_ = 0;
1780 #ifdef Xyce_DEBUG_TIME
1781  if (tiaParams.debugLevel > 1)
1782  {
1783  Xyce::dout() << std::endl;
1784  Xyce::dout() << Xyce::section_divider << std::endl;
1785  Xyce::dout() <<
1786  " N_TIA_OneStep::initialize" << std::endl;
1787  Xyce::dout() << "\n xHistory: \n" << std::endl;
1788  (ds.xHistory[0])->printPetraObject(Xyce::dout());
1789  Xyce::dout() << std::endl;
1790  (ds.xHistory[1])->printPetraObject(Xyce::dout());
1791  Xyce::dout() << std::endl;
1792  Xyce::dout() << "\n qHistory: \n" << std::endl;
1793  (ds.qHistory[0])->printPetraObject(Xyce::dout());
1794  Xyce::dout() << std::endl;
1795  (ds.qHistory[1])->printPetraObject(Xyce::dout());
1796  Xyce::dout() << std::endl;
1797  Xyce::dout() << "\n sHistory: \n" << std::endl;
1798  (ds.sHistory[0])->printPetraObject(Xyce::dout());
1799  Xyce::dout() << std::endl;
1800  (ds.sHistory[1])->printPetraObject(Xyce::dout());
1801  Xyce::dout() << std::endl;
1802  Xyce::dout() << "\n" << "currentTimeStep = " << currentTimeStep << "\n" << std::endl;
1803  Xyce::dout() << "\n" << "time_to_stop = " << time_to_stop << "\n" << std::endl;
1804  Xyce::dout() << Xyce::section_divider << std::endl;
1805  }
1806 #endif // Xyce_DEBUG_TIME
1807 
1809 }
1810 
1811 //-----------------------------------------------------------------------------
1812 // Function : N_TIA_OneStep::initializeSensitivities
1813 // Purpose :
1814 // Special Notes :
1815 // Scope : public
1816 // Creator : Eric Keiter, SNL
1817 // Creation Date :
1818 //-----------------------------------------------------------------------------
1820 {
1821  int numParams = ds.sensRHSPtrVector.size();
1822  for (int ip=0; ip<numParams;++ip)
1823  {
1824  std::vector<N_LAS_Vector*> & dqdpHistory = ds.dqdpHistory[ip];
1825  std::vector<N_LAS_Vector*> & dXdpHistory = ds.dXdpHistory[ip];
1826 
1827  N_LAS_Vector * currDxdpPtr = ds.currDXdpPtrVector[ip];
1828  N_LAS_Vector * currDqdpPtr = ds.currDqdpPtrVector[ip];
1829 
1830  // next or curr?
1831  N_LAS_Vector & dfdpVec = *(ds.currDfdpPtrVector[ip]);
1832  N_LAS_Vector & dqdpVec = *(ds.currDqdpPtrVector[ip]);
1833  N_LAS_Vector & dbdpVec = *(ds.currDbdpPtrVector[ip]);
1834 
1835  // dXdp history
1836  *(dXdpHistory[0]) = *(currDxdpPtr);
1837  (dXdpHistory[1])->putScalar(0.0);
1838 
1839  // dqdp history
1840  *(dqdpHistory[0]) = dqdpVec;
1841 #ifdef SEPARATE_F_AND_B
1842  (dqdpHistory[1])->linearCombo(1.0, *(dfdpVec), -1.0, *(dbdpVec));
1843 #else
1844  *(dqdpHistory[1]) = dfdpVec;
1845 #endif
1846  (dqdpHistory[1])->scale(-sec.currentTimeStep);
1847  }
1848 }
1849 
1850 //-----------------------------------------------------------------------------
1851 // Function : N_TIA_OneStep::setTwoLevelTimeInfo
1852 // Purpose :
1853 // Special Notes :
1854 // Scope : public
1855 // Creator : Eric Keiter, SNL
1856 // Creation Date : 03/01/07
1857 //-----------------------------------------------------------------------------
1859 (const N_TIA_TimeIntInfo & tiInfo)
1860 {
1861  // set initial step-size
1862  double time_to_stop = sec.stopTime - sec.currentTime;
1863 
1864 
1865  // x history
1866  *(ds.xHistory[0]) = *(ds.currSolutionPtr);
1867  (ds.xHistory[1])->putScalar(0.0); // no need to multiply by dt here
1868 
1869  // q history
1870  *(ds.qHistory[0]) = *(ds.daeQVectorPtr);
1871 #ifdef SEPARATE_F_AND_B
1872  (ds.qHistory[1])->linearCombo(1.0, *(ds.daeFVectorPtr), -1.0, *(ds.daeBVectorPtr));
1873 #else
1874  *(ds.qHistory[1]) = *(ds.daeFVectorPtr);
1875 #endif
1876  (ds.qHistory[1])->scale(-sec.currentTimeStep);
1877 
1878  // state history
1879  *(ds.sHistory[0]) = *(ds.nextStatePtr);
1880  (ds.sHistory[1])->putScalar(0.0);
1881 
1882  // Coefficient initialization
1883  sec.numberOfSteps_ = 0; // number of total time integration steps taken
1884  sec.usedOrder_ = 1;
1885  sec.psi_[0] = sec.currentTimeStep;
1886  sec.cj_ = 1/sec.psi_[0];
1887  sec.nscsco_ = 0;
1888 }
1889 
1890 //-----------------------------------------------------------------------------
1891 // Function : N_TIA_OneStep::checkReduceOrder()
1892 // Purpose : check whether to reduce order independent of local error test
1893 // Special Notes :
1894 // Scope : public
1895 // Creator : Ting Mei, SNL
1896 // Creation Date : 10/21/07
1897 //-----------------------------------------------------------------------------
1899 {
1900 
1901 }
1902 
1903 //-----------------------------------------------------------------------------
1904 // Function : N_TIA_OneStep::rejectStep()
1905 // Purpose : code to restore history, choose new order/step-size
1906 // Special Notes :
1907 // Scope : public
1908 // Creator : Ting Mei, SNL
1909 // Creation Date : 10/21/07
1910 //-----------------------------------------------------------------------------
1912 {
1913  // This routine puts its output in newTimeStep_ and sec.newOrder_
1914 
1915  // This routine changes the following variables:
1916  // lastAttemptedTimeStep, sec.initialPhase_, sec.nef_, sec.psi_, newTimeStep_,
1917  // sec.newOrder_, sec.currentOrder_, currentTimeStep_, ds.xHistory,
1918  // ds.qHistory, nextTimePt, nextTime, currentTimeStepRatio,
1919  // currentTimeStepSum, nextTimePt
1920 
1921  // This routine reades but does not change the following variables:
1922  // stepAttemptStatus, sec.r_factor_, sec.r_safety_, sec.Est_, sec.r_fudge_, sec.r_min_, sec.r_max_,
1923  // minTimeStep, maxTimeStep, currentTime, stopTime, lastTimeStep
1924 
1926 
1927 #ifdef Xyce_DEBUG_TIME
1928  if (tiaParams.debugLevel > 0)
1929  {
1930  Xyce::dout() << std::endl;
1931  Xyce::dout() << Xyce::section_divider << std::endl;
1932  Xyce::dout() << " N_TIA_OneStep::rejectStep" << std::endl;
1933  }
1934 #endif
1935 
1936  // Only update the time step if we are NOT running constant stepsize.
1937  bool adjustStep = (!(tiaParams.constantStepSize) );
1938 
1940 
1941 
1942  double newTimeStep_ = sec.currentTimeStep;
1943  double rr = 1.0; // step size ratio = new step / old step
1944  if ((sec.stepAttemptStatus == false) && (adjustStep))
1945  {
1946  if (tiaParams.errorAnalysisOption == 1)
1947  {
1948  newTimeStep_ = sec.currentTimeStep/8;
1949  }
1950  else
1951  {
1952  sec.initialPhase_ = false;
1953  sec.nef_++;
1954  restoreHistory();
1955  // restore sec.psi_
1956  // for (int i=1;i<=sec.currentOrder_;++i)
1957  // sec.psi_[i-1] = sec.psi_[i] - sec.currentTimeStep;
1958 
1959  if (sec.nef_ >= sec.max_LET_fail_)
1960  {
1961  std::string msg = "N_TIA_OneStep::rejectStep: ";
1962  msg += " Maximum number of local error test failures. ";
1963  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL_0, msg);
1964  }
1965 
1966  if ((sec.newtonConvergenceStatus <= 0))
1967  {
1968  /// 11/11/05 erkeite: If the Newton solver fails, don't
1969  // rely on the error estimate - it may be full of Nan's.
1970  // rr = sec.r_min_;
1971 
1972  newTimeStep_ = sec.currentTimeStep/8;
1973  // sec.currentOrder_ = 1;
1975  // if (sec.nef_ > 2) sec.newOrder_ = 1;//consistent with block below.
1976  }
1977  else
1978  {
1979  // 03/11/04 tscoffe: Here is the block for choosing order &
1980  // step-size when the local error test FAILS (but Newton
1981  // succeeded).
1982  if (sec.nef_ == 1) // first local error test failure
1983  {
1984  // sec.estOverTol_
1985  sec.Est_ = sec.estOverTol_;
1986 #ifndef Xyce_USE_Q_NORM
1987  rr = sec.r_factor_*pow(sec.r_safety_*(sec.Est_+sec.r_fudge_),-1.0/(sec.newOrder_+1.0));
1988  rr = Xycemax(sec.r_min_,Xycemin(sec.r_max_,rr));
1989 #else
1990  // rr = sec.r_factor_*pow(sec.r_safety_*sec.Est_+sec.r_fudge_,-1.0/(sec.newOrder_+1.0));
1991 
1992  rr = sec.tolAimFac_/(sec.estOverTol_ + 0.0001);
1993  rr = pow(rr, 1.0/(sec.currentOrder_+1.0));
1994  rr = Xycemax(sec.r_min_,Xycemin(sec.r_max_,rr));
1995 #endif
1996  newTimeStep_ = rr * sec.currentTimeStep;
1997 
1998  // sec.currentOrder_ = 1;
2000  }
2001  else // if (sec.nef_ == 2) // second Dae failure
2002  {
2003  rr = sec.r_min_;
2004  newTimeStep_ = rr * sec.currentTimeStep;
2005 
2006  // sec.currentOrder_ = 1;
2008  }
2009  }
2010 #ifdef Xyce_DEBUG_TIME
2011  if (tiaParams.debugLevel > 0)
2012  {
2013  Xyce::dout() << " currentTimeStep = " << sec.currentTimeStep << std::endl;
2014  Xyce::dout() << " numberOfSteps_ = " << sec.numberOfSteps_ << std::endl;
2015  Xyce::dout() << " currentOrder_ = " << sec.currentOrder_ << std::endl;
2016  Xyce::dout() << " nscsco_ = " << sec.nscsco_ << std::endl;
2017  Xyce::dout() << " alpha_[0] = " << sec.alpha_[0] << std::endl;
2018  Xyce::dout() << " alpha_[1] = " << sec.alpha_[1] << std::endl;
2019  Xyce::dout() << " alpha_[2] = " << sec.alpha_[2] << std::endl;
2020  Xyce::dout() << " alpha_[3] = " << sec.alpha_[3] << std::endl;
2021  Xyce::dout() << " alpha_[4] = " << sec.alpha_[4] << std::endl;
2022  Xyce::dout() << " psi_[0] = " << sec.psi_[0] << std::endl;
2023  Xyce::dout() << " psi_[1] = " << sec.psi_[1] << std::endl;
2024  Xyce::dout() << " psi_[2] = " << sec.psi_[2] << std::endl;
2025  Xyce::dout() << " psi_[3] = " << sec.psi_[3] << std::endl;
2026  Xyce::dout() << " psi_[4] = " << sec.psi_[4] << std::endl;
2027  Xyce::dout() << " sigma_[0] = " << sec.sigma_[0] << std::endl;
2028  Xyce::dout() << " sigma_[1] = " << sec.sigma_[1] << std::endl;
2029  Xyce::dout() << " sigma_[2] = " << sec.sigma_[2] << std::endl;
2030  Xyce::dout() << " sigma_[3] = " << sec.sigma_[3] << std::endl;
2031  Xyce::dout() << " sigma_[4] = " << sec.sigma_[4] << std::endl;
2032  Xyce::dout() << " rr = " << rr << std::endl;
2033  Xyce::dout() << " r_factor_ = " << sec.r_factor_ << std::endl;
2034  Xyce::dout() << " r_safety_ = " << sec.r_safety_ << std::endl;
2035  Xyce::dout() << " Est_ = " << sec.Est_ << std::endl;
2036  Xyce::dout() << " nef_ = " << sec.nef_ << std::endl;
2037  Xyce::dout() << " r_fudge_ = " << sec.r_fudge_ << std::endl;
2038  Xyce::dout() << " newOrder_ = " << sec.newOrder_ << std::endl;
2039  Xyce::dout() << " currentTimeStep = " << sec.currentTimeStep << std::endl;
2040  Xyce::dout() << " newTimeStep_ = " << newTimeStep_ << std::endl;
2041  }
2042 #endif // Xyce_DEBUG_TIME
2043  }
2044  }
2045  else if ((sec.stepAttemptStatus == false) & (!adjustStep))
2046  {
2047  std::string tmp = " OneStep:rejectStep: Warning: Local error test failed with constant step-size.\n";
2048  Xyce::dout() << tmp << std::endl;
2049  }
2050 
2051  // If the step needs to be adjusted:
2052  if (adjustStep)
2053  {
2054  newTimeStep_ = Xycemax(newTimeStep_, sec.minTimeStep);
2055  newTimeStep_ = Xycemin(newTimeStep_, sec.maxTimeStep);
2056 
2057  double nextTimePt = sec.currentTime + newTimeStep_;
2058 
2059  if (nextTimePt > sec.stopTime)
2060  {
2061  nextTimePt = sec.stopTime;
2062  newTimeStep_ = sec.stopTime - sec.currentTime;
2064  }
2065 
2066  sec.nextTime = nextTimePt;
2067 
2068  sec.currentTimeStepRatio = newTimeStep_/sec.lastTimeStep;
2069  sec.currentTimeStepSum = newTimeStep_ + sec.lastTimeStep;
2070 
2071 #ifdef Xyce_DEBUG_TIME
2072  if (tiaParams.debugLevel >0)
2073  {
2074  Xyce::dout() << " newTimeStep_ = " << newTimeStep_ << std::endl;
2075  Xyce::dout() << " nextTime = " << sec.nextTime << std::endl;
2076  }
2077 #endif // Xyce_DEBUG_TIME
2078 
2079  sec.currentTimeStep = newTimeStep_;
2080  }
2081  else // if time step is constant for this step:
2082  {
2083  double nextTimePt = sec.currentTime + sec.currentTimeStep;
2084 
2085  if (nextTimePt > sec.stopTime)
2086  {
2087  nextTimePt = sec.stopTime;
2089  }
2090 
2093 
2094  sec.nextTime = nextTimePt;
2095  }
2096 #ifdef Xyce_DEBUG_TIME
2097  if (tiaParams.debugLevel >0)
2098  {
2099  Xyce::dout() << Xyce::section_divider << std::endl;
2100  }
2101 #endif // Xyce_DEBUG_TIME
2102 }
2103 
2104 //-----------------------------------------------------------------------------
2105 // Function : N_TIA_OneStep::rejectStepForHabanero
2106 // Purpose : step rejection, but from an outside program (Habanero API)
2107 // Special Notes :
2108 // Scope : public
2109 // Creator : Eric Keiter, SNL
2110 // Creation Date : 08/11/09
2111 //-----------------------------------------------------------------------------
2113 {
2114  restoreHistory();
2116 }
2117 
2118 //-----------------------------------------------------------------------------
2119 // Function : N_TIA_OneStep::completeStep()
2120 // Purpose : code to update history, choose new order/step-size
2121 // Special Notes :
2122 // Scope : public
2123 // Creator : Ting Mei, SNL
2124 // Creation Date : 10/21/07
2125 //-----------------------------------------------------------------------------
2127 {
2129 
2130  sec.numberOfSteps_ ++;
2131  sec.nef_ = 0;
2134 
2135 #ifdef Xyce_DEBUG_TIME
2136  if (tiaParams.debugLevel > 0)
2137  {
2138  Xyce::dout() << std::endl;
2139  Xyce::dout() << Xyce::section_divider << std::endl;
2140 
2141  Xyce::dout() << " N_TIA_OneStep::completeStep" << std::endl;
2142  }
2143 #endif
2144 
2145  // Only update the time step if we are NOT running constant stepsize.
2146  bool adjustStep = (!(tiaParams.constantStepSize) );
2147 
2149 
2150  double newTimeStep_ = sec.currentTimeStep;
2151  double rr = 1.0; // step size ratio = new step / old step
2152  // 03/11/04 tscoffe: Here is the block for choosing order & step-size when
2153  // the local error test PASSES (and Newton succeeded).
2154 
2155  // need the lastTimeStep for output interpolation. Save it before
2156  // overwriting it with currentTimeStep;
2157  // This history for this step is saved in updateHistory(), called below
2159 
2161  sec.lastTimeStepRatio = sec.currentTimeStepRatio; // copied from calcTStep1
2162  sec.lastTimeStepSum = sec.currentTimeStepSum; // copied from calcTStep1
2163  int orderDiff = sec.currentOrder_ - sec.usedOrder_;
2166 
2167 #ifdef Xyce_DEBUG_TIME
2168  if (tiaParams.debugLevel > 0)
2169  {
2170  Xyce::dout() << " initialPhase_ = " << sec.initialPhase_ << std::endl;
2171  Xyce::dout() << " rr = " << rr << std::endl;
2172  Xyce::dout() << " currentTimeStep = " << sec.currentTimeStep << std::endl;
2173  Xyce::dout() << " currentTime = " << sec.currentTime << std::endl;
2174  Xyce::dout() << " nextTime = " << sec.nextTime << std::endl;
2175  Xyce::dout() << " newTimeStep_ = " << newTimeStep_ << std::endl;
2176  Xyce::dout() << " minTimeStep = " << sec.minTimeStep << std::endl;
2177  Xyce::dout() << " maxTimeStep = " << sec.maxTimeStep << std::endl;
2178  }
2179 #endif
2180  // 03/22/04 tscoffe: Note that updating the history has nothing to do with
2181  // the step-size and everything to do with the newton correction vectors.
2182 
2183  /* if (sec.numberOfSteps_ >= 2)
2184  {
2185  sec.currentOrder_ = 2;
2186  // (ds.relErrTolPtr)->putScalar(1e-2);
2187  }
2188  */
2189  if (tiaParams.errorAnalysisOption == 1)
2190  {
2191  if (sec.numberOfSteps_ >= 2 && sec.maxOrder_ == 2)
2192  {
2193  sec.currentOrder_ = 2;
2194  }
2195 
2196  rr = 1;
2197 
2198  if (sec.nIterations <= tiaParams.NLmin)
2199  rr = 2;
2200 
2202  rr = 1.0/8;
2203 
2204  newTimeStep_ = rr*sec.currentTimeStep;
2205  }
2206  else
2207  {
2208  rr = sec.tolAimFac_/(sec.estOverTol_ + 0.0001);
2209  rr = pow(rr, 1.0/(sec.currentOrder_+1.0));
2210 
2211  if (sec.numberOfSteps_ >= 2 && sec.maxOrder_ == 2)
2212  {
2213  if (sec.currentOrder_ == 1)
2214  {
2215  sec.currentOrder_ = 2;
2216  rr = sec.tolAimFac_/(sec.estOverTol_ + 0.0001);
2217  rr = pow(rr, 1.0/(sec.currentOrder_+1.0));
2218 
2219  if (rr <= 1.05)
2220  {
2221  // sec.currentOrder_ = 1;
2223  }
2224  }
2225  }
2226 #ifdef Xyce_DEBUG_TIME
2227  if (tiaParams.debugLevel > 0)
2228  {
2229  Xyce::dout() << " currentOrder_ = " << sec.currentOrder_ << std::endl;
2230  Xyce::dout() << " r_safety = " << sec.r_safety_ << std::endl;
2231  Xyce::dout() << " r_fudge_ = " << sec.r_fudge_ << std::endl;
2232  Xyce::dout() << " r_hincr_ = " << sec.r_hincr_ << std::endl;
2233  Xyce::dout() << " r_hincr_test_ = " << sec.r_hincr_test_ << std::endl;
2234  Xyce::dout() << " Est = " << sec.Est_ << std::endl;
2235  Xyce::dout() << " raw rr = " << rr << std::endl;
2236  }
2237 #endif
2238 
2239  if (rr >= sec.r_hincr_test_)
2240  {
2241  rr = sec.r_hincr_;
2242  newTimeStep_ = rr*sec.currentTimeStep;
2243  }
2244  else if (rr <= 1)
2245  {
2246  rr = Xycemax(sec.r_min_,Xycemin(sec.r_max_,rr));
2247  newTimeStep_ = rr*sec.currentTimeStep;
2248  }
2249  }
2250 
2251  updateHistory();
2252 
2253  newTimeStep_ = Xycemax(newTimeStep_, sec.minTimeStep);
2254  newTimeStep_ = Xycemin(newTimeStep_, sec.maxTimeStep);
2255 
2256  // 12/01/05 tscoffe: This is necessary to avoid currentTimeStep == 0 right
2257  // before a breakpoint. So I'm checking to see if currentTime is identically
2258  // equal to stopTime, in which case we are right before a breakpoint and we
2259  // should not adjust currentStepSize because that would result in
2260  // currentStepSize == 0.
2262  {
2263  // If the step needs to be adjusted:
2264  if (adjustStep)
2265  {
2266  // newTimeStep_ = Xycemax(newTimeStep_, sec.minTimeStep);
2267  // newTimeStep_ = Xycemin(newTimeStep_, sec.maxTimeStep);
2268 
2269  double nextTimePt = sec.currentTime + newTimeStep_;
2270  if (nextTimePt > sec.stopTime)
2271  {
2272 
2273  sec.savedTimeStep = newTimeStep_;
2274 
2275  nextTimePt = sec.stopTime;
2276 
2277  newTimeStep_ = sec.stopTime - sec.currentTime;
2279  }
2280 
2281  sec.nextTime = nextTimePt;
2282 
2283  sec.currentTimeStepRatio = newTimeStep_/sec.lastTimeStep;
2284  sec.currentTimeStepSum = newTimeStep_ + sec.lastTimeStep;
2285 
2286 #ifdef Xyce_DEBUG_TIME
2287  if (tiaParams.debugLevel >0)
2288  {
2289  Xyce::dout() << " nextTime = " << sec.nextTime << std::endl;
2290  Xyce::dout() << " newTimeStep_ = " << newTimeStep_ << std::endl;
2291  }
2292 #endif
2293 
2294  // sec.currentTimeStep = newTimeStep_;
2295  }
2296  else // if time step is constant for this step:
2297  {
2298  double nextTimePt = sec.currentTime + sec.currentTimeStep;
2299 
2300  if (nextTimePt > sec.stopTime)
2301  {
2303 
2304  nextTimePt = sec.stopTime;
2306  }
2307 
2310 
2311  sec.nextTime = nextTimePt;
2312  }
2313  }
2314 
2315 #ifdef Xyce_DEBUG_TIME
2316  if (tiaParams.debugLevel >0)
2317  {
2318  Xyce::dout() << Xyce::section_divider << std::endl;
2319  }
2320 #endif
2321 
2322  sec.currentTimeStep = newTimeStep_;
2323 }
2324 
2325 //-----------------------------------------------------------------------------
2326 // Function : N_TIA_OneStep::updateStateDeriv
2327 // Purpose :
2328 // Special Notes :
2329 // Scope : public
2330 // Creator : Ting Mei, SNL
2331 // Creation Date : 11/17/05
2332 //-----------------------------------------------------------------------------
2334 {
2335  // dS/dt = spn0 - (sec.alpha_/hn)(S(x)-sn0)
2336  ds.nextStateDerivPtr->linearCombo(1.0,*ds.nextStatePtr,-1.0,*ds.sn0Ptr);
2337  // ds.nextStateDerivPtr->linearCombo(1.0,*ds.spn0Ptr,-sec.alphas_/sec.currentTimeStep,*ds.nextStateDerivPtr);
2338 
2339  if (sec.currentOrder_ == 1)
2340  {
2342  }
2343  else
2344  {
2346  }
2347 
2348 #ifdef Xyce_DEBUG_TIME
2349  if (tiaParams.debugLevel > 1)
2350  {
2351  Xyce::dout() << "\n next state Ptr: \n" << std::endl;
2352  ds.nextStatePtr->printPetraObject(Xyce::dout());
2353  Xyce::dout() << std::endl;
2354 
2355  Xyce::dout() << "\n sn0: \n" << std::endl;
2356  ds.sn0Ptr->printPetraObject(Xyce::dout());
2357  Xyce::dout() << std::endl;
2358 
2359  Xyce::dout() << "\n next State Deriv: \n" << std::endl;
2360  ds.nextStateDerivPtr->printPetraObject(Xyce::dout());
2361  Xyce::dout() << std::endl;
2362  }
2363 #endif
2364 }
2365 
2366 //-----------------------------------------------------------------------------
2367 // Function : N_TIA_OneStep::updateLeadCurrent
2368 // Purpose : calculates lead currents in store vector with
2369 // the storeLeadCurrQVec.
2370 // Special Notes :
2371 // Scope : public
2372 // Creator : Richard Schiek, Electrical Systems Modeling, SNL
2373 // Creation Date : 03/22/2013
2374 //-----------------------------------------------------------------------------
2376 {
2378  if (sec.currentOrder_ == 1)
2379  {
2381  }
2382  else
2383  {
2385  linearCombo(-1.0,*(ds.currStoreLeadCurrQDerivPtr),
2387  }
2388  ds.nextStorePtr->linearCombo(1.0,*ds.nextStorePtr,1.0,*ds.nextStoreLeadCurrQDerivPtr);
2389 
2390 #ifdef Xyce_DEBUG_TIME
2391  if (tiaParams.debugLevel > 1)
2392  {
2393  Xyce::dout() << "\n next store Ptr: \n" << std::endl;
2394  ds.nextStorePtr->printPetraObject(Xyce::dout());
2395  Xyce::dout() << std::endl;
2396  }
2397 #endif
2398 
2399 }
2400 
2401 //-----------------------------------------------------------------------------
2402 // Function : N_TIA_OneStep::getInitialQnorm
2403 // Purpose : Needed by 2-level solves.
2404 // Special Notes :
2405 // Scope : public
2406 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
2407 // Creation Date : 03/18/07
2408 //-----------------------------------------------------------------------------
2410 {
2411  tle.q1HistorySum = ds.partialSum_q1();
2412 }
2413 
2414 //-----------------------------------------------------------------------------
2415 // Function : N_TIA_OneStep::setupTwoLevelError
2416 // Purpose : Needed by 2-level solves.
2417 // Special Notes :
2418 // Scope : public
2419 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
2420 // Creation Date : 03/15/07
2421 //-----------------------------------------------------------------------------
2423 {
2424  tle.xErrorSum = ds.partialErrorNormSum ();
2429  tle.innerSize = ds.globalLength ();
2430 
2431 #ifdef Xyce_DEBUG_TIME
2432  Xyce::dout() << tle;
2433 #endif
2434 }
2435