Xyce  6.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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-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 //-----------------------------------------------------------------------------
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.67.2.1 $
36 // Revision Date : $Date: 2014/08/26 22:31:06 $
37 // Current Owner : $Author: dgbaur $
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_IO_CmdParse.h>
48 #include <N_IO_OutputMgr.h>
49 #include <N_LAS_BlockMatrix.h>
50 #include <N_LAS_BlockSystemHelpers.h>
51 #include <N_LAS_BlockVector.h>
52 #include <N_LAS_Builder.h>
53 #include <N_LAS_Matrix.h>
54 #include <N_LAS_MultiVector.h>
55 #include <N_LAS_System.h>
56 #include <N_LOA_Loader.h>
57 #include <N_MPDE_Manager.h>
58 #include <N_TIA_DataStore.h>
59 #include <N_TIA_StepErrorControl.h>
61 #include <N_UTL_LogStream.h>
62 #include <N_UTL_Misc.h>
63 #include <N_UTL_Timer.h>
64 
65 #include <N_PDS_ParMap.h>
66 #ifdef Xyce_PARALLEL_MPI
67 #include <N_PDS_ParComm.h>
68 #include <mpi.h>
69 #else
70 #include <N_PDS_SerialComm.h>
71 #endif
72 
73 #include <Epetra_SerialComm.h>
74 #include <Epetra_Map.h>
75 #include <Epetra_BlockMap.h>
76 #include <Epetra_CrsMatrix.h>
77 #include <Epetra_CrsGraph.h>
78 #include <Epetra_MultiVector.h>
79 #include <Epetra_Vector.h>
80 #include <Epetra_Export.h>
81 #include <Epetra_LinearProblem.h>
82 #include <Amesos.h>
83 
84 #include<N_UTL_ExpressionData.h>
85 #include<N_NLS_ReturnCodes.h>
86 
87 namespace Xyce {
88 namespace Analysis {
89 
90 //-----------------------------------------------------------------------------
91 // Function : AC::AC( AnalysisManager * )
92 // Purpose : constructor
93 // Special Notes :
94 // Scope : public
95 // Creator : Ting Mei
96 // Creation Date : 5/11
97 //-----------------------------------------------------------------------------
98 AC::AC( AnalysisManager &analysis_manager )
99  : AnalysisBase(analysis_manager),
100  Util::ListenerAutoSubscribe<StepEvent>(&analysis_manager),
101  bVecRealPtr(0),
102  bVecImagPtr(0),
103  acLoopSize_(0),
104  stepMult_(0.0),
105  fstep_(0.0),
106  currentFreq_(0.0)
107 {
108  bVecRealPtr = linearSystem_.builder().createVector();
109  bVecImagPtr = linearSystem_.builder().createVector();
110 
111  bVecRealPtr->putScalar(0.0);
112  bVecImagPtr->putScalar(0.0);
113 
114 }
115 
116 //-----------------------------------------------------------------------------
117 // Function : AC::~AC()
118 // Purpose : destructor
119 // Special Notes :
120 // Scope : public
121 // Creator : Ting Mei
122 // Creation Date : 5/11
123 //-----------------------------------------------------------------------------
125 {
126 
127  if (bVecRealPtr) { delete bVecRealPtr; bVecRealPtr=0; }
128 
129  if (bVecImagPtr) { delete bVecImagPtr; bVecImagPtr=0; }
130 
131 }
132 
133 void AC::notify(const StepEvent &event)
134 {
135  if (event.state_ == StepEvent::STEP_STARTED)
136  {
138 
139  bVecRealPtr->putScalar(0.0);
140  bVecImagPtr->putScalar(0.0);
141  }
142 }
143 
144 
145 //-----------------------------------------------------------------------------
146 // Function : AC::setAnalysisParams
147 // Purpose :
148 // Special Notes : These are from the .AC statement.
149 // Scope : public
150 // Creator : Ting Mei, SNL
151 // Creation Date : 6/11
152 //-----------------------------------------------------------------------------
153 bool AC::setAnalysisParams(const N_UTL_OptionBlock & paramsBlock)
154 {
155  Util::ParameterList::const_iterator it_tp;
156  Util::ParameterList::const_iterator first = paramsBlock.getParams().begin();
157  Util::ParameterList::const_iterator last = paramsBlock.getParams().end();
158  for (it_tp = first; it_tp != last; ++it_tp)
159  {
160  if (it_tp->uTag() == "TYPE")
161  {
162  tiaParams_.type = it_tp->stringValue();
163 // tiaParams_.tStartGiven = true;
164  }
165  else if (it_tp->uTag() == "NP")
166  {
167  tiaParams_.np = it_tp->getImmutableValue<double>();
168  }
169  else if (it_tp->uTag() == "FSTART")
170  {
171  tiaParams_.fStart = it_tp->getImmutableValue<double>();
172  }
173  else if (it_tp->uTag() == "FSTOP")
174  {
175  tiaParams_.fStop = it_tp->getImmutableValue<double>();
176  }
177  }
178 
179 // tiaParams_.debugLevel = 1;
180  if (DEBUG_ANALYSIS && tiaParams_.debugLevel > 0)
181  {
182  dout() << section_divider << std::endl
183  << "AC simulation parameters"
184  //<< Util::push << std::endl
185  << std::endl
186  << "number of points = " << tiaParams_.np << std::endl
187  << "starting frequency = " << tiaParams_.fStart << std::endl
188  << "stop frequency = " << tiaParams_.fStop
189  //<< Util::pop << std::endl;
190  << std::endl;
191  }
192 
193  return true;
194 }
195 
196 //-----------------------------------------------------------------------------
197 // Function : AC::getDCOPFlag()
198 // Purpose :
199 // Special Notes :
200 // Scope : public
201 // Creator : Eric Keiter, SNL
202 // Creation Date : 3/24/2014
203 //-----------------------------------------------------------------------------
205 {
206  return ((getIntegrationMethod())==TIAMethod_NONE);
207 }
208 
209 //-----------------------------------------------------------------------------
210 // Function : AC::run()
211 // Purpose :
212 // Special Notes :
213 // Scope : public
214 // Creator : Ting Mei, SNL
215 // Creation Date : 5/11
216 //-----------------------------------------------------------------------------
217 bool AC::run()
218 {
219  bool bsuccess = true;
220 
221  bsuccess = bsuccess & init();
222  bsuccess = bsuccess & loopProcess();
223 
224  // if processing the loop failed,
225  // then skip finish step
226  if( bsuccess )
227  {
228  bsuccess = bsuccess & finish();
229  }
230 
231  return bsuccess;
232 }
233 
234 //-----------------------------------------------------------------------------
235 // Function : AC::init()
236 // Purpose :
237 // Special Notes :
238 // Scope : public
239 // Creator : Ting Mei, SNL
240 // Creation Date : 5/11
241 //-----------------------------------------------------------------------------
242 bool AC::init()
243 {
244  bool bsuccess = true;
245 
247 
249 
250  // Get set to do the operating point.
253 
254  stepNumber = 0;
255  doubleDCOPFlag_ = loader_.getDoubleDCOPFlag ();
257 
258  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::INITIALIZE, AnalysisEvent::AC_IC));
259 
260  // set initial guess, if there is one to be set.
261  // this setInitialGuess call is to up an initial guess in the
262  // devices that have them (usually PDE devices). This is different than
263  // the "intializeProblem" call, which sets IC's. (initial conditions are
264  // different than initial guesses.
266 
267  // If available, set initial solution (.IC, .NODESET, etc).
269 
270  // Set a constant history for operating point calculation
274 
275  // solving for DC op
276  handlePredictor();
277  loader_.updateSources();
282 
284  {
285  Report::UserError() << "Solving for DC operating point failed! Cannot continue AC analysis";
286  return false;
287  }
288 
289  // only output DC op if the .op was specified in the netlist
290  // or if this AC analysis is called from a .step loop
292  {
294  stepNumber,
298  objectiveVec_,
301  }
302 
303  // This for saving the data from the DC op. different from the above where we are
304  // concerned with generating normal output.
306 
307  loader_.loadBVectorsforAC (bVecRealPtr, bVecImagPtr);
308 
309  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::FINISH, AnalysisEvent::AC_IC));
310 
311  return bsuccess;
312 }
313 
314 
315 //-----------------------------------------------------------------------------
316 // Function : AC::loopProcess()
317 // Purpose : Conduct the stepping loop.
318 // Special Notes :
319 // Scope : public
320 // Creator : Ting Mei, SNL
321 // Creation Date : 5/11
322 //-----------------------------------------------------------------------------
324 {
327 
332 
333 // integrationMethod_ = 6;
334 
335  loader_.updateState
345  );
346 
347  loader_.loadDAEVectors
364 
369 
372 
373  if (tiaParams_.debugLevel > 0)
374  {
375  Xyce::dout() << "dQdxMatrixPtr:" << std::endl;
376  analysisManager_.getTIADataStore()->dQdxMatrixPtr->printPetraObject( Xyce::dout() );
377 
378  Xyce::dout() << "dFdxMatrixPtr:" << std::endl;
379  analysisManager_.getTIADataStore()->dFdxMatrixPtr->printPetraObject( Xyce::dout() );
380 
381  Xyce::dout() << std::endl;
382  }
383 
385 
386  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::INITIALIZE, AnalysisEvent::AC));
387 
388  for (int currentStep = 0; currentStep < acLoopSize_; ++currentStep)
389  {
390  updateCurrentFreq_(currentStep);
391 
393 
394  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::STEP_STARTED, AnalysisEvent::AC, currentFreq_, currentStep));
395 
396  bool stepAttemptStatus = solveLinearSystem_();
397 
398  if (stepAttemptStatus)
399  {
400  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::STEP_SUCCESSFUL, AnalysisEvent::AC, currentFreq_, currentStep));
402  }
403  else // stepAttemptStatus (ie do this if the step FAILED)
404  {
405  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::STEP_FAILED, AnalysisEvent::AC, currentFreq_, currentStep));
407  }
408  }
409 
410  static_cast<Xyce::Util::Notifier<AnalysisEvent> &>(analysisManager_).publish(AnalysisEvent(AnalysisEvent::FINISH, AnalysisEvent::AC));
411 
412  return true;
413 }
414 
415 
416 //-----------------------------------------------------------------------------
417 // Function : AC::createLinearSystem_()
418 // Purpose :
419 // Special Notes :
420 // Scope : public
421 // Creator : Ting Mei, Heidi Thornquist, SNL
422 // Creation Date : 6/20/2011
423 //-----------------------------------------------------------------------------
425 {
426  bool bsuccess = true;
427 
428  RCP<N_PDS_Manager> pdsMgrPtr_;
429  pdsMgrPtr_ = rcp(linearSystem_.getPDSManager(), false);
430 
431  RCP<N_PDS_ParMap> baseMap;
432  RCP<Epetra_CrsGraph> BaseFullGraph_;
433 
434  baseMap = rcp(pdsMgrPtr_->getParallelMap( "SOLUTION" ), false);
435  BaseFullGraph_ = rcp( pdsMgrPtr_->getMatrixGraph("JACOBIAN"), false );
436 
437  int numBlocks = 2;
438 
439  RCP<N_PDS_ParMap> blockMap = createBlockParMap(numBlocks, *baseMap);
440 
441  // Create a block vector
442  BPtr_ = rcp ( new N_LAS_BlockVector ( numBlocks, blockMap, baseMap ) );
443 
444  // -----------------------------------------------------
445  // Now test block graphs.
446  // -----------------------------------------------------
447 
448  std::vector<std::vector<int> > blockPattern(2);
449  blockPattern[0].resize(2);
450  blockPattern[0][0] = 0; blockPattern[0][1] = 1;
451  blockPattern[1].resize(2);
452  blockPattern[1][0] = 0; blockPattern[1][1] = 1;
453 
454  int offset = generateOffset( *baseMap );
455 
456  RCP<Epetra_CrsGraph> blockGraph = createBlockGraph( offset, blockPattern, *blockMap, *BaseFullGraph_);
457 
458  ACMatrixPtr_ = rcp ( new N_LAS_BlockMatrix( numBlocks, offset, blockPattern, *blockGraph, *BaseFullGraph_) );
459 
460  // First diagonal block
461  ACMatrixPtr_->put( 0.0 ); // Zero out whole matrix
462  ACMatrixPtr_->block( 0, 0 ).add(*GPtr_);
463  // Second diagonal block
464  ACMatrixPtr_->block( 1, 1 ).add(*GPtr_);
465 
466  BPtr_->putScalar( 0.0 );
467  BPtr_->block( 0 ).addVec( 1.0, *bVecRealPtr);
468  BPtr_->block( 1 ).addVec( 1.0, *bVecImagPtr);
469 
470  Amesos amesosFactory;
471 
472  XPtr_ = rcp ( new N_LAS_BlockVector (numBlocks, blockMap, baseMap) );
473  XPtr_->putScalar( 0.0 );
474 
475  blockProblem = rcp(new Epetra_LinearProblem(&ACMatrixPtr_->epetraObj(), &XPtr_->epetraObj(), &BPtr_->epetraObj() ) );
476 
477  blockSolver = rcp( amesosFactory.Create( "Klu", *blockProblem ) );
478 
479  // Need to reindex the linear system because of the noncontiguous block map.
480  Teuchos::ParameterList params;
481  params.set( "Reindex", true );
482  blockSolver->SetParameters( params );
483 
484  // Call symbolic factorization without syncronizing the values, since they are not necessary here.
485  int linearStatus = blockSolver->SymbolicFactorization();
486 
487  if (linearStatus != 0)
488  {
489  Xyce::dout() << "Amesos symbolic factorization exited with error: " << linearStatus;
490  bsuccess = false;
491  }
492 
493  return bsuccess;
494 
495 }
496 
497 //-----------------------------------------------------------------------------
498 // Function : AC::updateLinearSystemFreq_()
499 // Purpose :
500 // Special Notes :
501 // Scope : public
502 // Creator : Ting Mei, Heidi Thornquist, SNL
503 // Creation Date : 6/20/2011
504 //-----------------------------------------------------------------------------
505 
507 {
508  double omega;
509 
510  omega = 2.0 * M_PI * currentFreq_;
511 
512  ACMatrixPtr_->block( 0, 1).put( 0.0);
513  ACMatrixPtr_->block( 0, 1).add(*CPtr_);
514  ACMatrixPtr_->block( 0, 1).scale(-omega);
515 
516  ACMatrixPtr_->block(1, 0).put( 0.0);
517  ACMatrixPtr_->block(1, 0).add(*CPtr_);
518  ACMatrixPtr_->block(1, 0).scale(omega);
519 
520  // Copy the values loaded into the blocks into the global matrix for the solve.
521  ACMatrixPtr_->assembleGlobalMatrix();
522 
523  return true;
524 }
525 
526 //-----------------------------------------------------------------------------
527 // Function : AC::solveLinearSystem_()
528 // Purpose :
529 // Special Notes :
530 // Scope : public
531 // Creator : Ting Mei, Heidi Thornquist, SNL
532 // Creation Date : 6/2011
533 //-----------------------------------------------------------------------------
534 
536 {
537 
538  bool bsuccess = true;
539  // Solve the block problem
540 
541  int linearStatus = blockSolver->NumericFactorization();
542 
543  if (linearStatus != 0)
544  {
545  Xyce::dout() << "Amesos numeric factorization exited with error: " << linearStatus;
546  bsuccess = false;
547  }
548 
549  linearStatus = blockSolver->Solve();
550  if (linearStatus != 0)
551  {
552  Xyce::dout() << "Amesos solve exited with error: " << linearStatus;
553  bsuccess = false;
554  }
555 
556  return bsuccess;
557 }
558 
559 //-----------------------------------------------------------------------------
560 // Function : AC::processSuccessfulStep()
561 // Purpose :
562 // Special Notes :
563 // Scope : public
564 // Creator : Ting Mei, SNL
565 // Creation Date :
566 //-----------------------------------------------------------------------------
568 {
569  bool bsuccess = true;
570 
571 // Output x.
572  outputManagerAdapter_.outputAC (currentFreq_, XPtr_->block(0), XPtr_-> block(1));
573 
574  if ( !firstDoubleDCOPStep_() )
575  {
576  stepNumber += 1;
579  }
580 
581  // This output call is for device-specific output (like from a PDE device,
582  // outputting mesh-based tecplot files). It will only work in parallel if on
583  // a machine where all processors have I/O capability.
584  loader_.output();
585 
586  return bsuccess;
587 }
588 
589 
590 //-----------------------------------------------------------------------------
591 // Function : AC::processFailedStep
592 // Purpose :
593 // Special Notes :
594 // Scope : public
595 // Creator : Ting Mei, SNL
596 // Creation Date :
597 //-----------------------------------------------------------------------------
599 {
600  bool bsuccess = true;
601 
602  stepNumber += 1;
603  acSweepFailures_.push_back(stepNumber);
606 
607  return bsuccess;
608 }
609 
610 //-----------------------------------------------------------------------------
611 // Function : AC::finish
612 // Purpose :
613 // Special Notes :
614 // Scope : public
615 // Creator : Ting Mei, SNL
616 // Creation Date :
617 //-----------------------------------------------------------------------------
619 {
620  bool bsuccess = true;
621 
622 #ifdef Xyce_DEBUG_ANALYSIS
623  Xyce::dout() << "Calling AC::finish() outputs!" << std::endl;
624 #endif
626  if (!(acSweepFailures_.empty()))
627  {
628  bsuccess = false;
629  }
630 
631  return bsuccess;
632 }
633 
634 
635 //-----------------------------------------------------------------------------
636 // Function : AC::handlePredictor
637 // Purpose :
638 // Special Notes :
639 // Scope : private
640 // Creator : Eric Keiter, SNL
641 // Creation Date : 06/24/2013
642 //-----------------------------------------------------------------------------
644 {
648 
649  // In case this is the upper level of a 2-level sim, tell the
650  // inner solve to do its prediction:
651  loader_.startTimeStep ();
652 
653  return true;
654 }
655 
656 //-----------------------------------------------------------------------------
657 // Function : AC::updateCurrentFreq_(int stepNumber )
658 // Purpose :
659 // Special Notes : Used for AC analysis classes.
660 // Scope : public
661 // Creator : Ting Mei, SNL.
662 // Creation Date :
663 //-----------------------------------------------------------------------------
664 bool AC::updateCurrentFreq_(int stepNumber)
665 {
666  double fstart;
667  fstart = tiaParams_.fStart;
668 
669 // tiaParams_.debugLevel = 1;
670 
671  if (tiaParams_.type=="LIN")
672  {
673 
674  currentFreq_ = fstart + static_cast<double>(stepNumber)*fstep_;
675  }
676  else if(tiaParams_.type=="DEC" || tiaParams_.type=="OCT")
677  {
678 
679  currentFreq_ = fstart*pow(stepMult_, static_cast<double>(stepNumber) );
680  }
681  else
682  {
683  Report::DevelFatal().in("AC::updateCurrentFreq_") << "AC::updateCurrentFreq_: unsupported STEP type";
684  }
685 
687  {
688  Xyce::dout() << "currentFreq_ = " << currentFreq_ << std::endl;
689  }
690 
691  return true;
692 }
693 
694 //-----------------------------------------------------------------------------
695 // Function : AC::setupSweepParam_
696 // Purpose : Processes sweep parameters.
697 // Special Notes : Used for AC analysis classes.
698 // Scope : public
699 // Creator : Ting Mei, SNL.
700 // Creation Date :
701 //-----------------------------------------------------------------------------
703 {
704  double fstart, fstop;
705  double fcount = 0.0;
706 
707  fstart = tiaParams_.fStart;
708  fstop = tiaParams_.fStop;
709 
710 #ifdef Xyce_DEBUG_ANALYSIS
712  {
713  Xyce::dout() << std::endl << std::endl;
714  Xyce::dout() << section_divider << std::endl;
715  Xyce::dout() << "AC::setupSweepParam_" << std::endl;
716  }
717 #endif
718 
719  if (tiaParams_.type=="LIN")
720  {
721  int np = static_cast<int> (tiaParams_.np);
722 
723  if ( np == 1)
724  fstep_ = 0;
725  else
726  fstep_ = (fstop - fstart)/(tiaParams_.np - 1);
727 
728  fcount = tiaParams_.np;
729 #ifdef Xyce_DEBUG_ANALYSIS
731  {
732  Xyce::dout() << "fstep = " << fstep_ << std::endl;
733  }
734 #endif
735  }
736  else if(tiaParams_.type=="DEC")
737  {
738  stepMult_ = pow(10,(1/tiaParams_.np));
739  fcount = floor(fabs(log10(fstart) - log10(fstop)) * tiaParams_.np + 1);
740 #ifdef Xyce_DEBUG_ANALYSIS
742  {
743  Xyce::dout() << "stepMult_ = " << stepMult_ << std::endl;
744  }
745 #endif
746  }
747  else if(tiaParams_.type=="OCT")
748  {
749  stepMult_ = pow(2,1/(tiaParams_.np));
750 
751  // changed to remove dependence on "log2" function, which apparently
752  // doesn't exist in the math libraries of FreeBSD or the mingw
753  // cross-compilation suite. Log_2(x)=log_e(x)/log_e(2.0)
754  double ln2=log(2.0);
755  fcount = floor(fabs(log(fstart) - log(fstop))/ln2 * tiaParams_.np + 1);
756 #ifdef Xyce_DEBUG_ANALYSIS
758  {
759  Xyce::dout() << "stepMult_ = " << stepMult_ << std::endl;
760  }
761 #endif
762  }
763  else
764  {
765  Report::DevelFatal().in("AC::setupSweepParam") << "Unsupported type";
766  }
767 
768  // At this point, pinterval equals the total number of steps
769  // for the step loop.
770  return static_cast<int> (fcount);
771 }
772 
773 } // namespace Analysis
774 } // namespace Xyce