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