Xyce  6.1
N_ANP_AC.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-2015 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 //-----------------------------------------------------------------------------
27 // Filename : $RCSfile: N_ANP_AC.C,v $
28 // Purpose : AC analysis functions.
29 // Special Notes :
30 // Creator : Ting Mei
31 // Creation Date : 7/11
32 //
33 // Revision Information:
34 // ---------------------
35 // Revision Number: $Revision: 1.136 $
36 // Revision Date : $Date: 2015/09/03 22:11:39 $
37 // Current Owner : $Author: rlschie $
38 //-----------------------------------------------------------------------------
39 
40 #include <Xyce_config.h>
41 
42 #include <iomanip>
43 
44 #include <N_ANP_AC.h>
45 
46 #include <N_ANP_AnalysisManager.h>
47 #include <N_ANP_OutputMgrAdapter.h>
48 #include <N_ERH_Message.h>
49 #include <N_IO_CircuitBlock.h>
50 #include <N_IO_CmdParse.h>
51 #include <N_IO_InitialConditions.h>
52 #include <N_IO_OptionBlock.h>
53 #include <N_IO_PkgOptionsMgr.h>
54 #include <N_IO_SpiceSeparatedFieldTool.h>
55 #include <N_LAS_BlockMatrix.h>
56 #include <N_LAS_BlockSystemHelpers.h>
57 #include <N_LAS_BlockVector.h>
58 #include <N_LAS_Builder.h>
59 #include <N_LAS_Matrix.h>
60 #include <N_LAS_MultiVector.h>
61 #include <N_LAS_System.h>
62 #include <N_LOA_Loader.h>
63 #include <N_NLS_Manager.h>
64 #include <N_TIA_DataStore.h>
65 #include <N_TIA_StepErrorControl.h>
68 #include <N_UTL_Diagnostic.h>
69 #include <N_TOP_Topology.h>
70 #include <N_UTL_Factory.h>
71 #include <N_UTL_FeatureTest.h>
72 #include <N_UTL_LogStream.h>
73 #include <N_UTL_Math.h>
74 #include <N_UTL_Timer.h>
75 
76 #include <N_PDS_ParMap.h>
77 
78 #include <Epetra_SerialComm.h>
79 #include <Epetra_Map.h>
80 #include <Epetra_BlockMap.h>
81 #include <Epetra_CrsMatrix.h>
82 #include <Epetra_CrsGraph.h>
83 #include <Epetra_MultiVector.h>
84 #include <Epetra_Vector.h>
85 #include <Epetra_Export.h>
86 #include <Epetra_LinearProblem.h>
87 #include <Amesos.h>
88 
89 #include<N_UTL_ExtendedString.h>
90 #include<N_UTL_ExpressionData.h>
91 #include<N_NLS_ReturnCodes.h>
92 
93 #include <Teuchos_RCP.hpp>
94 using Teuchos::RCP;
95 using Teuchos::rcp;
96 
97 namespace Xyce {
98 namespace Analysis {
99 
100 //-----------------------------------------------------------------------------
101 // Function : AnalysisManager::setTimeIntegratorOptions
102 // Purpose :
103 // Special Notes : These are from '.options timeint'
104 // Scope : public
105 // Creator : Eric R. Keiter, SNL, Parallel Computational Sciences
106 // Creation Date : 04/18/02
107 //-----------------------------------------------------------------------------
109  const Util::OptionBlock & option_block)
110 {
111  for (Util::ParamList::const_iterator it = option_block.begin(), end = option_block.end(); it != end; ++it)
112  {
113  const Util::Param &param = *it;
114 
115  if (param.uTag() == "DEBUGLEVEL" )
116  IO::setTimeIntegratorDebugLevel(analysisManager_.getCommandLine(), param.getImmutableValue<int>());
117  else if (nonlinearManager_.setReturnCodeOption(param))
118  ;
119  else if (tiaParams_.setTimeIntegratorOption(param))
120  ;
121  else if (setDCOPOption(param))
122  ;
123  else
124  Report::UserError() << param.uTag() << " is not a recognized time integration option";
125  }
126 
127  return true;
128 }
129 
130 
131 //-----------------------------------------------------------------------------
132 // Function : AC::AC
133 // Purpose : constructor
134 // Special Notes :
135 // Scope : public
136 // Creator : Ting Mei
137 // Creation Date : 5/11
138 //-----------------------------------------------------------------------------
140  AnalysisManager & analysis_manager,
141  Linear::System & linear_system,
142  Nonlinear::Manager & nonlinear_manager,
143  Loader::Loader & loader,
144  Topo::Topology & topology,
145  IO::InitialConditionsManager & initial_conditions_manager)
146  : AnalysisBase(analysis_manager, "AC"),
147  StepEventListener(&analysis_manager),
148  analysisManager_(analysis_manager),
149  loader_(loader),
150  linearSystem_(linear_system),
151  nonlinearManager_(nonlinear_manager),
152  topology_(topology),
153  initialConditionsManager_(initial_conditions_manager),
154  outputManagerAdapter_(analysis_manager.getOutputManagerAdapter()),
155  bVecRealPtr(linearSystem_.builder().createVector()),
156  bVecImagPtr(linearSystem_.builder().createVector()),
157  ACMatrix_(0),
158  B_(0),
159  X_(0),
160  blockSolver_(0),
161  blockProblem_(0),
162  acLoopSize_(0),
163  stepMult_(0.0),
164  fstep_(0.0),
165  type_("DEC"),
166  stepFlag_(false),
167  np_(10.0),
168  fStart_(1.0),
169  fStop_(1.0),
170  currentFreq_(0.0)
171 {
172  bVecRealPtr->putScalar(0.0);
173  bVecImagPtr->putScalar(0.0);
174 }
175 
176 //-----------------------------------------------------------------------------
177 // Function : AC::~AC()
178 // Purpose : destructor
179 // Special Notes :
180 // Scope : public
181 // Creator : Ting Mei
182 // Creation Date : 5/11
183 //-----------------------------------------------------------------------------
185 {
186  delete bVecRealPtr;
187  delete bVecImagPtr;
188  delete ACMatrix_;
189  delete B_;
190  delete X_;
191  delete blockSolver_;
192  delete blockProblem_;
193 }
194 
195 //-----------------------------------------------------------------------------
196 // Function : AC::notify
197 // Purpose :
198 // Special Notes :
199 // Scope : public
200 // Creator :
201 // Creation Date :
202 //-----------------------------------------------------------------------------
203 void AC::notify(const StepEvent &event)
204 {
205  if (event.state_ == StepEvent::STEP_STARTED)
206  {
208 
209  stepFlag_ = true;
211 
212  bVecRealPtr->putScalar(0.0);
213  bVecImagPtr->putScalar(0.0);
214  }
215 }
216 
217 //-----------------------------------------------------------------------------
218 // Function : AC::setAnalysisParams
219 // Purpose :
220 // Special Notes : These are from the .AC statement.
221 // Scope : public
222 // Creator : Ting Mei, SNL
223 // Creation Date : 6/11
224 //-----------------------------------------------------------------------------
225 bool
227  const Util::OptionBlock & paramsBlock)
228 {
229  for (Util::ParamList::const_iterator it = paramsBlock.begin(), end = paramsBlock.end(); it != end; ++it)
230  {
231  if ((*it).uTag() == "TYPE")
232  {
233  type_ = (*it).stringValue();
234  }
235  else if ((*it).uTag() == "NP")
236  {
237  np_ = (*it).getImmutableValue<double>();
238  }
239  else if ((*it).uTag() == "FSTART")
240  {
241  fStart_ = (*it).getImmutableValue<double>();
242  }
243  else if ((*it).uTag() == "FSTOP")
244  {
245  fStop_ = (*it).getImmutableValue<double>();
246  }
247  }
248 
249  if (DEBUG_ANALYSIS && isActive(Diag::TIME_PARAMETERS))
250  {
251  dout() << section_divider << std::endl
252  << "AC simulation parameters"
253  << std::endl
254  << "number of points = " << np_ << std::endl
255  << "starting frequency = " << fStart_ << std::endl
256  << "stop frequency = " << fStop_
257  << std::endl;
258  }
259 
260  return true;
261 }
262 
263 //-----------------------------------------------------------------------------
264 // Function : AC::getDCOPFlag()
265 // Purpose :
266 // Special Notes :
267 // Scope : public
268 // Creator : Eric Keiter, SNL
269 // Creation Date : 3/24/2014
270 //-----------------------------------------------------------------------------
271 bool AC::getDCOPFlag() const
272 {
274 }
275 
276 //-----------------------------------------------------------------------------
277 // Function : AC::run()
278 // Purpose :
279 // Special Notes :
280 // Scope : public
281 // Creator : Ting Mei, SNL
282 // Creation Date : 5/11
283 //-----------------------------------------------------------------------------
284 bool AC::doRun()
285 {
286  return doInit() && doLoopProcess() && doFinish();
287 }
288 
289 //-----------------------------------------------------------------------------
290 // Function : AC::init()
291 // Purpose :
292 // Special Notes :
293 // Scope : public
294 // Creator : Ting Mei, SNL
295 // Creation Date : 5/11
296 //-----------------------------------------------------------------------------
298 {
299  bool bsuccess = true;
300 
302 
303  // Get set to do the operating point.
306 
307  stepNumber = 0;
309 
310  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::INITIALIZE, AnalysisEvent::AC_IC));
311 
312  // set initial guess, if there is one to be set.
313  // this setInitialGuess call is to up an initial guess in the
314  // devices that have them (usually PDE devices). This is different than
315  // the "intializeProblem" call, which sets IC's. (initial conditions are
316  // different than initial guesses.
318 
319  // If available, set initial solution (.IC, .NODESET, etc).
320  // setInputOPFlag(
321  // outputManagerAdapter_.setupInitialConditions(
322  // *analysisManager_.getDataStore()->nextSolutionPtr,
323  // *analysisManager_.getDataStore()->flagSolutionPtr));
324 
327  topology_.getSolutionNodeNameMap(),
330 
331  // Set a constant history for operating point calculation
335 
336  // solving for DC op
343 
345  {
346  Report::UserError() << "Solving for DC operating point failed! Cannot continue AC analysis";
347  return false;
348  }
349 
350  // only output DC op if the .op was specified in the netlist
351  // or if this AC analysis is called from a .step loop
353  {
355  stepNumber,
362  objectiveVec_,
365  }
366 
367  // This for saving the data from the DC op. different from the above where we are
368  // concerned with generating normal output.
369  // outputManagerAdapter_.outputDCOP(*analysisManager_.getDataStore()->nextSolutionPtr);
371 
373 
374  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::FINISH, AnalysisEvent::AC_IC));
375 
376  return bsuccess;
377 }
378 
379 
380 //-----------------------------------------------------------------------------
381 // Function : AC::loopProcess()
382 // Purpose : Conduct the stepping loop.
383 // Special Notes :
384 // Scope : public
385 // Creator : Ting Mei, SNL
386 // Creation Date : 5/11
387 //-----------------------------------------------------------------------------
389 {
390  analysisManager_.getDataStore()->daeQVectorPtr->putScalar(0.0);
391  analysisManager_.getDataStore()->daeFVectorPtr->putScalar(0.0);
392 
397 
408  );
409 
434 
439 
442 
443  if (DEBUG_TIME && isActive(Diag::TIME_PARAMETERS))
444  {
445  Xyce::dout() << "dQdxMatrixPtr:" << std::endl;
446  analysisManager_.getDataStore()->dQdxMatrixPtr->printPetraObject( Xyce::dout() );
447 
448  Xyce::dout() << "dFdxMatrixPtr:" << std::endl;
449  analysisManager_.getDataStore()->dFdxMatrixPtr->printPetraObject( Xyce::dout() );
450 
451  Xyce::dout() << std::endl;
452  }
453 
455 
456  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::INITIALIZE, AnalysisEvent::AC));
457 
458  for (int currentStep = 0; currentStep < acLoopSize_; ++currentStep)
459  {
460  updateCurrentFreq_(currentStep);
461 
463 
464  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::STEP_STARTED, AnalysisEvent::AC, currentFreq_, currentStep));
465 
466  bool stepAttemptStatus;
467  {
468  Stats::StatTop _nonlinearStat("Nonlinear Solve");
469  Stats::TimeBlock _nonlinearTimer(_nonlinearStat);
470 
471  stepAttemptStatus = solveLinearSystem_();
472  }
473 
474  if (stepAttemptStatus)
475  {
476  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::STEP_SUCCESSFUL, AnalysisEvent::AC, currentFreq_, currentStep));
478  }
479  else // stepAttemptStatus (ie do this if the step FAILED)
480  {
481  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::STEP_FAILED, AnalysisEvent::AC, currentFreq_, currentStep));
483  }
484  }
485 
486  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::FINISH, AnalysisEvent::AC));
487 
488  return true;
489 }
490 
491 
492 //-----------------------------------------------------------------------------
493 // Function : AC::createLinearSystem_()
494 // Purpose :
495 // Special Notes :
496 // Scope : public
497 // Creator : Ting Mei, Heidi Thornquist, SNL
498 // Creation Date : 6/20/2011
499 //-----------------------------------------------------------------------------
501 {
502  bool bsuccess = true;
503 
504  N_PDS_Manager &pds_manager = *analysisManager_.getPDSManager();
505 
506  RCP<N_PDS_ParMap> baseMap = rcp(pds_manager.getParallelMap( Parallel::SOLUTION ), false);
507  Epetra_CrsGraph &baseFullGraph = *pds_manager.getMatrixGraph(Parallel::JACOBIAN);
508 
509  int numBlocks = 2;
510 
511  RCP<N_PDS_ParMap> blockMap = Linear::createBlockParMap(numBlocks, *baseMap);
512 
513  // Create a block vector
514  delete B_;
515  B_ = new Linear::BlockVector(numBlocks, blockMap, baseMap);
516 
517  // -----------------------------------------------------
518  // Now test block graphs.
519  // -----------------------------------------------------
520 
521  std::vector<std::vector<int> > blockPattern(2);
522  blockPattern[0].resize(2);
523  blockPattern[0][0] = 0; blockPattern[0][1] = 1;
524  blockPattern[1].resize(2);
525  blockPattern[1][0] = 0; blockPattern[1][1] = 1;
526 
527  int offset = Linear::generateOffset( *baseMap );
528 
529  RCP<Epetra_CrsGraph> blockGraph = Linear::createBlockGraph( offset, blockPattern, *blockMap, baseFullGraph);
530 
531  delete ACMatrix_;
532  ACMatrix_ = new Linear::BlockMatrix( numBlocks, offset, blockPattern, *blockGraph, baseFullGraph);
533 
534  // First diagonal block
535  ACMatrix_->put( 0.0 ); // Zero out whole matrix
536  ACMatrix_->block( 0, 0 ).add(*G_);
537 
538  // Second diagonal block
539  ACMatrix_->block( 1, 1 ).add(*G_);
540 
541  B_->putScalar( 0.0 );
542  B_->block( 0 ).addVec( 1.0, *bVecRealPtr);
543  B_->block( 1 ).addVec( 1.0, *bVecImagPtr);
544 
545  Amesos amesosFactory;
546 
547  delete X_;
548  X_ = new Linear::BlockVector (numBlocks, blockMap, baseMap);
549  X_->putScalar( 0.0 );
550 
551  delete blockProblem_;
552  blockProblem_ = new Epetra_LinearProblem(&ACMatrix_->epetraObj(), &X_->epetraObj(), &B_->epetraObj() );
553 
554  delete blockSolver_;
555  blockSolver_ = amesosFactory.Create( "Klu", *blockProblem_ );
556 
557  // Need to reindex the linear system because of the noncontiguous block map.
558  Teuchos::ParameterList params;
559  params.set( "Reindex", true );
560  blockSolver_->SetParameters( params );
561 
562  // Call symbolic factorization without syncronizing the values, since they are not necessary here.
563  int linearStatus = blockSolver_->SymbolicFactorization();
564 
565  if (linearStatus != 0)
566  {
567  Xyce::dout() << "Amesos symbolic factorization exited with error: " << linearStatus;
568  bsuccess = false;
569  }
570 
571  return bsuccess;
572 }
573 
574 //-----------------------------------------------------------------------------
575 // Function : AC::updateLinearSystemFreq_()
576 // Purpose :
577 // Special Notes :
578 // Scope : public
579 // Creator : Ting Mei, Heidi Thornquist, SNL
580 // Creation Date : 6/20/2011
581 //-----------------------------------------------------------------------------
583 {
584  double omega = 2.0 * M_PI * currentFreq_;
585 
586  ACMatrix_->block( 0, 1).put( 0.0);
587  ACMatrix_->block( 0, 1).add(*C_);
588  ACMatrix_->block( 0, 1).scale(-omega);
589 
590  ACMatrix_->block(1, 0).put( 0.0);
591  ACMatrix_->block(1, 0).add(*C_);
592  ACMatrix_->block(1, 0).scale(omega);
593 
594  // Copy the values loaded into the blocks into the global matrix for the solve.
595  ACMatrix_->assembleGlobalMatrix();
596 
597  return true;
598 }
599 
600 //-----------------------------------------------------------------------------
601 // Function : AC::solveLinearSystem_()
602 // Purpose :
603 // Special Notes :
604 // Scope : public
605 // Creator : Ting Mei, Heidi Thornquist, SNL
606 // Creation Date : 6/2011
607 //-----------------------------------------------------------------------------
609 {
610  bool bsuccess = true;
611  // Solve the block problem
612 
613  int linearStatus = blockSolver_->NumericFactorization();
614 
615  if (linearStatus != 0)
616  {
617  Xyce::dout() << "Amesos numeric factorization exited with error: " << linearStatus;
618  bsuccess = false;
619  }
620 
621  linearStatus = blockSolver_->Solve();
622  if (linearStatus != 0)
623  {
624  Xyce::dout() << "Amesos solve exited with error: " << linearStatus;
625  bsuccess = false;
626  }
627 
628  return bsuccess;
629 }
630 
631 //-----------------------------------------------------------------------------
632 // Function : AC::processSuccessfulStep()
633 // Purpose :
634 // Special Notes :
635 // Scope : public
636 // Creator : Ting Mei, SNL
637 // Creation Date :
638 //-----------------------------------------------------------------------------
640 {
641  bool bsuccess = true;
642 
643 // Output x.
644  outputManagerAdapter_.outputAC (currentFreq_, X_->block(0), X_-> block(1));
645 
646  if ( !firstDoubleDCOPStep() )
647  {
648  stepNumber += 1;
651  }
652 
653  // This output call is for device-specific output (like from a PDE device,
654  // outputting mesh-based tecplot files). It will only work in parallel if on
655  // a machine where all processors have I/O capability.
657 
658  return bsuccess;
659 }
660 
661 
662 //-----------------------------------------------------------------------------
663 // Function : AC::processFailedStep
664 // Purpose :
665 // Special Notes :
666 // Scope : public
667 // Creator : Ting Mei, SNL
668 // Creation Date :
669 //-----------------------------------------------------------------------------
671 {
672  bool bsuccess = true;
673 
674  stepNumber += 1;
675  acSweepFailures_.push_back(stepNumber);
678 
679  return bsuccess;
680 }
681 
682 //-----------------------------------------------------------------------------
683 // Function : AC::finish
684 // Purpose :
685 // Special Notes :
686 // Scope : public
687 // Creator : Ting Mei, SNL
688 // Creation Date :
689 //-----------------------------------------------------------------------------
691 {
692  bool bsuccess = true;
693 
694  if (DEBUG_ANALYSIS)
695  {
696  Xyce::dout() << "Calling AC::doFinish() outputs!" << std::endl;
697  }
698 
700  if (!(acSweepFailures_.empty()))
701  {
702  bsuccess = false;
703  }
704 
705  return bsuccess;
706 }
707 
708 
709 //-----------------------------------------------------------------------------
710 // Function : AC::doHandlePredictor
711 // Purpose :
712 // Special Notes :
713 // Scope : private
714 // Creator : Eric Keiter, SNL
715 // Creation Date : 06/24/2013
716 //-----------------------------------------------------------------------------
718 {
722 
723  // In case this is the upper level of a 2-level sim, tell the
724  // inner solve to do its prediction:
725  bool beginIntegrationFlag = analysisManager_.getBeginningIntegrationFlag(); // system_state.beginIntegrationFlag;
726  double nextTimeStep = analysisManager_.getStepErrorControl().currentTimeStep; // system_state.nextTimeStep;
727  double nextTime = analysisManager_.getStepErrorControl().nextTime; // system_state.nextTime;
728  int currentOrder = analysisManager_.getWorkingIntegrationMethod().getOrder(); // system_state.currentOrder;
729 
730  loader_.startTimeStep(beginIntegrationFlag, nextTimeStep, nextTime, currentOrder);
731 
732  return true;
733 }
734 
735 //-----------------------------------------------------------------------------
736 // Function : AC::updateCurrentFreq_(int stepNumber )
737 // Purpose :
738 // Special Notes : Used for AC analysis classes.
739 // Scope : public
740 // Creator : Ting Mei, SNL.
741 // Creation Date :
742 //-----------------------------------------------------------------------------
743 bool AC::updateCurrentFreq_(int stepNumber)
744 {
745  if (type_ == "LIN")
746  {
747 
748  currentFreq_ = fStart_ + static_cast<double>(stepNumber)*fstep_;
749  }
750  else if(type_ == "DEC" || type_ == "OCT")
751  {
752 
753  currentFreq_ = fStart_*pow(stepMult_, static_cast<double>(stepNumber) );
754  }
755  else
756  {
757  Report::DevelFatal().in("AC::updateCurrentFreq_") << "AC::updateCurrentFreq_: unsupported STEP type";
758  }
759 
760  return true;
761 }
762 
763 //-----------------------------------------------------------------------------
764 // Function : AC::setupSweepParam_
765 // Purpose : Processes sweep parameters.
766 // Special Notes : Used for AC analysis classes.
767 // Scope : public
768 // Creator : Ting Mei, SNL.
769 // Creation Date :
770 //-----------------------------------------------------------------------------
772 {
773  double fstart, fstop;
774  double fcount = 0.0;
775 
776  fstart = fStart_;
777  fstop = fStop_;
778 
779  if (DEBUG_ANALYSIS && isActive(Diag::TIME_PARAMETERS))
780  {
781  Xyce::dout() << std::endl << std::endl;
782  Xyce::dout() << section_divider << std::endl;
783  Xyce::dout() << "AC::setupSweepParam_" << std::endl;
784  }
785 
786  if (type_ == "LIN")
787  {
788  int np = static_cast<int>(np_);
789 
790  if ( np == 1)
791  fstep_ = 0;
792  else
793  fstep_ = (fstop - fstart)/(np_ - 1.0);
794 
795  fcount = np_;
796  if (DEBUG_ANALYSIS && isActive(Diag::TIME_PARAMETERS))
797  {
798  Xyce::dout() << "fstep = " << fstep_ << std::endl;
799  }
800  }
801  else if (type_ == "DEC")
802  {
803  stepMult_ = pow(10.0, 1.0/np_);
804  fcount = floor(fabs(log10(fstart) - log10(fstop)) * np_ + 1.0);
805  if (DEBUG_ANALYSIS && isActive(Diag::TIME_PARAMETERS))
806  {
807  Xyce::dout() << "stepMult_ = " << stepMult_ << std::endl;
808  }
809  }
810  else if (type_ == "OCT")
811  {
812  stepMult_ = pow(2.0, 1.0/np_);
813 
814  // changed to remove dependence on "log2" function, which apparently
815  // doesn't exist in the math libraries of FreeBSD or the mingw
816  // cross-compilation suite. Log_2(x)=log_e(x)/log_e(2.0)
817  double ln2 = log(2.0);
818  fcount = floor(fabs(log(fstart) - log(fstop))/ln2 * np_ + 1.0);
819  if (DEBUG_ANALYSIS && isActive(Diag::TIME_PARAMETERS))
820  {
821  Xyce::dout() << "stepMult_ = " << stepMult_ << std::endl;
822  }
823  }
824  else
825  {
826  Report::DevelFatal().in("AC::setupSweepParam") << "Unsupported type";
827  }
828 
829  // At this point, pinterval equals the total number of steps
830  // for the step loop.
831  return static_cast<int> (fcount);
832 }
833 
834 namespace {
835 
836 typedef Util::Factory<AnalysisBase, AC> ACFactoryBase;
837 
838 //-----------------------------------------------------------------------------
839 // Class : ACFactory
840 // Purpose :
841 // Special Notes :
842 // Scope : public
843 // Creator : David G. Baur Raytheon Sandia National Laboratories 1355
844 // Creation Date : Thu Jan 29 12:53:02 2015
845 //-----------------------------------------------------------------------------
846 ///
847 /// Factory for parsing AC parameters from the netlist and creating AC analysis.
848 ///
849 class ACFactory : public ACFactoryBase
850 {
851 public:
852  //-----------------------------------------------------------------------------
853  // Function : ACFactory
854  // Purpose :
855  // Special Notes :
856  // Scope : public
857  // Creator : David G. Baur Raytheon Sandia National Laboratories 1355
858  // Creation Date : Thu Jan 29 12:54:09 2015
859  //-----------------------------------------------------------------------------
860  ///
861  /// Constructs the AC analysis factory
862  ///
863  /// @invariant Stores the results of parsing, so if more than one of the analysis and
864  /// associated lines are parsed, the second options simply overwrite the previously parsed
865  /// values.
866  ///
867  /// @invariant The existence of the parameters specified in the constructor cannot
868  /// change.
869  ///
870  /// @param analysis_manager
871  /// @param linear_system
872  /// @param nonlinear_manager
873  /// @param topology
874  ///
875  ACFactory(
876  Analysis::AnalysisManager & analysis_manager,
877  Linear::System & linear_system,
878  Nonlinear::Manager & nonlinear_manager,
879  Loader::Loader & loader,
880  Topo::Topology & topology,
881  IO::InitialConditionsManager & initial_conditions_manager)
882  : ACFactoryBase(),
883  analysisManager_(analysis_manager),
884  linearSystem_(linear_system),
885  nonlinearManager_(nonlinear_manager),
886  loader_(loader),
887  topology_(topology),
888  initialConditionsManager_(initial_conditions_manager)
889  {}
890 
891  virtual ~ACFactory()
892  {}
893 
894  //-----------------------------------------------------------------------------
895  // Function : create
896  // Purpose :
897  // Special Notes :
898  // Scope : public
899  // Creator : David G. Baur Raytheon Sandia National Laboratories 1355
900  // Creation Date : Thu Jan 29 12:59:00 2015
901  //-----------------------------------------------------------------------------
902  ///
903  /// Create a new AC analysis and applies the analysis and time integrator option blocks.
904  ///
905  /// @return new AC analysis object
906  ///
907  AC *create() const
908  {
909  analysisManager_.setAnalysisMode(ANP_MODE_AC);
911  ac->setAnalysisParams(acAnalysisOptionBlock_);
912  ac->setTimeIntegratorOptions(timeIntegratorOptionBlock_);
913 
914  return ac;
915  }
916 
917  //-----------------------------------------------------------------------------
918  // Function : setACAnalysisOptionBlock
919  // Purpose :
920  // Special Notes :
921  // Scope : public
922  // Creator : David G. Baur Raytheon Sandia National Laboratories 1355
923  // Creation Date : Thu Jan 29 13:00:14 2015
924  //-----------------------------------------------------------------------------
925  ///
926  /// Saves the analysis parsed options block in the factory.
927  ///
928  /// @invariant Overwrites any previously specified analysis option block.
929  ///
930  /// @param option_block parsed option block
931  ///
932  void setACAnalysisOptionBlock(const Util::OptionBlock &option_block)
933  {
934  acAnalysisOptionBlock_ = option_block;
935  }
936 
937  //-----------------------------------------------------------------------------
938  // Function : setTimeIntegratorOptionBlock
939  // Purpose :
940  // Special Notes :
941  // Scope : public
942  // Creator : David G. Baur Raytheon Sandia National Laboratories 1355
943  // Creation Date : Thu Jan 29 13:01:27 2015
944  //-----------------------------------------------------------------------------
945  ///
946  /// Saves the time integrator parsed option block.
947  ///
948  /// @invariant Overwrites any previously specified time integrator option block.
949  ///
950  /// @param option_block parsed option block
951  ///
952  bool setTimeIntegratorOptionBlock(const Util::OptionBlock &option_block)
953  {
954  timeIntegratorOptionBlock_ = option_block;
955 
956  return true;
957  }
958 
959 public:
960  AnalysisManager & analysisManager_;
961  Linear::System & linearSystem_;
962  Nonlinear::Manager & nonlinearManager_;
963  Loader::Loader & loader_;
964  Topo::Topology & topology_;
965  IO::InitialConditionsManager & initialConditionsManager_;
966 
967 private:
968  Util::OptionBlock acAnalysisOptionBlock_;
969  Util::OptionBlock timeIntegratorOptionBlock_;
970 };
971 
972 // .AC
973 struct ACAnalysisReg : public IO::PkgOptionsReg
974 {
975  ACAnalysisReg(
976  ACFactory & factory )
977  : factory_(factory)
978  {}
979 
980  bool operator()(const Util::OptionBlock &option_block)
981  {
982  factory_.setACAnalysisOptionBlock(option_block);
983 
984  factory_.analysisManager_.addAnalysis(&factory_);
985 
986  return true;
987  }
988 
989  ACFactory & factory_;
990 };
991 
992 //-----------------------------------------------------------------------------
993 // Function : extractACData
994 // Purpose : Extract the parameters from a netlist .AC line held in
995 // parsed_line.
996 // Special Notes :
997 // Scope : public
998 // Creator :
999 // Creation Date : 7/23/08
1000 //-----------------------------------------------------------------------------
1001 bool
1002 extractACData(
1003  IO::PkgOptionsMgr & options_manager,
1004  IO::CircuitBlock & circuit_block,
1005  const std::string & netlist_filename,
1006  const IO::TokenVector & parsed_line)
1007 {
1008  Util::OptionBlock option_block("AC", Util::OptionBlock::NO_EXPRESSIONS, netlist_filename, parsed_line[0].lineNumber_);
1009 
1010  int numFields = parsed_line.size();
1011 
1012  // Check that the minimum required number of fields are on the line.
1013  if (numFields != 5)
1014  {
1015  Report::UserError0().at(netlist_filename, parsed_line[0].lineNumber_)
1016  << ".AC line has an unexpected number of fields";
1017  return false;
1018  }
1019 
1020  int linePosition = 1; // Start of parameters on .param line.
1021 
1022  Util::Param parameter("", "");
1023 
1024  // type is required
1025  parameter.setTag( "TYPE" );
1026 // parameter.setVal( parsed_line[linePosition].string_ );
1027  ExtendedString stringVal ( parsed_line[linePosition].string_ );
1028  stringVal.toUpper();
1029  parameter.setVal(std::string(stringVal));
1030  option_block.addParam( parameter );
1031 
1032  ++linePosition; // Advance to next parameter.
1033 
1034  // np is required
1035  parameter.setTag( "NP" );
1036  parameter.setVal( parsed_line[linePosition].string_ );
1037  option_block.addParam( parameter );
1038  ++linePosition; // Advance to next parameter.
1039 
1040  // fstart is required
1041  parameter.setTag( "FSTART" );
1042  parameter.setVal( parsed_line[linePosition].string_ );
1043  option_block.addParam( parameter );
1044  ++linePosition; // Advance to next parameter.
1045 
1046 
1047  // fstop is required
1048  parameter.setTag( "FSTOP" );
1049  parameter.setVal( parsed_line[linePosition].string_ );
1050  option_block.addParam( parameter );
1051 // ++linePosition; // Advance to next parameter.
1052 
1053  circuit_block.addOptions(option_block);
1054 
1055  return true;
1056 }
1057 
1058 } // namespace <unnamed>
1059 
1060 
1061 bool
1063  FactoryBlock & factory_block)
1064 {
1065  ACFactory *factory = new ACFactory(factory_block.analysisManager_, factory_block.linearSystem_, factory_block.nonlinearManager_, factory_block.loader_, factory_block.topology_, factory_block.initialConditionsManager_);
1066 
1067  addAnalysisFactory(factory_block, factory);
1068 
1069  factory_block.optionsManager_.addCommandParser(".AC", extractACData);
1070 
1071  factory_block.optionsManager_.addCommandProcessor("AC", new ACAnalysisReg(*factory));
1072 
1073  factory_block.optionsManager_.addOptionsProcessor("TIMEINT", IO::createRegistrationOptions(*factory, &ACFactory::setTimeIntegratorOptionBlock));
1074 
1075  return true;
1076 }
1077 
1078 } // namespace Analysis
1079 } // namespace Xyce
bool doProcessSuccessfulStep()
Definition: N_ANP_AC.C:639
Linear::Vector * lastSolutionPtr
bool doHandlePredictor()
Definition: N_ANP_AC.C:717
bool registerACFactory(FactoryBlock &factory_block)
Definition: N_ANP_AC.C:1062
bool getDCOPFlag() const
Definition: N_ANP_AC.C:271
bool createLinearSystem_()
Definition: N_ANP_AC.C:500
unsigned int successfulStepsTaken_
Number of consecutive successful time-integration steps.
virtual bool setInitialGuess(Linear::Vector *solVectorPtr)
Definition: N_LOA_Loader.h:232
Linear::Vector * lastLeadDeltaVPtr
bool setAnalysisParams(const Util::OptionBlock &paramsBlock)
Definition: N_ANP_AC.C:226
bool doLoopProcess()
Definition: N_ANP_AC.C:388
Linear::Matrix * C_
Definition: N_ANP_AC.h:163
IO::InitialConditionsManager & initialConditionsManager_
Definition: N_ANP_AC.C:965
bool setTimeIntegratorOption(const Util::Param &param)
ACFactory & factory_
Definition: N_ANP_AC.C:989
void outputAC(double freq, const Linear::Vector &solnVecRealPtr, const Linear::Vector &solnVecImaginaryPtr)
void evaluateStepError(const Loader::Loader &loader, const TIAParams &tia_params)
Pure virtual class to augment a linear system.
TimeIntg::TIAParams tiaParams_
Definition: N_ANP_AC.h:144
Util::OptionBlock timeIntegratorOptionBlock_
Definition: N_ANP_AC.C:969
virtual bool loadDAEVectors(Linear::Vector *nextSolVectorPtr, Linear::Vector *currSolVectorPtr, Linear::Vector *lastSolVectorPtr, Linear::Vector *nextStaVectorPtr, Linear::Vector *currStaVectorPtr, Linear::Vector *lastStaVectorPtr, Linear::Vector *StaDerivVectorPtr, Linear::Vector *nextStoVectorPtr, Linear::Vector *currStoVectorPtr, Linear::Vector *lastStoVectorPtr, Linear::Vector *stoLeadCurrQVectorPtr, Linear::Vector *nextLeadFVectorPtr, Linear::Vector *currLeadFVectorPtr, Linear::Vector *lastLeadFVectorPtr, Linear::Vector *nextLeadQVectorPtr, Linear::Vector *nextJunctionVVectorPtr, Linear::Vector *currJunctionVVectorPtr, Linear::Vector *lastJunctionVVectorPtr, Linear::Vector *QVectorPtr, Linear::Vector *FVectorPtr, Linear::Vector *BVectorPtr, Linear::Vector *dFdxdVpVectorPtr, Linear::Vector *dQdxdVpVectorPtr)
Definition: N_LOA_Loader.h:122
virtual bool outputPlotFiles() const
Definition: N_LOA_Loader.h:261
unsigned int stepNumber
Time-integration step number counter.
virtual bool loadBVectorsforAC(Linear::Vector *bVecRealPtr, Linear::Vector *bVecImagPtr)
Definition: N_LOA_Loader.h:203
bool doProcessFailedStep()
Definition: N_ANP_AC.C:670
Linear::Vector * currStorePtr
Linear::Vector * lastStatePtr
Topo::Topology & topology_
Definition: N_ANP_AC.C:964
void gatherStepStatistics(StatCounts &stats, Nonlinear::NonLinearSolver &nonlinear_solver, int newton_convergence_status)
double currentFreq_
Definition: N_ANP_AC.h:161
virtual ~AC()
Definition: N_ANP_AC.C:184
void createTimeIntegratorMethod(const TimeIntg::TIAParams &tia_params, const unsigned int integration_method)
bool resetAll(const TIAParams &tia_params)
std::vector< int > acSweepFailures_
Definition: N_ANP_AC.h:151
Util::ListenerAutoSubscribe< StepEvent > StepEventListener
Linear::Vector * nextLeadDeltaVPtr
Linear::Vector * currStatePtr
unsigned int failedStepsAttempted_
Total number of failed time-integration steps.
Parallel::Manager * getPDSManager() const
bool setDCOPOption(const Util::Param &param)
TimeIntg::StepErrorControl & getStepErrorControl()
bool setTimeIntegratorOptions(const Util::OptionBlock &option_block)
Definition: N_ANP_AC.C:108
Linear::Vector * bVecRealPtr
Definition: N_ANP_AC.h:146
Linear::BlockVector * B_
Definition: N_ANP_AC.h:166
int setupSweepParam_()
Definition: N_ANP_AC.C:771
Topo::Topology & topology_
Definition: N_ANP_AC.h:141
NonLinearSolver & getNonlinearSolver()
std::vector< double > scaled_dOdpVec_
Definition: N_ANP_AC.h:175
Linear::Vector * daeQVectorPtr
std::string type_
Definition: N_ANP_AC.h:154
Linear::Vector * lastLeadCurrentPtr
virtual bool startTimeStep(bool beginIntegrationFlag, double nextTimeStep, double nextTime, int currentOrder)
Definition: N_LOA_Loader.h:305
void setErrorWtVector(const TIAParams &tia_params)
virtual bool isPDESystem() const
Definition: N_LOA_Loader.h:256
std::vector< double > dOdpAdjVec_
Definition: N_ANP_AC.h:174
Amesos_BaseSolver * blockSolver_
Definition: N_ANP_AC.h:169
Linear::Vector * daeFVectorPtr
Linear::Vector * currLeadCurrentQPtr
Nonlinear::Manager & nonlinearManager_
Definition: N_ANP_AC.h:140
AC(AnalysisManager &analysis_manager, Linear::System &linear_system, Nonlinear::Manager &nonlinear_manager, Loader::Loader &loader, Topo::Topology &topology, IO::InitialConditionsManager &initial_conditions_manager)
Definition: N_ANP_AC.C:139
OutputMgrAdapter & outputManagerAdapter_
Definition: N_ANP_AC.h:143
The FactoryBlock contains parameters needed by the analysis creation functions.
IO::InitialConditionsManager & initialConditionsManager_
Linear::Matrix * dFdxMatrixPtr
virtual bool updateSources()
Definition: N_LOA_Loader.h:244
Linear::Vector * nextStatePtr
Linear::Vector * dFdxdVpVectorPtr
virtual bool loadDAEMatrices(Linear::Vector *tmpSolVectorPtr, Linear::Vector *tmpStaVectorPtr, Linear::Vector *tmpStaDerivVectorPtr, Linear::Vector *tmpStoVectorPtr, Linear::Matrix *tmpdQdxMatrixPtr, Linear::Matrix *tmpdFdxMatrixPtr)
Definition: N_LOA_Loader.h:110
virtual bool updateState(Linear::Vector *nextSolVectorPtr, Linear::Vector *currSolVectorPtr, Linear::Vector *lastSolVectorPtr, Linear::Vector *nextStaVectorPtr, Linear::Vector *currStaVectorPtr, Linear::Vector *lastStaVectorPtr, Linear::Vector *nextStoVectorPtr, Linear::Vector *currStoVectorPtr, Linear::Vector *lastStoVectorPtr)
Definition: N_LOA_Loader.h:189
Linear::Vector * currLeadDeltaVPtr
AnalysisManager & analysisManager_
Definition: N_ANP_AC.C:960
void setDoubleDCOPEnabled(bool enable)
Nonlinear::Manager & nonlinearManager_
bool updateCurrentFreq_(int stepNumber)
Definition: N_ANP_AC.C:743
std::vector< double > scaled_dOdpAdjVec_
Definition: N_ANP_AC.h:176
IO::InitialConditionsManager & initialConditionsManager_
Definition: N_ANP_AC.h:142
Linear::Vector * nextSolutionPtr
Linear::BlockVector * X_
Definition: N_ANP_AC.h:167
void dcOutput(int dcStepNumber, Linear::Vector &currSolutionPtr, Linear::Vector &stateVecPtr, Linear::Vector &storeVecPtr, Linear::Vector &lead_current_vector, Linear::Vector &junction_voltage_vector, Linear::Vector &lead_current_dqdt_vector, std::vector< double > &objectiveVec_, std::vector< double > &dOdpVec_, std::vector< double > &dOdpAdjVec_, std::vector< double > &scaled_dOdpVec_, std::vector< double > &scaled_dOdpAdjVec_)
#define M_PI
Linear::Vector * nextStorePtr
bool updateLinearSystemFreq_()
Definition: N_ANP_AC.C:582
void addAnalysisFactory(FactoryBlock &factory_block, Util::Factory< AnalysisBase, void > *factory)
Linear::Vector * dQdxdVpVectorPtr
Util::OptionBlock acAnalysisOptionBlock_
Definition: N_ANP_AC.C:968
Linear::Vector * currLeadCurrentPtr
Linear::Vector * lastStorePtr
void setInputOPFlag(bool initial_conditions_loaded)
TimeIntg::WorkingIntegrationMethod & getWorkingIntegrationMethod()
const IO::CmdParse & getCommandLine() const
AnalysisManager & analysisManager_
Definition: N_ANP_AC.h:137
std::vector< double > objectiveVec_
Definition: N_ANP_AC.h:172
Linear::Vector * bVecImagPtr
Definition: N_ANP_AC.h:147
Linear::Vector * nextLeadCurrentQDerivPtr
Linear::Vector * currSolutionPtr
TimeIntg::DataStore * getDataStore()
Linear::BlockMatrix * ACMatrix_
Definition: N_ANP_AC.h:165
unsigned int baseIntegrationMethod_
Current time-integration method flag.
Linear::Vector * nextStoreLeadCurrQPtr
Linear::Vector * nextStateDerivPtr
Loader::Loader & loader_
Definition: N_ANP_AC.C:963
std::vector< double > dOdpVec_
Definition: N_ANP_AC.h:173
Nonlinear::Manager & nonlinearManager_
Definition: N_ANP_AC.C:962
Linear::Matrix * dQdxMatrixPtr
bool solveLinearSystem_()
Definition: N_ANP_AC.C:608
Epetra_LinearProblem * blockProblem_
Definition: N_ANP_AC.h:170
Linear::Vector * flagSolutionPtr
void notify(const StepEvent &event)
Definition: N_ANP_AC.C:203
Loader::Loader & loader_
Definition: N_ANP_AC.h:138
Linear::Vector * daeBVectorPtr
Linear::Matrix * G_
Definition: N_ANP_AC.h:164
bool setReturnCodeOption(const Util::Param &param)
Linear::System & linearSystem_
Definition: N_ANP_AC.C:961
Linear::Vector * nextLeadCurrentPtr