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