Xyce  6.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
N_ANP_HB.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------------
2 // Copyright Notice
3 //
4 // Copyright 2002 Sandia Corporation. Under the terms
5 // of Contract DE-AC04-94AL85000 with Sandia Corporation, the U.S.
6 // Government retains certain rights in this software.
7 //
8 // Xyce(TM) Parallel Electrical Simulator
9 // Copyright (C) 2002-2014 Sandia Corporation
10 //
11 // This program is free software: you can redistribute it and/or modify
12 // it under the terms of the GNU General Public License as published by
13 // the Free Software Foundation, either version 3 of the License, or
14 // (at your option) any later version.
15 //
16 // This program is distributed in the hope that it will be useful,
17 // but WITHOUT ANY WARRANTY; without even the implied warranty of
18 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 // GNU General Public License for more details.
20 //
21 // You should have received a copy of the GNU General Public License
22 // along with this program. If not, see <http://www.gnu.org/licenses/>.
23 //-----------------------------------------------------------------------------
24 
25 //-----------------------------------------------------------------------------
26 // Filename : $RCSfile: N_ANP_HB.C,v $
27 // Purpose : HB analysis functions.
28 // Special Notes :
29 // Creator : Todd Coffey, 1414, Ting Mei, 1437
30 // Creation Date : 07/23/08
31 //
32 // Revision Information:
33 // ---------------------
34 // Revision Number: $Revision: 1.94.2.2 $
35 // Revision Date : $Date: 2014/08/28 21:00:43 $
36 // Current Owner : $Author: dgbaur $
37 //-----------------------------------------------------------------------------
38 #include <Xyce_config.h>
39 
40 #include <N_UTL_APFT.h>
41 
42 #include <N_ANP_AnalysisManager.h>
43 #include <N_ANP_HB.h>
44 #include <N_ANP_Transient.h>
45 #include <N_ANP_DCSweep.h>
46 #include <N_ANP_Report.h>
47 #include <N_MPDE_Manager.h>
48 #include <N_MPDE_Discretization.h>
49 
50 #include <N_TIA_StepErrorControl.h>
52 #include <N_TIA_DataStore.h>
53 
54 #include <N_LOA_HBLoader.h>
55 #include <N_LOA_NonlinearEquationLoader.h>
56 #include <N_LAS_HBBuilder.h>
57 #include <N_LAS_HBPrecondFactory.h>
58 #include <N_LAS_PrecondFactory.h>
59 #include <N_LAS_System.h>
60 #include <N_LAS_BlockSystemHelpers.h>
61 
62 #include <N_NLS_Manager.h>
63 
64 // #include <N_IO_OutputMgr.h>
65 #include <N_IO_ActiveOutput.h>
66 
67 #include <N_UTL_FFTInterface.hpp>
68 #include <N_UTL_Timer.h>
69 
70 #include <Teuchos_BLAS.hpp>
71 #include <Teuchos_Utils.hpp>
72 #include <Teuchos_ScalarTraits.hpp>
73 #include <Teuchos_SerialDenseMatrix.hpp>
74 #include <Teuchos_SerialDenseVector.hpp>
75 #include <Teuchos_SerialDenseHelpers.hpp>
76 #include <Teuchos_SerialDenseSolver.hpp>
77 
78 using Teuchos::rcp;
79 using Teuchos::RCP;
80 using Teuchos::rcp_dynamic_cast;
81 
82 namespace Xyce {
83 namespace Analysis {
84 
85 //-----------------------------------------------------------------------------
86 // Function : HB::HB( AnalysisManager * )
87 // Purpose :
88 // Special Notes :
89 // Scope : public
90 // Creator : Todd Coffey, SNL
91 // Creation Date : 09/18/2008
92 //-----------------------------------------------------------------------------
93 HB::HB( AnalysisManager &analysis_manager )
94  : AnalysisBase(analysis_manager),
95  Util::ListenerAutoSubscribe<StepEvent>(&analysis_manager),
96  debugLevel(0),
97  isPaused(false),
98  startDCOPtime(0.0),
99  endTRANtime(0.0),
100  isTransient_(false),
101  isDCSweep_(false),
102  test_(false),
103  size_(21),
104  period_(1.0),
105  startUpPeriods_(0),
106  startUpPeriodsGiven_(false),
107  startUpPeriodsFinished_(false),
108  saveIcData_(false),
109  transTiaParams_( AnalysisBase::tiaParams_),
110  voltLimFlag_(1),
111  intmodMax_(0),
112  method_("APFT"),
113  taHB_(1),
114  intmodMaxGiven_(false),
115  fastTimeDisc_(0),
116  fastTimeDiscOrder_(1),
117  resetForStepCalledBefore_(false),
118  hbLoaderPtr_(0)
119 {
124 }
125 
127 {
128  delete hbLoaderPtr_;
129 }
130 
131 void HB::notify(const StepEvent &event)
132 {
133  if (event.state_ == StepEvent::STEP_STARTED)
134  {
136 
138  {
139  goodSolutionVec_.clear();
140  goodStateVec_.clear();
141  goodQVec_.clear();
142  goodStoreVec_.clear();
143 
145 
149 
155 
157 
158  // un-set the fast source flag in the devices
159  std::vector<std::string> srcVec;
161 
163 
165  devInterfacePtr_->setMPDEFlag( false );
166 
168 
170 
171  analysisManager_.getXyceTranTimer().resetStartTime();
172  }
173 
175  }
176 }
177 
178 //-----------------------------------------------------------------------------
179 // Function : HB::getStepNumber()
180 // Purpose :
181 // Special Notes :
182 // Scope : public
183 // Creator : Heidi Thornquist
184 // Creation Date : 5/20/13
185 //-----------------------------------------------------------------------------
187 {
188  if ( !Teuchos::is_null( analysisObject_ ) )
189  {
190  return analysisObject_->getStepNumber();
191  }
192  return 0;
193 }
194 
195 //-----------------------------------------------------------------------------
196 // Function : HB::setStepNumber()
197 // Purpose :
198 // Special Notes :
199 // Scope : public
200 // Creator : Heidi Thornquist
201 // Creation Date : 5/20/13
202 //-----------------------------------------------------------------------------
203 void HB::setStepNumber (int step)
204 {
205  if ( !Teuchos::is_null( analysisObject_ ) )
206  {
207  analysisObject_->setStepNumber( step );
208  }
209 }
210 
211 //-----------------------------------------------------------------------------
212 // Function : HB::setBeginningIntegrationFlag()
213 // Purpose :
214 // Special Notes :
215 // Scope : public
216 // Creator : Heidi Thornquist
217 // Creation Date : 5/20/13
218 //-----------------------------------------------------------------------------
220 {
221  if ( !Teuchos::is_null( analysisObject_ ) )
222  {
223  analysisObject_->setBeginningIntegrationFlag( bif );
224  }
225 }
226 
227 //-----------------------------------------------------------------------------
228 // Function : HB::getBeginningIntegrationFlag()
229 // Purpose :
230 // Special Notes :
231 // Scope : public
232 // Creator : Heidi Thornquist
233 // Creation Date : 5/20/13
234 //-----------------------------------------------------------------------------
236 {
237  if ( !Teuchos::is_null( analysisObject_ ) )
238  {
239  return analysisObject_->getBeginningIntegrationFlag();
240  }
241  return true;
242 }
243 
244 //-----------------------------------------------------------------------------
245 // Function : HB::setIntegrationMethod
246 // Purpose :
247 // Special Notes :
248 // Scope : public
249 // Creator : Heidi Thornquist
250 // Creation Date : 5/20/13
251 //-----------------------------------------------------------------------------
253 {
254  if ( !Teuchos::is_null( analysisObject_ ) )
255  {
256  analysisObject_->setIntegrationMethod( im );
257  }
258 }
259 
260 //-----------------------------------------------------------------------------
261 // Function : HB::getIntegrationMethod
262 // Purpose :
263 // Special Notes :
264 // Scope : public
265 // Creator : Heidi Thornquist
266 // Creation Date : 5/20/13
267 //-----------------------------------------------------------------------------
269 {
270  if ( !Teuchos::is_null( analysisObject_ ) )
271  {
272  return analysisObject_->getIntegrationMethod ();
273  }
274  return TIAMethod_NONE;
275 }
276 
277 
278 //-----------------------------------------------------------------------------
279 // Function : HB::getDoubleDCOPStep
280 // Purpose :
281 // Special Notes :
282 // Scope : public
283 // Creator : Eric Keiter
284 // Creation Date : 3/24/2014
285 //-----------------------------------------------------------------------------
287 {
288  if ( !Teuchos::is_null( analysisObject_ ) )
289  {
290  return analysisObject_->getDoubleDCOPStep();
291  }
292  else
293  {
294  return 0;
295  }
296 }
297 
298 //-----------------------------------------------------------------------------
299 // Function : HB::getDCOPFlag
300 // Purpose :
301 // Special Notes :
302 // Scope : public
303 // Creator : Eric Keiter
304 // Creation Date : 3/24/2014
305 //-----------------------------------------------------------------------------
307 {
308  if ( !Teuchos::is_null( analysisObject_ ) )
309  {
310  if ( isTransient_ )
311  {
312  return analysisObject_->getDCOPFlag ();
313  }
314  // DCSweep is a special case in the HB analysis type. The HB calculation is
315  // performed through a DC analysis object. However, this function (getDCOPFlag)
316  // is called (ultimately) from the device package to determine a bunch of
317  // state-dependent load decisions. For an HB calculation, it should NOT
318  // do a DCOP load. It needs to do a full transient load.
319  else
320  {
321  return false;
322  }
323  }
324  else
325  {
326  // the assumption here is that if there is no analysis object, then
327  // the simulation hasn't started yet. When it *does* start, the first
328  // analysis will be DCOP.
329  //return true;
330  return false;
331  }
332 }
333 
334 //-----------------------------------------------------------------------------
335 // Function : HB::run()
336 // Purpose :
337 // Special Notes :
338 // Scope : public
339 // Creator : Todd Coffey, SNL
340 // Creation Date : 09/18/2008
341 //-----------------------------------------------------------------------------
342 bool HB::run()
343 {
344 
345  // initializeAll_();
346  // get TIAParams from AnalysisManager
347  // create HBBuilder
348  // generateMaps
349  // generateStateMaps
350  // create vectors for IC
351  //
352  // if (test_) {
353  // runTests_();
354  // } else {
355  //
356  // computeInitialCondition_();
357  // run startup periods
358  // integrate one period
359  // interpolate to evenly spaced points
360  // set this as IC for HB nonlinear problem
361  //
362  // setupHBProblem_();
363  //
364  // runHBProblem_();
365  //
366  // }
367 
368  bool bsuccess = true;
369 
370  bsuccess = bsuccess & init();
371  bsuccess = bsuccess & loopProcess();
372 
373  // if processing the loop failed,
374  // then skip finish step
375  if( bsuccess )
376  {
377  bsuccess = bsuccess & finish();
378  }
379 
380  return bsuccess;
381 }
382 
383 //-----------------------------------------------------------------------------
384 // Function : HB::init()
385 // Purpose :
386 // Special Notes :
387 // Scope : public
388 // Creator : Todd Coffey, SNL
389 // Creation Date : 09/18/2008
390 //-----------------------------------------------------------------------------
391 bool HB::init()
392 {
393  bool returnValue=true;
394 
395  Xyce::lout() << " ***** Running HB initial conditions....\n" << std::endl;
396 #ifdef Xyce_DEBUG_HB
397  Xyce::dout() << std::endl
398  << section_divider << std::endl
399  << " HB::init()" << std::endl;
400 #endif // Xyce_DEBUG_HB
401 
402  //Store copy of transient TIAParams for HB run
404  doubleDCOPFlag_ = loader_.getDoubleDCOPFlag();
406 
407  setFreqPoints_();
408  //
409  // If it was requested, advance the solution a fixed number of startup periods
410 
411  period_ = 1.0/freqPoints_[ (size_-1)/2 + 1 ];
412 
413 #ifdef Xyce_DEBUG_HB
414  Xyce::dout() << "HB period =" << period_ << std::endl;
415 #endif
416 
417  if (transTiaParams_.freqs.size() == 1)
419  else
420  {
421  setTimePoints_();
422 
423  fastTimes_.resize(size_+1);
424  goodTimePoints_.resize(size_+1);
425 
427 
428  }
429 
430  // now that we have size_, continue with the initialization of objects for HB
431  mpdeDiscPtr_ = rcp(new N_MPDE_Discretization(
432  static_cast<N_MPDE_Discretization::Type>(fastTimeDisc_), fastTimeDiscOrder_ ));
433 
434  hbBuilderPtr_ = rcp(new N_LAS_HBBuilder( size_, mpdeDiscPtr_ ));
435 
436 #ifdef Xyce_DEBUG_HB
437  Xyce::dout() << "HB::init(): Generate Maps\n";
438 #endif // Xyce_DEBUG_HB
439  hbBuilderPtr_->generateMaps( rcp(pdsMgrPtr_->getParallelMap( "SOLUTION" ), false),
440  rcp(pdsMgrPtr_->getParallelMap( "SOLUTION_OVERLAP_GND" ), false) );
441  hbBuilderPtr_->generateStateMaps( rcp(pdsMgrPtr_->getParallelMap( "STATE" ),false) );
442  hbBuilderPtr_->generateStoreMaps( rcp(pdsMgrPtr_->getParallelMap( "STORE" ),false) );
443  hbBuilderPtr_->generateGraphs( *pdsMgrPtr_->getMatrixGraph( "JACOBIAN" ));
444 
445  HBICVectorPtr_ = hbBuilderPtr_->createTimeDomainBlockVector();
446  HBICStateVectorPtr_ = hbBuilderPtr_->createTimeDomainStateBlockVector();
447  HBICStoreVectorPtr_ = hbBuilderPtr_->createTimeDomainStoreBlockVector();
448  HBICQVectorPtr_ = hbBuilderPtr_->createTimeDomainBlockVector();
449 
450  HBICVectorFreqPtr_ = hbBuilderPtr_->createExpandedRealFormTransposeBlockVector();
451 
452 // HBICStateVectorFreqPtr_ = hbBuilderPtr_->createExpandedRealFormTransposeStateBlockVector();
453 // HBICQVectorFreqPtr_ = hbBuilderPtr_->createExpandedRealFormTransposeBlockVector();
454 // HBICStoreVectorFreqPtr_ = hbBuilderPtr_->createExpandedRealFormTransposeStoreBlockVector();
455 
458 
459  // set the fast source flag in the devices
460  std::vector<std::string> srcVec;
462 
463  //analysisManager_.anaIntPtr->getnlHBOptions(saved_nlHBOB_);
465 
466  // Create HB Loader.
467  delete hbLoaderPtr_;
468  hbLoaderPtr_ = new N_LOA_HBLoader( mpdeState_, mpdeDiscPtr_ );
469  hbLoaderPtr_->registerHBBuilder(hbBuilderPtr_);
470  hbLoaderPtr_->registerAppBuilder(appBuilderPtr_);
471 
472  // Create DFT for HB Loader
473  // NOTE: For single-tone HB the DFT will probably be a FFT, for multi-tone a specialized
474  // implementation of the Util::DFTInterfaceDecl will need to be made and registered with
475  // the HB loader.
476 
477  ftInData_.resize( size_ );
478  ftOutData_.resize( size_ +1 );
479  iftInData_.resize( size_ +1 );
480  iftOutData_.resize( size_ );
481  if (transTiaParams_.freqs.size() == 1)
482  {
483  if (ftInterface_ == Teuchos::null)
484  {
485  ftInterface_ = Teuchos::rcp( new N_UTL_FFTInterface<std::vector<double> >( size_ ) );
486  ftInterface_->registerVectors( ftInData_, &ftOutData_, iftInData_, &iftOutData_ );
487  }
488  else if (ftInterface_->getFFTInterface()->getSignalLength() != size_)
489  {
490  ftInterface_ = Teuchos::rcp( new N_UTL_FFTInterface<std::vector<double> >( size_ ) );
491  ftInterface_->registerVectors( ftInData_, &ftOutData_, iftInData_, &iftOutData_ );
492  }
493  hbLoaderPtr_->registerDFTInterface( ftInterface_->getFFTInterface() );
494  }
495  else
496  {
497  createFT_();
498 
499  dftInterface_ = Teuchos::rcp( new N_UTL_APFT<std::vector<double> >( idftMatrix_, dftMatrix_ ) );
500  dftInterface_->registerVectors( Teuchos::rcp( &ftInData_, false ), Teuchos::rcp( &ftOutData_, false ), Teuchos::rcp( &iftInData_, false ), Teuchos::rcp( &iftOutData_, false ) );
501 // ftInterface_-> registerMatrices(idftMatrix_, dftMatrix_);
502 
503  hbLoaderPtr_->registerDFTInterface( dftInterface_ );
504  }
505 
506  if (taHB_==1)
507  {
508  // Pick up IC data from the initial transient.
509  for (int i=0 ; i<size_ ; ++i)
510  {
511 #ifdef Xyce_DEBUG_HB
512  if( debugLevel > 0 )
513  {
514  Xyce::dout() << "HB::init(): Loading initial condition data from time: fastTimes_["
515  << i << "] = " << fastTimes_[i] << std::endl;
516  }
517 #endif // Xyce_DEBUG_HB
518 
519  HBICVectorPtr_->block(i) = *(goodSolutionVec_[i]);
520  HBICStateVectorPtr_->block(i) = *(goodStateVec_[i]);
521  HBICQVectorPtr_->block(i) = *(goodQVec_[i]);
522  HBICStoreVectorPtr_->block(i) = *(goodStoreVec_[i]);
523  }
524 
526 
527  }
528 #ifdef Xyce_DEBUG_HB
529  if ( debugLevel > 1 )
530  {
531  Xyce::dout() << "HB Initial Condition Solution!\n";
532  HBICVectorPtr_->printPetraObject(std::cout);
533  Xyce::dout() << "HB Initial Condition State Vector!\n";
534  HBICStateVectorPtr_->printPetraObject(std::cout);
535  Xyce::dout() << "HB Initial Condition Store Vector!\n";
536  HBICStoreVectorPtr_->printPetraObject(std::cout);
537  }
538 #endif // Xyce_DEBUG_HB
539 
540  //Destroy Solvers, etc. from IC phase and prepare for HB
541  //-----------------------------------------
542 
543 // if ( voltLimFlag_ == 0 )
544 // devInterfacePtr_->setVoltageLimiterFlag (false);
546 
547  devInterfacePtr_->setMPDEFlag( true );
549 
550  //-----------------------------------------
551 
552  //Finish setup of HB Loader
553  //-----------------------------------------
554  hbLoaderPtr_->registerAppLoader( rcp(&loader_, false) );
555  hbLoaderPtr_->registerDeviceInterface( devInterfacePtr_ );
556  goodTimePoints_.resize(size_+1);
558 
559  for( int i = 0; i < size_; ++i )
560  {
562  }
563 
564  hbLoaderPtr_->setFastTimes(goodTimePoints_);
565 
566  //-----------------------------------------
567  //Construct Solvers, etc. for HB Phase
568  //-----------------------------------------
569  lasHBSysPtr_ = rcp(new N_LAS_System());
570  //-----------------------------------------
571 
574 
575  //hack needed by TIA initialization currently
576  hbBuilderPtr_->registerPDSManager( pdsMgrPtr_ );
577 
578  lasHBSysPtr_->registerAnalysisManager(&analysisManager_);
579  lasHBSysPtr_->registerPDSManager( pdsMgrPtr_ );
580  lasHBSysPtr_->registerBuilder( &*hbBuilderPtr_ );
581 
582  //need to cut out unnecessary stuff from this call for new dae
583  lasHBSysPtr_->initializeSystem();
584 
585  // Give NLS Manager the same old nonlinearEquationLoader as it just calls the TIA loader in newDAE
589 
590  // Let the HB loader know that the application of the operator is matrix free
591  hbLoaderPtr_->setMatrixFreeFlag( true );
592 
593  if (Teuchos::is_null( precFactory_ ))
594  {
595  // Generate the HB preconditioner factory.
596  precFactory_ = rcp( new N_LAS_HBPrecondFactory( saved_lsHBOB_ ) );
597  }
598 
599  // Register application loader with preconditioner factory
600  RCP<N_LAS_HBPrecondFactory> tmpPrecFactory
601  = rcp_dynamic_cast<N_LAS_HBPrecondFactory>( precFactory_ );
602 
603  tmpPrecFactory->registerAppBuilder( appBuilderPtr_ );
604  tmpPrecFactory->registerHBLoader( rcp(hbLoaderPtr_, false) );
605  tmpPrecFactory->registerHBBuilder( hbBuilderPtr_ );
606  tmpPrecFactory->setFastTimes( goodTimePoints_ );
607  tmpPrecFactory->setTimeSteps( timeSteps_ );
608 
610  //-----------------------------------------
611 
612  //Initialization of Solvers, etc. for HB Phase
613  //-----------------------------------------
614  //Dummy call to setup time integrator for transient
618 
621 
622 #ifdef Xyce_DEBUG_HB
623  Xyce::dout() << section_divider << std::endl;
624 #endif // Xyce_DEBUG_HB
625 
626  return returnValue;
627 }
628 
629 //-----------------------------------------------------------------------------
630 // Function : HB::loopProcess()
631 // Purpose :
632 // Special Notes :
633 // Scope : public
634 // Creator : Todd Coffey, SNL
635 // Creation Date : 09/18/2008
636 //-----------------------------------------------------------------------------
638 {
639  bool returnValue = true;
640 
641  Xyce::lout() << " ***** Beginning full HB simulation....\n" << std::endl;
642 
643 #ifdef Xyce_DEBUG_HB
644  Xyce::dout() << std::endl
645  << section_divider << std::endl
646  << " HB::loopProcess" << std::endl;
647 #endif // Xyce_DEBUG_HB
648 
650  *(dsPtr->nextSolutionPtr) = *(HBICVectorFreqPtr_.get());
651  *(dsPtr->nextStatePtr) = *(HBICStateVectorPtr_.get());
652  *(dsPtr->nextStorePtr) = *(HBICStoreVectorPtr_.get());
653 
654  // try to run the problem
655  analysisObject_ = Teuchos::rcp(new DCSweep(analysisManager_));
656  returnValue = analysisObject_->run();
657 
658  // Add in simulation times
660 
661  // print out analysis info
662  Xyce::lout() << " ***** Harmonic Balance Computation Summary *****" << std::endl;
663  analysisObject_->printLoopInfo( 0, 0 );
664 
665 #ifdef Xyce_DEBUG_HB
666  dout() << section_divider << std::endl;
667 #endif // Xyce_DEBUG_HB
668 
669  return returnValue;
670 }
671 
672 //-----------------------------------------------------------------------------
673 // Function : HB::processSuccessfulDCOP()
674 // Purpose :
675 // Special Notes :
676 // Scope : public
677 // Creator : Todd Coffey, SNL
678 // Creation Date : 09/18/2008
679 //-----------------------------------------------------------------------------
681 {
682  return false;
683 }
684 
685 //-----------------------------------------------------------------------------
686 // Function : HB::processSuccessfulStep()
687 // Purpose :
688 // Special Notes :
689 // Scope : public
690 // Creator : Todd Coffey, SNL
691 // Creation Date : 09/18/2008
692 //-----------------------------------------------------------------------------
694 {
695  return false;
696 }
697 
698 //-----------------------------------------------------------------------------
699 // Function : HB::processFailedStep
700 // Purpose :
701 // Special Notes :
702 // Scope : public
703 // Creator : Todd Coffey, SNL
704 // Creation Date : 09/18/2008
705 //-----------------------------------------------------------------------------
707 {
708  return false;
709 }
710 
711 //-----------------------------------------------------------------------------
712 // Function : HB::processFailedDCOP
713 // Purpose :
714 // Special Notes :
715 // Scope : public
716 // Creator : Todd Coffey, SNL
717 // Creation Date : 09/18/2008
718 //-----------------------------------------------------------------------------
720 {
721  return false;
722 }
723 
724 //-----------------------------------------------------------------------------
725 // Function : HB::finish
726 // Purpose :
727 // Special Notes :
728 // Scope : public
729 // Creator : Todd Coffey, SNL
730 // Creation Date : 09/18/2008
731 //-----------------------------------------------------------------------------
733 {
735 
736  return true;
737 }
738 
740 {
741  return true;
742 }
743 
744 //-----------------------------------------------------------------------------
745 // Function : HB::finalVerboseOutput
746 // Purpose :
747 // Special Notes :
748 // Scope : public
749 // Creator : Todd Coffey, SNL
750 // Creation Date : 09/18/2008
751 //-----------------------------------------------------------------------------
753 {
754  return false;
755 }
756 
757 //-----------------------------------------------------------------------------
758 // Function : HB::setHBOptions
759 // Purpose :
760 // Special Notes :
761 // Scope : public
762 // Creator : Heidi Thornquist
763 // Creation Date : 05/13/13
764 //-----------------------------------------------------------------------------
765 bool HB::setHBOptions(const Util::OptionBlock & OB)
766 {
767  Util::ParameterList::const_iterator iterPL = OB.getParams().begin();
768  Util::ParameterList::const_iterator endPL = OB.getParams().end();
769 
770  for( ; iterPL != endPL; ++iterPL )
771  {
772  ExtendedString tag = iterPL->tag();
773  tag.toUpper();
774 
775  if (std::string(tag,0,7) == "NUMFREQ" )
776  {
777  size_ = iterPL->getImmutableValue<int>();
778  numFreqs_.push_back(size_);
779  }
780  else if ( tag == "STARTUPPERIODS" )
781  {
782  startUpPeriods_ = iterPL->getImmutableValue<int>();
783 
784  if (startUpPeriods_ > 0)
785  startUpPeriodsGiven_ = true;
786  }
787  else if( tag == "SAVEICDATA" )
788  {
789  saveIcData_ = true;
790  }
791  else if( tag == "TEST" )
792  {
793  test_ = static_cast<bool> (iterPL->getImmutableValue<int>());
794  }
795  else if (tag == "DEBUGLEVEL" )
796  {
797  debugLevel = iterPL->getImmutableValue<int>();
798  }
799  else if ( tag == "TAHB" )
800  {
801  taHB_ = iterPL->getImmutableValue<int>();
802  }
803  else if ( tag == "VOLTLIM" )
804  {
805  voltLimFlag_ = static_cast<bool> (iterPL->getImmutableValue<int>());
806  }
807  else if ( tag == "INTMODMAX" )
808  {
809  intmodMax_ = iterPL->getImmutableValue<int>();
810 
811  if ( intmodMax_ > 0)
812  intmodMaxGiven_ = true;
813  }
814  else if ( tag == "METHOD" )
815  {
816  ExtendedString stringVal ( iterPL->stringValue() );
817  stringVal.toUpper();
818  method_ = stringVal;
819  }
820  else
821  {
822  UserWarning(*this) << "Unrecognized HBINT option " << tag;
823  }
824  }
825 
826  if (numFreqs_.size() != 0 && numFreqs_.size() != transTiaParams_.freqs.size() )
827  {
828  Report::UserError() << "The size of numFreq does not match the number of tones in .hb!";
829  }
830 
831 
832  if (numFreqs_.size() == 0 )
833  {
834 
835  numFreqs_.resize( transTiaParams_.freqs.size() ) ;
836  for (int i=0; i < transTiaParams_.freqs.size(); i++ )
837  {
838  numFreqs_[i] = size_;
839  }
840 
841  }
842 
843  return true;
844 }
845 
846 
847 //-----------------------------------------------------------------------------
848 // Function : HB::setLinSol
849 // Purpose : this is needed for .STEP to work with HB
850 // Special Notes :
851 // Scope : public
852 // Creator : Eric R. Keiter
853 // Creation Date : 7/12/2013
854 //-----------------------------------------------------------------------------
855 bool HB::setLinSol(const Util::OptionBlock & OB)
856 {
857  // Save the non-HB linear solver option block
858  saved_lsOB_ = OB;
859 
860  return true;
861 }
862 
863 //-----------------------------------------------------------------------------
864 // Function : HB::setHBLinSol
865 // Purpose :
866 // Special Notes :
867 // Scope : public
868 // Creator : Heidi Thornquist
869 // Creation Date : 5/13/13
870 //-----------------------------------------------------------------------------
871 bool HB::setHBLinSol(const Util::OptionBlock & OB)
872 {
873  // Save the HB linear solver option block
874  saved_lsHBOB_ = OB;
875 
876  // Generate the HB preconditioner factory.
877  precFactory_ = rcp( new N_LAS_HBPrecondFactory( OB ) );
878 
879  return true;
880 }
881 
882 //-----------------------------------------------------------------------------
883 // Function : HB::isAnalysis
884 // Purpose :
885 // Special Notes :
886 // Scope : public
887 // Creator : Heidi Thornquist
888 // Creation Date : 5/13/13
889 // Notes : Alternatively, we could try to cast the analysis object
890 // : However, this method is called a lot.
891 //-----------------------------------------------------------------------------
892 bool HB::isAnalysis( int analysis_type )
893 {
894  bool returnValue = false;
895 
896  if ( analysis_type == ANP_MODE_TRANSIENT )
897  {
898  returnValue = isTransient_;
899  }
900  if ( analysis_type == ANP_MODE_DC_SWEEP )
901  {
902  returnValue = isDCSweep_;
903  }
904  return returnValue;
905 }
906 
907 //-----------------------------------------------------------------------------
908 // Function : HB::prepareHBOutput
909 // Purpose :
910 // Special Notes :
911 // Scope : public
912 // Creator : Eric Keiter, 9233, Computational Sciences
913 // Creation Date : 08/20/07
914 //-----------------------------------------------------------------------------
916  N_LAS_Vector & solnVecPtr,
917  std::vector<double> & timePoints,
918  std::vector<double> & freqPoints,
919  RCP<N_LAS_BlockVector> & timeDomainSolnVec,
920  RCP<N_LAS_BlockVector> & freqDomainSolnVecReal,
921  RCP<N_LAS_BlockVector> & freqDomainSolnVecImag,
922  RCP<N_LAS_BlockVector> & timeDomainStoreVec,
923  RCP<N_LAS_BlockVector> & freqDomainStoreVecReal,
924  RCP<N_LAS_BlockVector> & freqDomainStoreVecImag
925  ) const
926 {
927  N_LAS_BlockVector & blockSolVecPtr = dynamic_cast<N_LAS_BlockVector &>(solnVecPtr);
928 
929  Teuchos::RCP<N_LAS_BlockVector> bStoreVecFreqPtr_ = hbLoaderPtr_->getStoreVecFreqPtr();
930 
931  timeDomainStoreVec = hbBuilderPtr_->createTimeDomainStoreBlockVector();
932 
933  if (bStoreVecFreqPtr_->blockCount() > 0 )
934  {
935  hbLoaderPtr_->permutedIFT(*bStoreVecFreqPtr_, &*timeDomainStoreVec);
936  }
937 
938  //TD solution
939  timeDomainSolnVec = hbBuilderPtr_->createTimeDomainBlockVector();
940  int blockCount = timeDomainSolnVec->blockCount();
941  int N = timeDomainSolnVec->block(0).globalLength();
942 
943  timePoints.resize(size_);
944 
945  for( int i = 0; i < size_; ++i )
946  {
947  timePoints[i] = fastTimes_[i] - transTiaParams_.initialTime;
948  }
949 
950  freqPoints = freqPoints_;
951 
952  // Create individual block vectors to store the real and imaginary parts separately.
953  Teuchos::RCP<N_PDS_ParMap> baseMap = Teuchos::rcp_const_cast<N_PDS_ParMap>( hbBuilderPtr_->getBaseSolutionMap() );
954  Teuchos::RCP<N_PDS_ParMap> globalMap = createBlockParMap( blockCount, *baseMap );
955  freqDomainSolnVecReal = Teuchos::rcp( new N_LAS_BlockVector( blockCount, globalMap, baseMap ) );
956  freqDomainSolnVecImag = Teuchos::rcp( new N_LAS_BlockVector( blockCount, globalMap, baseMap ) );
957 
958  hbLoaderPtr_->permutedIFT(blockSolVecPtr, &*timeDomainSolnVec);
959 
960  // Now copy over the frequency domain solution, real and imaginary parts separately, into the output vectors.
961  for (int j=0; j<N; j++)
962  {
963  // See if this time-domain solution variable is owned by the local processor.
964  // If so, this processor owns the entire j-th block of the blockSolVecPtr vector,
965  // and the j-th entry of every block in the freqDomainSolnVec[Real/Imag] vector.
966  int lid = baseMap->globalToLocalIndex( j );
967  N_LAS_Vector& solBlock = blockSolVecPtr.block( j );
968 
969  N_LAS_Vector& realVecRef = freqDomainSolnVecReal->block((blockCount-1)/2);
970  N_LAS_Vector& imagVecRef = freqDomainSolnVecImag->block((blockCount-1)/2);
971 
972  if (lid >= 0)
973  {
974  realVecRef[lid] = solBlock[0];
975  imagVecRef[lid] = solBlock[1];
976  }
977 
978  for (int i=1; i <= (blockCount-1)/2; ++i)
979  {
980  N_LAS_Vector& realVecRef_neg = freqDomainSolnVecReal->block((blockCount-1)/2 - i);
981  N_LAS_Vector& imagVecRef_neg = freqDomainSolnVecImag->block((blockCount-1)/2 - i);
982  N_LAS_Vector& realVecRef_pos = freqDomainSolnVecReal->block((blockCount-1)/2 + i);
983  N_LAS_Vector& imagVecRef_pos = freqDomainSolnVecImag->block((blockCount-1)/2 + i);
984 
985  if (lid >= 0)
986  {
987  realVecRef_neg[lid] = solBlock[ 2*(blockCount-i) ];
988  imagVecRef_neg[lid] = solBlock[ 2*(blockCount-i) + 1 ];
989  realVecRef_pos[lid] = solBlock[ 2*i ];
990  imagVecRef_pos[lid] = solBlock[ 2*i+1 ];
991  }
992  }
993  }
994 
995  // proceed to store variables
996  Teuchos::RCP<N_PDS_ParMap> baseStoreMap = Teuchos::rcp_const_cast<N_PDS_ParMap>( hbBuilderPtr_->getBaseStoreMap() );
997  Teuchos::RCP<N_PDS_ParMap> globalStoreMap = createBlockParMap( blockCount, *baseStoreMap );
998  freqDomainStoreVecReal = Teuchos::rcp( new N_LAS_BlockVector( blockCount, globalStoreMap, baseStoreMap ) );
999  freqDomainStoreVecImag = Teuchos::rcp( new N_LAS_BlockVector( blockCount, globalStoreMap, baseStoreMap ) );
1000 
1001  N = timeDomainStoreVec->block(0).globalLength();
1002 
1003  for (int j=0; j<N; j++)
1004  {
1005  // See if this time-domain solution variable is owned by the local processor.
1006  // If so, this processor owns the entire j-th block of the blockSolVecPtr vector,
1007  // and the j-th entry of every block in the freqDomainSolnVec[Real/Imag] vector.
1008  int lid = baseStoreMap->globalToLocalIndex( j );
1009  N_LAS_Vector& storeBlock = bStoreVecFreqPtr_->block( j );
1010 
1011  N_LAS_Vector& realVecRef = freqDomainStoreVecReal->block((blockCount-1)/2);
1012  N_LAS_Vector& imagVecRef = freqDomainStoreVecImag->block((blockCount-1)/2);
1013 
1014  if (lid >= 0)
1015  {
1016  realVecRef[lid] = storeBlock[0];
1017  imagVecRef[lid] = storeBlock[1];
1018  }
1019 
1020  for (int i=1; i <= (blockCount-1)/2; ++i)
1021  {
1022  N_LAS_Vector& realVecRef_neg = freqDomainStoreVecReal->block((blockCount-1)/2 - i);
1023  N_LAS_Vector& imagVecRef_neg = freqDomainStoreVecImag->block((blockCount-1)/2 - i);
1024  N_LAS_Vector& realVecRef_pos = freqDomainStoreVecReal->block((blockCount-1)/2 + i);
1025  N_LAS_Vector& imagVecRef_pos = freqDomainStoreVecImag->block((blockCount-1)/2 + i);
1026 
1027  if (lid >= 0)
1028  {
1029  realVecRef_neg[lid] = storeBlock[ 2*(blockCount-i) ];
1030  imagVecRef_neg[lid] = storeBlock[ 2*(blockCount-i) + 1 ];
1031  realVecRef_pos[lid] = storeBlock[ 2*i ];
1032  imagVecRef_pos[lid] = storeBlock[ 2*i+1 ];
1033  }
1034  }
1035  }
1036 
1037 #ifdef Xyce_DEBUG_HB
1038 // Xyce::dout() << "HB X Vector FD" << std::endl;
1039  freqDomainSolnVecReal->printPetraObject(std::cout);
1040  freqDomainSolnVecImag->printPetraObject(std::cout);
1041 
1042  Xyce::dout() << "HB Store Vector FD" << std::endl;
1043 
1044  freqDomainStoreVecReal->printPetraObject(std::cout);
1045  freqDomainStoreVecImag->printPetraObject(std::cout);
1046 #endif // Xyce_DEBUG_HB
1047 
1048 
1049 }
1050 
1051 
1052 //-----------------------------------------------------------------------------
1053 // Function : HB::accumulateStatistics()
1054 // Purpose : Add in the statistics from the current analysis object
1055 // Special Notes :
1056 // Scope : private
1057 // Creator : Heidi Thornquist, 1355, Electrical Models & Simulation
1058 // Creation Date : 05/29/13
1059 //-----------------------------------------------------------------------------
1061 {
1062  hbStatCounts_ += analysisObject_->stats_;
1063 }
1064 
1065 
1066 //-----------------------------------------------------------------------------
1067 // Function : HB::runTol_
1068 // Purpose : Conducts transient run to determine right tolerance
1069 // parameters for IC calculation
1070 // Special Notes :
1071 // Scope : private
1072 // Creator : T. Mei, 1437, Electrical and Micro Modeling
1073 // Creation Date : 02/23/09
1074 //-----------------------------------------------------------------------------
1076 {
1077  bool returnValue = true;
1078 
1079  Xyce::lout() << " ***** Computing tolerance parameters for HB IC calculation....\n" << std::endl;
1080 
1082  N_TIA_TIAParams tiaParamsSave = tiaParams;
1083 
1084  // now try the real thing
1085  tiaParams.initialTime = 0;
1086  tiaParams.finalTime = period_;
1087  tiaParams.pauseTime = tiaParams.finalTime;
1088  tiaParams.resume = false;
1089  tiaParams.maxOrder = 1;
1090 
1091  // If start up periods are not being run then we can use this transient to compute HB ICs.
1092  if (!startUpPeriodsGiven_)
1093  {
1094  tiaParams.saveTimeStepsFlag = true;
1095  }
1096 
1097  // register new tiaParams with time integrator
1098  analysisManager_.registerTIAParams (tiaParams);
1099 
1100  // Create a transient analysis object for this section.
1101  isTransient_ = true;
1102  analysisObject_ = Teuchos::rcp(new Transient(analysisManager_));
1103  analysisObject_->setAnalysisParams( Util::OptionBlock() );
1104  Teuchos::rcp_dynamic_cast<Transient>(analysisObject_)->resetForHB();
1105  returnValue = analysisObject_->run();
1106 
1107  if (!returnValue)
1108  {
1109  Report::UserError() << "Calculation of tolerance parameters failed for relErrorTol = " << tiaParams.relErrorTol;
1110  return false;
1111  }
1112 
1113  int numPoints = analysisManager_.getStepNumber();
1114 
1115  while((numPoints < (1.2*size_)) && (tiaParams.relErrorTol>= 1e-6))
1116  {
1117  {
1118  Report::UserWarning() << "Tolerance parameters refined, re-running with relErrorTol = " << tiaParams.relErrorTol/10;
1119  }
1120 
1121  if (!startUpPeriodsGiven_)
1122  {
1123  // Clear the fast time data storage before performing the next transient
1125  dsPtr->resetFastTimeData();
1126  }
1127 
1128  tiaParams.relErrorTol = tiaParams.relErrorTol/10;
1129  analysisManager_.registerTIAParams (tiaParams);
1130 
1131  // Create a transient analysis object for this section.
1132  analysisObject_ = Teuchos::rcp(new Transient(analysisManager_));
1133  analysisObject_->setAnalysisParams( Util::OptionBlock() );
1134  Teuchos::rcp_dynamic_cast<Transient>(analysisObject_)->resetForHB();
1135  bool retV = analysisObject_->run();
1136 
1137  if (!retV)
1138  {
1139  Report::UserError() << "Calculation of tolerance parameters failed for relErrorTol = " << tiaParams.relErrorTol;
1140  return false;
1141  }
1142  returnValue = retV && returnValue;
1143 
1144  numPoints = analysisManager_.getStepNumber();
1145  }
1146 
1147  // Add in simulation times
1149 
1150  // Reset parameters
1151  tiaParamsSave.relErrorTol = tiaParams.relErrorTol;
1152  analysisManager_.registerTIAParams (tiaParamsSave);
1153 
1154  // Reset transient flag before exiting.
1155  isTransient_ = false;
1156 
1157  return returnValue;
1158 }
1159 
1160 //-----------------------------------------------------------------------------
1161 // Function : HB::runStartupPeriods_()
1162 // Purpose : Runs normal transient problem through the requested
1163 // number of startup periods
1164 // Special Notes :
1165 // Scope : private
1166 // Creator : Richard Schiek, 1437, Electrical and Micro Modeling
1167 // Creation Date : 09/17/07
1168 //-----------------------------------------------------------------------------
1170 {
1171  bool returnValue = true;
1172 
1173  Xyce::lout() << " ***** Computing " << startUpPeriods_ << " start up periods for HB IC calculation...." << std::endl;
1174 
1175  // need to advance time by startUpPeriods_ * period_
1177  N_TIA_TIAParams tiaParamsSave = tiaParams;
1178 
1179  // set DAE initial time = 0.0
1180  tiaParams.initialTime = 0.0;
1181  tiaParams.finalTime = startUpPeriods_ * period_;
1182  tiaParams.pauseTime = tiaParams.finalTime;
1183 
1184 #ifdef Xyce_DEBUG_HB
1185  Xyce::dout() << "HB::runStartupPeriods_(): Advancing time through "
1186  << startUpPeriods_ << " startup periods"
1187  << " initialTime = " << tiaParams.initialTime
1188  << " finalTime = " << tiaParams.finalTime << std::endl;
1189 
1190  Xyce::dout() << "HB::runStartupPeriods_(): Double DCOP tiaParams:"
1191  << " doubleDCOPStep = " << tiaParams.doubleDCOPStep
1192  << " firstDCOPStep = " << tiaParams.firstDCOPStep
1193  << " lastDCOPStep = " << tiaParams.lastDCOPStep << std::endl;
1194 
1195  Xyce::dout() << "HB::runStartupPeriods_(): Double DCOP tiaParamsSave:"
1196  << " doubleDCOPStep = " << tiaParamsSave.doubleDCOPStep
1197  << " firstDCOPStep = " << tiaParamsSave.firstDCOPStep
1198  << " lastDCOPStep = " << tiaParamsSave.lastDCOPStep << std::endl;
1199 
1200 #endif // Xyce_DEBUG_HB
1201 
1202  {
1203  Xyce::IO::ActiveOutput x(*analysisManager_.getOutputManager());
1204  x.add(Xyce::IO::PrintType::HB_STARTUP, ANP_MODE_HB);
1205 
1206  // register new tiaParams with time integrator
1207  analysisManager_.registerTIAParams (tiaParams);
1208 
1209  // Create a transient analysis object for this section.
1210  isTransient_ = true;
1211  analysisObject_ = Teuchos::rcp(new Transient(analysisManager_));
1212  analysisObject_->setAnalysisParams( Util::OptionBlock() );
1213  Teuchos::rcp_dynamic_cast<Transient>(analysisObject_)->resetForHB();
1214  returnValue = analysisObject_->run();
1215  isTransient_ = false;
1216 
1217  // Add in simulation times
1219 
1221  }
1222 
1223  // reset the output filename suffix
1224  // analysisManager_.outMgrPtr->setOutputFilenameSuffix("");
1225 
1226  // put the dsPtr->currentSolutionPtr into dcOpSol and State Vec so that it
1227  // is used as our initial condition for the pending fast time scale runs
1229  dcOpSolVecPtr_ = rcp( new N_LAS_Vector( *(dsPtr->currSolutionPtr) ));
1230  dcOpStateVecPtr_ = rcp( new N_LAS_Vector( *(dsPtr->currStatePtr) ));
1231  dcOpQVecPtr_ = rcp( new N_LAS_Vector( *(dsPtr->daeQVectorPtr) ));
1232  dcOpStoreVecPtr_ = rcp( new N_LAS_Vector( *(dsPtr->currStorePtr) ));
1233 
1234  // tell HB to start after this startup period
1235  tiaParamsSave.initialTime = startUpPeriods_ * period_;
1237  analysisManager_.registerTIAParams (tiaParamsSave);
1238  startUpPeriodsFinished_ = true;
1239 
1240  return returnValue;
1241 }
1242 
1243 
1244 // bool setfreqPoints_();
1245 
1246 //-----------------------------------------------------------------------------
1247 // Function : HB::setFreqPoints_
1248 // Purpose : Set frequency spectrum for HB analysis.
1249 // Special Notes :
1250 // Scope : private
1251 // Creator : Ting Mei, SNL
1252 // Creation Date : 03/03/2014
1253 //-----------------------------------------------------------------------------
1255 {
1256  if ( !intmodMaxGiven_)
1257  {
1258 
1259  int maxValue = 0;
1260  if (numFreqs_.size() != 0 )
1261  // find the max of numFreqs
1262  {
1263  maxValue = (numFreqs_[0] - 1)/2;
1264  for (int i=1; i<numFreqs_.size(); ++i)
1265  {
1266  if ((numFreqs_[i] - 1)/2 > maxValue)
1267  maxValue = (numFreqs_[i] - 1)/2 ;
1268 
1269  }
1270  }
1271  else
1272  {
1273  maxValue = (size_ - 1)/2;
1274  }
1275 
1276  intmodMax_ = maxValue;
1277 
1278  }
1279 
1280 
1281 // std::vector<int> numPosFreqs;
1282 
1283  std::vector<int> k;
1284 
1285  int numAnalysisFreqs = transTiaParams_.freqs.size();
1286 
1287  numPosFreqs.resize(numAnalysisFreqs);
1288 
1289  k.resize(numAnalysisFreqs);
1290 
1291  k[0] = 1;
1292 
1293  if (numFreqs_.size() != 0 )
1294  numPosFreqs[0] = (numFreqs_[0] - 1)/2;
1295  else
1296  numPosFreqs[0]= (size_ - 1)/2;
1297 
1298  int numTotalFrequencies;
1299 
1300  int numExtraFreqs = 0;
1301 
1302  if (numPosFreqs[0] > intmodMax_)
1303  {
1304  numExtraFreqs += ( numPosFreqs[0] - intmodMax_ );
1305  numFreqs_[0] = (intmodMax_*2 + 1);
1306  }
1307 
1308  numTotalFrequencies = numFreqs_[0];
1309 
1310  for (int i=1; i < numAnalysisFreqs; i++)
1311  {
1312  numPosFreqs[i] = (numFreqs_[i] - 1)/2;
1313 
1314  if (numPosFreqs[i] > intmodMax_)
1315  {
1316  numExtraFreqs += ( numPosFreqs[i] - intmodMax_ );
1317  numFreqs_[i] = (intmodMax_*2 + 1);
1318  }
1319 
1320  k[i] = k[i-1] * numFreqs_[i -1 ];
1321 
1322  numTotalFrequencies *= numFreqs_[i];
1323  }
1324 
1325 #ifdef Xyce_DEBUG_HB
1326  for (int i=0; i< numAnalysisFreqs; i++)
1327  {
1328  Xyce::dout() << "HB index i" << i << std::endl;
1329  Xyce::dout() << "HB numPosFreqs =" << numPosFreqs[i] << std::endl;
1330  Xyce::dout() << "HB k =" << k[i] << std::endl;
1331  }
1332  Xyce::dout() << "HB numTotalFrequencies =" << numTotalFrequencies<< std::endl;
1333 
1334  Xyce::dout() << "HB numextrafreqs =" << numExtraFreqs << std::endl;
1335 
1336 #endif
1337 
1338  int numIndex = numTotalFrequencies;
1339 
1340  Teuchos::SerialDenseMatrix<int,double> indexMatrix(numAnalysisFreqs, numTotalFrequencies);
1341 // Teuchos::SerialDenseMatrix<int,int> indexMatrix(numAnalysisFreqs, numTotalFrequencies);
1342 
1343 
1344 #ifdef Xyce_DEBUG_HB
1345  Xyce::dout() << "HB intmodMax =" << intmodMax_ << std::endl;
1346 #endif
1347 
1348  int nextIndex;
1349 
1350  int idxMod, idxValues;
1351  int sumIndex;
1352 
1353  std::vector<int> goodIndex;
1354  for (int i=0; i < numIndex; i++) // column
1355  {
1356  nextIndex = i;
1357  sumIndex = 0;
1358 
1359  for (int j= (numAnalysisFreqs - 1); j >= 0; j-- ) // row
1360  {
1361  idxMod = nextIndex%k[j];
1362  idxValues = (nextIndex - idxMod)/k[j];
1363 
1364  indexMatrix (j, i) = static_cast<double>(idxValues - (numFreqs_[j] - 1)/2 );
1365  nextIndex = idxMod;
1366  sumIndex += abs(idxValues - (numFreqs_[j] - 1)/2 );
1367 
1368  }
1369 
1370  if( sumIndex <= intmodMax_)
1371  goodIndex.push_back(i);
1372 
1373  }
1374 
1375  int diaindexSize = goodIndex.size();
1376 
1377  Teuchos::SerialDenseMatrix<int,double> diaindexMatrix( numAnalysisFreqs, (diaindexSize + numExtraFreqs) );
1378  diaindexMatrix.putScalar(0.0);
1379 
1380  for (int i=0; i < diaindexSize; i++)
1381  {
1382  for (int j= (numAnalysisFreqs - 1); j >= 0; j-- )
1383  diaindexMatrix (j, i) = indexMatrix (j, goodIndex[i]);
1384  }
1385 
1386 #ifdef Xyce_DEBUG_HB
1387  for (int i=0; i< diaindexSize; i++)
1388  {
1389  std::cout << "good index i = " << i << goodIndex[i] << std::endl;
1390  }
1391 
1392 
1393  std::cout << " checking diamond indexMatrix" << std::endl;
1394  diaindexMatrix.print(std::cout);
1395 
1396  std::cout << " checking indexMatrix" << std::endl;
1397  indexMatrix.print(std::cout);
1398 
1399 #endif
1400 
1401  int extraIndexPos = diaindexSize;
1402  for (int i=0; i < numAnalysisFreqs ; i++) // column
1403  {
1404 
1405  if (numPosFreqs[i] > intmodMax_)
1406  {
1407 
1408  for (int j=0; j < (numPosFreqs[i] - intmodMax_ ) ; j++)
1409  diaindexMatrix (i, extraIndexPos + j ) = static_cast<double>(intmodMax_ + j + 1 );
1410 
1411  extraIndexPos += (numPosFreqs[i] - intmodMax_ );
1412 
1413  }
1414 
1415  }
1416 
1417 
1418 #ifdef Xyce_DEBUG_HB
1419  std::cout << " checking diamond indexMatrix after axis" << std::endl;
1420  diaindexMatrix.print(std::cout);
1421 #endif
1422 
1423 // get the positive frequencies
1424 
1425  std::vector<double> posfreqPoints_;
1426 
1427  int posindexSize = (diaindexSize - 1)/2;
1428  posfreqPoints_.resize(posindexSize + numExtraFreqs);
1429 
1430  Teuchos::SerialDenseMatrix<int,double> currindexMatrix( Teuchos::View, diaindexMatrix, numAnalysisFreqs, (posindexSize + numExtraFreqs), 0, posindexSize+1 );
1431  Teuchos::SerialDenseVector<int,double> currfreqPoints( Teuchos::View, &posfreqPoints_[0], (posindexSize + numExtraFreqs ) );
1432 
1433  Teuchos::SerialDenseVector<int,double> hbFreqs( Teuchos::View, &transTiaParams_.freqs[0], numAnalysisFreqs);
1434 // Teuchos::SerialDenseVector<int,double> currWeightVector( Teuchos::View, &weightVector[i], oversampleRate*size_-(i+1) );
1435  currfreqPoints.multiply( Teuchos::TRANS, Teuchos::NO_TRANS, 1.0, currindexMatrix, hbFreqs, 0.0 );
1436 
1437 
1438 #ifdef Xyce_DEBUG_HB
1439  std::cout << "checking positive frequencies" << std::endl;
1440  currindexMatrix.print(std::cout);
1441  hbFreqs.print(std::cout);
1442  currfreqPoints.print(std::cout);
1443 #endif
1444 
1445  for (int i=0; i < posindexSize; i++)
1446  {
1447 
1448  if (posfreqPoints_[i] < 0.0)
1449  posfreqPoints_[i] = fabs( posfreqPoints_[i]);
1450 
1451  }
1452 
1453 // size_ = (posindexSize + numExtraFreqs ) *2 + 1;
1454 
1455  std::sort(posfreqPoints_.begin(), posfreqPoints_.end() );
1456 
1457 
1458 #ifdef Xyce_DEBUG_HB
1459  for (int i=0; i< posfreqPoints_.size(); i++)
1460  {
1461  std::cout << "pos frequency point " << posfreqPoints_[i] << std::endl;
1462  }
1463 #endif
1464  posfreqPoints_.erase(std::unique(posfreqPoints_.begin(), posfreqPoints_.end() ), posfreqPoints_.end() );
1465 
1466 
1467  if (abs( posfreqPoints_[0]) < 2.0*Util::MachineDependentParams::MachinePrecision() )
1468  posfreqPoints_.erase( posfreqPoints_.begin());
1469 
1470  size_ = ( posfreqPoints_.size() ) *2 + 1;
1471 #ifdef Xyce_DEBUG_HB
1472  for (int i=0; i<posfreqPoints_.size(); i++)
1473  {
1474  std::cout << "pos frequency point after " << posfreqPoints_[i] << std::endl;
1475  }
1476 
1477  Xyce::dout() << "HB size =" << size_ << std::endl;
1478 #endif
1479 
1480  int i=0;
1481  freqPoints_.resize(size_);
1482 
1483  for( i = 0; i < size_; ++i )
1484  {
1485  if (i < (size_-1)/2)
1486  freqPoints_[i] = - posfreqPoints_[ (size_-1)/2 - i - 1 ];
1487  else if (i > (size_-1)/2)
1488  freqPoints_[i] = posfreqPoints_[ i - (size_-1)/2 - 1 ];
1489  else
1490  freqPoints_[i] = 0.0;
1491  }
1492 
1493 #ifdef Xyce_DEBUG_HB
1494  for (int i=0; i<freqPoints_.size(); i++)
1495  {
1496  std::cout << " frequency point " << freqPoints_[i] << std::endl;
1497  }
1498 #endif
1499 
1500 
1501  return true;
1502 }
1503 
1504 //-----------------------------------------------------------------------------
1505 // Function : HB::setTimePoints_
1506 // Purpose : Set time points for multi-tone HB analysis.
1507 // Special Notes :
1508 // Scope : private
1509 // Creator : Heidi Thornquist and Ting Mei, SNL
1510 // Creation Date : 03/05/2014
1511 //-----------------------------------------------------------------------------
1513 {
1514  // NOTE: Need to make this parallel safe.
1515 
1516  int posFreq = (size_-1)/2;
1517  int oversampleRate = 2;
1518  int periodSampleMultiplier = 1;
1519 // int periodSampleMultiplier = 1;
1520 
1521  Teuchos::BLAS<int,double> blas;
1522  std::vector<double> testPoints(oversampleRate*size_);
1523 
1524 // std::cout << "Period is" << period_ << std::endl;
1525  for (int i=0; i<oversampleRate*size_; ++i)
1526  {
1527  testPoints[i] = periodSampleMultiplier*period_*((Teuchos::ScalarTraits<double>::random()+1)/2);
1528  // testPoints[i] = periodSampleMultiplier*period_/(oversampleRate*size_)* static_cast<double>(i);
1529 // std::cout << "testPoints" << i << "=" << testPoints[i] << std::endl;
1530  }
1531 
1532  Teuchos::SerialDenseMatrix<int,double> testMatrix(size_,oversampleRate*size_);
1533  // NOTE: i represents frequency, j represents sample time
1534  // Set DC values first.
1535  for (int j=0; j<oversampleRate*size_; ++j)
1536  {
1537  testMatrix(0,j) = 1.0;
1538  }
1539  // Set rest of frequency values
1540  for (int i=1; i<=posFreq; i++)
1541  {
1542 // for (int j=0; j<oversampleRate*size_; j+2)
1543  for (int j=0; j<oversampleRate*size_; j++)
1544  {
1545  testMatrix(2*i-1,j) = cos(2*M_PI*freqPoints_[posFreq+i]*testPoints[j]);
1546  testMatrix(2*i,j) = sin(2*M_PI*freqPoints_[posFreq+i]*testPoints[j]);
1547  }
1548  }
1549  // Check matrix is correct
1550 // std::cout << "Checking testMatrix" << std::endl;
1551 // testMatrix.print(std::cout);
1552 
1553  // Now for orthogonalization (yipeee)
1554  std::vector<double> weightVector(oversampleRate*size_);
1555  for (int i=0; i<size_; ++i)
1556  {
1557  // Find column with largest norm, choose as next vector for orthogonalization
1558  int maxIndex = 0;
1559  double maxValue = 0.0;
1560  for (int j=i; j<oversampleRate*size_; ++j)
1561  {
1562  Teuchos::SerialDenseMatrix<int,double> tempVector( Teuchos::View, testMatrix, size_, 1, 0, j );
1563  weightVector[j] = tempVector.normFrobenius();
1564 
1565  if (weightVector[j] > maxValue)
1566  {
1567  maxIndex = j;
1568  maxValue = weightVector[j];
1569  }
1570  }
1571 
1572 
1573  // Swap time and vector.
1574  std::swap( testPoints[i], testPoints[maxIndex] );
1575  Teuchos::SerialDenseVector<int,double> newSwapVector2 = Teuchos::getCol<int,double>( Teuchos::Copy, testMatrix, maxIndex );
1576  Teuchos::SerialDenseVector<int,double> newSwapVector = Teuchos::getCol<int,double>( Teuchos::View, testMatrix, i );
1577  Teuchos::setCol<int,double>( newSwapVector, maxIndex, testMatrix );
1578  Teuchos::setCol<int,double>( newSwapVector2, i, testMatrix );
1579 
1580 // std::cout << "Checking testMatrix for i= " << i << std::endl;
1581 // testMatrix.print(std::cout);
1582 
1583  // Compute inner product with vector from last time point with rest of time points
1584  Teuchos::SerialDenseMatrix<int,double> currTestMatrix( Teuchos::View, testMatrix, size_, oversampleRate*size_-(i+1), 0, i+1 );
1585  Teuchos::SerialDenseVector<int,double> currWeightVector( Teuchos::View, &weightVector[i+1], oversampleRate*size_-(i+1) );
1586 // Teuchos::SerialDenseVector<int,double> currWeightVector( Teuchos::View, &weightVector[i], oversampleRate*size_-(i+1) );
1587  currWeightVector.multiply( Teuchos::TRANS, Teuchos::NO_TRANS, 1.0, currTestMatrix, newSwapVector, 0.0 );
1588 
1589 // std::cout << "The current norm is" << std::endl;
1590 // currWeightVector.print(std::cout);
1591 
1592  // Subtract off scaled vector from rest of time points.
1593  // NOTE: maxValue is the norm of the last time point.
1594  for (int j=i+1; j<oversampleRate*size_; ++j)
1595  {
1596  Teuchos::SerialDenseMatrix<int,double> currVector( Teuchos::View, testMatrix, size_, 1, 0, j );
1597  blas.AXPY( size_, -(currWeightVector[j-(i+1)]/(maxValue*maxValue) ), newSwapVector.values(), 1, currVector.values(), 1 );
1598  }
1599 
1600 // std::cout << "Checking testMatrix after orthogonalization for i = " << i << std::endl;
1601 // testMatrix.print(std::cout);
1602 
1603  }
1604 // std::cout << "Checking testMatrix after orthogonalization." << std::endl;
1605 // testMatrix.print(std::cout);
1606 
1607  // Sort the chosen test points and then copy them into the fastTimes_ vector.
1608  std::sort( testPoints.begin(), testPoints.begin()+size_ );
1609 
1610  fastTimes_.resize (size_);
1611 
1612  for (int i=0; i<size_; ++i)
1613  {
1614  fastTimes_[i] = testPoints[i];
1615  }
1616 
1617  return true;
1618 }
1619 
1620 //-----------------------------------------------------------------------------
1621 // Function : HB::createFT_
1622 // Purpose : Create the DFT and IFT matrices using the time points for multi-tone HB analysis.
1623 // Special Notes :
1624 // Scope : private
1625 // Creator : Heidi Thornquist and Ting Mei, SNL
1626 // Creation Date : 03/05/2014
1627 //-----------------------------------------------------------------------------
1629 {
1630  int posFreq = (size_-1)/2;
1631  idftMatrix_.reshape(size_,size_);
1632 
1633  // NOTE: i represents frequency, j represents sample time
1634  // Set DC values first.
1635  for (int i=0; i<size_; ++i)
1636  {
1637  idftMatrix_(i,0) = 1.0;
1638  }
1639  // Set rest of frequency values
1640  for (int i=0; i<size_; i++)
1641  {
1642  for (int j=1; j<=posFreq; j++)
1643  {
1644  idftMatrix_(i,2*j-1) = cos(2*M_PI*freqPoints_[posFreq+j]*fastTimes_[i]);
1645  idftMatrix_(i,2*j) = sin(2*M_PI*freqPoints_[posFreq+j]*fastTimes_[i]);
1646  }
1647  }
1648 
1649  std::cout << "Checking IDFTmatrix" << std::endl;
1650  idftMatrix_.print(std::cout);
1651 
1652  // Compute DFT matrix.
1654  Teuchos::SerialDenseSolver<int,double> ftSolver;
1655  ftSolver.setMatrix( Teuchos::rcp( &dftMatrix_, false ) );
1656  ftSolver.invert();
1657 
1658  std::cout << "Checking DFTmatrix" << std::endl;
1659  dftMatrix_.print(std::cout);
1660 
1661  std::cout << "Checking IDFTmatrix after inverse:" << std::endl;
1662  idftMatrix_.print(std::cout);
1663 
1664  return true;
1665 
1666 }
1667 
1668 //-----------------------------------------------------------------------------
1669 // Function : HB::setInitialGuess_
1670 // Purpose : Set initial guess for HB analysis.
1671 // Special Notes :
1672 // Scope : private
1673 // Creator : Ting Mei, SNL
1674 // Creation Date : 03/03/2014
1675 //-----------------------------------------------------------------------------
1677 {
1678 
1679  bool returnValue = true;
1680 
1681  if (taHB_ == 1)
1682  {
1683  bool retTol1 = runTol_(); returnValue = returnValue && retTol1;
1684 
1685  // Start up periods need to be run before the initial condition is computed, otherwise
1686  // just used the solution from the tolerance calculation.
1687  if( startUpPeriodsGiven_ )
1688  {
1689  bool startupPeriodsSuccess = runStartupPeriods_();
1690  if (!startupPeriodsSuccess)
1691  {
1692  Report::UserError() << "Failed to calculate the startup periods";
1693  return false;
1694  }
1695  returnValue = returnValue && startupPeriodsSuccess;
1696 
1697  bool icSuccess = runTransientIC_();
1698  if (!icSuccess)
1699  {
1700  Report::UserError() << "Initial HB Transient failed";
1701  return false;
1702  }
1703  returnValue = returnValue && icSuccess;
1704  }
1705 
1706  interpolateIC_();
1707 
1708  }
1709  else
1710  {
1711 
1712  double TimeStep = period_/static_cast<double>(size_);
1713  timeSteps_.push_back( TimeStep );
1714  fastTimes_.resize(size_+1);
1715 
1716  goodTimePoints_.resize(size_+1);
1717  for( int i = 0; i <= size_; ++i )
1718  {
1719  fastTimes_[i] = transTiaParams_.initialTime + static_cast<double>(i) * TimeStep;
1720  }
1722 
1723  }
1724 
1725  return returnValue;
1726 }
1727 
1728 //-----------------------------------------------------------------------------
1729 // Function : HB::runTransientIC_
1730 // Purpose : Conducts a regular transient run for HB initial conditions
1731 // Special Notes :
1732 // Scope : private
1733 // Creator : Ting Mei, SNL
1734 // Creation Date : 10/03/2008
1735 //-----------------------------------------------------------------------------
1737 {
1738  bool returnValue = true;
1739 
1740  Xyce::lout() << " ***** Running transient to compute HB initial condition....\n" << std::endl;
1741 
1742  // this prevents extra DC op data from being printed.
1743  devInterfacePtr_->setMPDEFlag( true );
1744 
1745  if(saveIcData_)
1746  {
1747  // Keep the initial condition data
1748  // analysisManager_.outMgrPtr->setOutputFilenameSuffix( ".hb_ic" );
1749  }
1750 
1751  // use an initial transient run to create a set of time points for the fast time scale
1753  N_TIA_TIAParams tiaParamsSave = tiaParams;
1754 
1755  //tiaParams.initialTime = 0; // should be tiaParams_.initialTime;
1757  //tiaParams.finalTime = period_; // should be tiaParams_.initialTime + period_;
1759  tiaParams.pauseTime = tiaParams.finalTime;
1760  tiaParams.saveTimeStepsFlag = true;
1761  tiaParams.maxOrder = 1;
1762 
1763 #ifdef Xyce_DEBUG_HB
1764  Xyce::dout() << "HB::runTransientIC_(): Advancing time from"
1765  << " initialTime = " << tiaParams.initialTime
1766  << " finalTime = " << tiaParams.finalTime << std::endl;
1767 #endif // Xyce_DEBUG_HB
1768 
1769  // Initial conditions will be set if startup periods were run.
1770  if ( startUpPeriodsGiven_ )
1771  {
1772  tiaParams.NOOP = true;
1773 
1775  *(dsPtr->nextSolutionPtr) = *(dcOpSolVecPtr_.get());
1776  *(dsPtr->nextStatePtr) = *(dcOpStateVecPtr_.get());
1777  *(dsPtr->daeQVectorPtr) = *(dcOpQVecPtr_.get());
1778  *(dsPtr->nextStorePtr) = *(dcOpStoreVecPtr_.get());
1779  }
1780 
1781  // register new tiaParams with time integrator
1782  analysisManager_.registerTIAParams (tiaParams);
1783 
1784  // Create a transient analysis object for this section.
1785  isTransient_ = true;
1786  analysisObject_ = Teuchos::rcp(new Transient(analysisManager_));
1787  analysisObject_->setAnalysisParams( Util::OptionBlock() );
1788  Teuchos::rcp_dynamic_cast<Transient>(analysisObject_)->resetForHB();
1789  returnValue = analysisObject_->run();
1790  isTransient_ = false;
1791 
1792  // Add in simulation times
1794 
1795  if(saveIcData_)
1796  {
1797  // reset suffix
1798  // analysisManager_.outMgrPtr->setOutputFilenameSuffix( "" );
1799  }
1800 
1801  tiaParamsSave.initialTime += period_; // start HB problem after this transient init.
1802  analysisManager_.registerTIAParams (tiaParamsSave);
1803  devInterfacePtr_->setMPDEFlag( false );
1804 
1805  return returnValue;
1806 }
1807 
1808 //-----------------------------------------------------------------------------
1809 // Function : HB::interpolateIC_()
1810 // Purpose : Tries to filter the fast time points from a transient run
1811 // so that points are not too close together
1812 // Special Notes :
1813 // Scope : private
1814 // Creator : Richard Schiek, 1437, Electrical and Micro Modeling
1815 // Creation Date : 09/17/07
1816 //-----------------------------------------------------------------------------
1818 {
1819  Xyce::lout() << " ***** Interpolating transient solution for IC calculation....\n" << std::endl;
1820 
1822  int numPoints = dsPtr->timeSteps.size();
1823 
1824 #ifdef Xyce_DEBUG_HB
1825  Xyce::dout() << "HB::interpolateIC_(): Initial transient run produced " << numPoints << " points." << std::endl;
1826 #endif
1827 
1828  std::vector<int> goodIndicies;
1829  goodTimePoints_.resize(size_);
1830 
1831  double TimeStep = period_/static_cast<double>(size_);
1832  timeSteps_.push_back( TimeStep );
1833  for( int i = 0; i < size_; ++i )
1834  {
1835  goodTimePoints_[i] = transTiaParams_.initialTime + static_cast<double>(i) * TimeStep;
1836  }
1838  fastTimes_.resize(size_+1);
1840 
1841  int breakpoints = 0; // need to keep track of how many breakpoints there are
1842  int startIndex = 0;
1843 
1844  // always keep first point
1845  goodIndicies.push_back(startIndex);
1846  int GoodTimePointIndex = startIndex + 1;
1847 
1848  for( int i=startIndex; i < numPoints - 1 ; i++ )
1849  {
1850  // count up breakpoints
1851  if( dsPtr->timeStepsBreakpointFlag[i] == true )
1852  {
1853  breakpoints++;
1854  }
1855 
1856 #ifdef Xyce_DEBUG_HB
1857  if( debugLevel > 0 )
1858  {
1859  Xyce::dout() << "\t\t timeStep[ " << i << " ] = " << dsPtr->timeSteps[i];
1860  if( dsPtr->timeStepsBreakpointFlag[i] == true )
1861  {
1862  Xyce::dout() << " Breakpoint";
1863  }
1864  Xyce::dout() << std::endl;
1865  }
1866 #endif
1867  while( ( GoodTimePointIndex < size_ ) && (dsPtr->timeSteps[i] <= goodTimePoints_[GoodTimePointIndex]) && (goodTimePoints_[GoodTimePointIndex] < dsPtr->timeSteps[i+1]))
1868  {
1869  // found a good point so save the index
1870  goodIndicies.push_back( i );
1871  GoodTimePointIndex = GoodTimePointIndex+1;
1872  }
1873  }
1874 
1875  for(int i=0; i<size_; i++ )
1876  {
1877  int currentIndex = goodIndicies[i];
1878  N_LAS_Vector * firstSolVecPtr = dsPtr->fastTimeSolutionVec[currentIndex];
1879  N_LAS_Vector * secondSolVecPtr = dsPtr->fastTimeSolutionVec[currentIndex+1];
1880 
1881  N_LAS_Vector * firstStateVecPtr = dsPtr->fastTimeStateVec[currentIndex];
1882  N_LAS_Vector * secondStateVecPtr = dsPtr->fastTimeStateVec[currentIndex+1];
1883 
1884  N_LAS_Vector * firstQVecPtr = dsPtr->fastTimeQVec[currentIndex];
1885  N_LAS_Vector * secondQVecPtr = dsPtr->fastTimeQVec[currentIndex+1];
1886 
1887  N_LAS_Vector * firstStoreVecPtr = dsPtr->fastTimeStoreVec[currentIndex];
1888  N_LAS_Vector * secondStoreVecPtr = dsPtr->fastTimeStoreVec[currentIndex+1];
1889 
1890  double fraction = (goodTimePoints_[i] - dsPtr->timeSteps[currentIndex])/(dsPtr->timeSteps[currentIndex+1] - dsPtr->timeSteps[currentIndex]);
1891 
1892  RCP<N_LAS_Vector> InterpICSolVecPtr = rcp( new N_LAS_Vector( *secondSolVecPtr ) );
1893  RCP<N_LAS_Vector> InterpICStateVecPtr = rcp( new N_LAS_Vector( *secondStateVecPtr ) );
1894  RCP<N_LAS_Vector> InterpICQVecPtr = rcp( new N_LAS_Vector( *secondQVecPtr ) );
1895  RCP<N_LAS_Vector> InterpICStoreVecPtr = rcp( new N_LAS_Vector( *secondStoreVecPtr ) );
1896 
1897  InterpICSolVecPtr->putScalar(0.0);
1898  InterpICStateVecPtr->putScalar(0.0);
1899  InterpICQVecPtr->putScalar(0.0);
1900  InterpICStoreVecPtr->putScalar(0.0);
1901 
1902  InterpICSolVecPtr->linearCombo(-1.0, *firstSolVecPtr, 1.0, *secondSolVecPtr );
1903  InterpICSolVecPtr->linearCombo(1.0, *firstSolVecPtr, fraction , *InterpICSolVecPtr);
1904 
1905  InterpICStateVecPtr->linearCombo(-1.0, *firstStateVecPtr, 1.0, *secondStateVecPtr );
1906  InterpICStateVecPtr->linearCombo(1.0, *firstStateVecPtr, fraction , *InterpICStateVecPtr);
1907 
1908  InterpICQVecPtr->linearCombo(-1.0, *firstQVecPtr, 1.0, *secondQVecPtr );
1909  InterpICQVecPtr->linearCombo(1.0, *firstQVecPtr, fraction , *InterpICQVecPtr);
1910 
1911  InterpICStoreVecPtr->linearCombo(-1.0, *firstStoreVecPtr, 1.0, *secondStoreVecPtr );
1912  InterpICStoreVecPtr->linearCombo(1.0, *firstStoreVecPtr, fraction , *InterpICStoreVecPtr);
1913 
1914  goodSolutionVec_.push_back(InterpICSolVecPtr);
1915  goodStateVec_.push_back(InterpICStateVecPtr);
1916  goodQVec_.push_back(InterpICQVecPtr);
1917  goodStoreVec_.push_back(InterpICStoreVecPtr);
1918  }
1919 
1920  // Clean up the fast time data since we are finished computing the initial condition.
1921  // The fast time data can take a considerable amount of memory for large problems.
1922  dsPtr->resetFastTimeData();
1923 
1924  return true;
1925 }
1926 
1927 } // namespace Analysis
1928 } // namespace Xyce
1929