Xyce  6.1
N_LOA_HBLoader.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 // Filename : $RCSfile: N_LOA_HBLoader.C,v $
27 // Purpose :
28 // Special Notes :
29 // Creator :
30 // Creation Date :
31 //
32 // Revision Information:
33 // ---------------------
34 // Revision Number: $Revision: 1.57 $
35 // Revision Date : $Date: 2015/09/24 00:13:11 $
36 // Current Owner : $Author: hkthorn $
37 //-----------------------------------------------------------------------------
38 
39 #include <Xyce_config.h>
40 
41 
42 // ---------- Standard Includes ----------
43 #include <iostream>
44 
45 // ---------- Xyce Includes ----------
46 
47 #include <N_LOA_HBLoader.h>
48 #include <N_MPDE_Discretization.h>
49 #include <N_ERH_ErrorMgr.h>
50 
51 #include <N_LAS_Builder.h>
52 #include <N_LAS_HBBuilder.h>
53 #include <N_LAS_Matrix.h>
54 #include <N_LAS_BlockVector.h>
55 #include <N_LAS_BlockMatrix.h>
56 #include <N_LAS_BlockSystemHelpers.h>
57 #include <N_LAS_SystemHelpers.h>
58 #include <N_LAS_MultiVector.h>
59 
60 #include <N_UTL_fwd.h>
61 #include <N_UTL_FeatureTest.h>
62 #include <N_UTL_Math.h>
63 
64 #include <N_PDS_ParMap.h>
65 #include <N_DEV_DeviceMgr.h>
66 
67 #include <Epetra_BlockMap.h>
68 #include <Epetra_CrsGraph.h>
69 #include <Epetra_CrsMatrix.h>
70 
71 using Teuchos::rcp;
72 using Teuchos::RCP;
73 using Teuchos::rcp_dynamic_cast;
74 
75 namespace Xyce {
76 namespace Loader {
77 
79  const Teuchos::RCP<const N_MPDE_Discretization> discPtr,
80  Device::DeviceMgr & device_manager,
81  Linear::Builder & builder)
82  : CktLoader(device_manager),
83  deviceManager_(device_manager),
84  builder_(builder),
85  fastTimeDiscPtr_(discPtr),
86  periodicTimesOffset_(0),
87  period_(0),
88  matrixFreeFlag_(false)
89 {
90  // Now initialize all the time domain working vectors.
91  appVecPtr_ = rcp(builder_.createVector());
92  appNextStaVecPtr_ = rcp(builder_.createStateVector());
93  appCurrStaVecPtr_ = rcp(builder_.createStateVector());
94  appLastStaVecPtr_ = rcp(builder_.createStateVector());
95 
96  appdQdxPtr_ = rcp(builder_.createMatrix());
97  appdFdxPtr_ = rcp(builder_.createMatrix());
98 
99  appNextStoVecPtr_ = rcp(builder_.createStoreVector());
100  appCurrStoVecPtr_ = rcp(builder_.createStoreVector());
101  appLastStoVecPtr_ = rcp(builder_.createStoreVector());
102  appStoLeadCurrQVecPtr_ = rcp(builder_.createStoreVector());
103 
104  appNextLeadFVecPtr_ = rcp(builder.createLeadCurrentVector());
105  appCurrLeadFVecPtr_ = rcp(builder.createLeadCurrentVector());
106  appLastLeadFVecPtr_ = rcp(builder.createLeadCurrentVector());
107  appLeadQVecPtr_ = rcp(builder.createLeadCurrentVector());
108  appNextJunctionVVecPtr_ = rcp(builder.createLeadCurrentVector());
109  appCurrJunctionVVecPtr_ = rcp(builder.createLeadCurrentVector());
110  appLastJunctionVVecPtr_ = rcp(builder.createLeadCurrentVector());
111 }
112 
113 //-----------------------------------------------------------------------------
114 // Function : HBLoader::registerHBBuilder
115 // Purpose : Registration method for the HB builder
116 // Special Notes :
117 // Scope : public
118 // Creator : Heidi Thornquist, Sandia Labs
119 // Creation Date : 11/08/13
120 //-----------------------------------------------------------------------------
121 void HBLoader::registerHBBuilder( Teuchos::RCP<Linear::HBBuilder> hbBuilderPtr )
122 {
123  hbBuilderPtr_ = hbBuilderPtr;
124 
125  // Now initialize all the frequency domain working vectors.
126  bXtPtr_ = hbBuilderPtr_->createTimeDomainBlockVector();
127  bVtPtr_ = hbBuilderPtr_->createTimeDomainBlockVector();
128  bmdQdxPtr_ = rcp_dynamic_cast<Linear::BlockMatrix>(rcp(hbBuilderPtr_->createMatrix()));
129  bmdFdxPtr_ = rcp_dynamic_cast<Linear::BlockMatrix>(rcp(hbBuilderPtr_->createMatrix()));
130 
131  // Vectors related to lead currents
132  bStoreVecFreqPtr_ = hbBuilderPtr_->createExpandedRealFormTransposeStoreBlockVector();
133  bStoreLeadCurrQVecFreqPtr_ = hbBuilderPtr_->createExpandedRealFormTransposeStoreBlockVector();
134 
135  bLeadCurrentVecFreqPtr_ = hbBuilderPtr_->createExpandedRealFormTransposeLeadCurrentBlockVector();
136  bLeadCurrentQVecFreqPtr_ = hbBuilderPtr_->createExpandedRealFormTransposeLeadCurrentBlockVector();
137 }
138 
139 //-----------------------------------------------------------------------------
140 // Function : HBLoader::setFreq
141 // Purpose : Assign times for fast time scale
142 // Special Notes :
143 // Scope : public
144 // Creator :
145 // Creation Date :
146 //-----------------------------------------------------------------------------
147 void HBLoader::setHBFreqs( const std::vector<double> & freqs)
148 {
149  freqs_ = freqs;
150 }
151 
152 
153 
154 //-----------------------------------------------------------------------------
155 // Function : HBLoader::setFastTimes
156 // Purpose : Assign times for fast time scale
157 // Special Notes :
158 // Scope : public
159 // Creator :
160 // Creation Date :
161 //-----------------------------------------------------------------------------
162 void HBLoader::setFastTimes( const std::vector<double> & times )
163 {
164  times_ = times;
165 }
166 
167 //-----------------------------------------------------------------------------
168 // Function : HBLoader::loadDAEMatrices
169 // Purpose :
170 // Special Notes :
171 // Scope : public
172 // Creator :
173 // Creation Date :
174 //-----------------------------------------------------------------------------
175 bool HBLoader::loadDAEMatrices( Linear::Vector * X,
176  Linear::Vector * S,
177  Linear::Vector * dSdt,
178  Linear::Vector * Store,
179  Linear::Matrix * dQdx,
180  Linear::Matrix * dFdx)
181 {
182  if (matrixFreeFlag_)
183  {
184  if (DEBUG_HB)
185  {
186  Xyce::dout() << std::endl
187  << Xyce::section_divider << std::endl
188  << " HBLoader::loadDAEMatrices: matrixFree case" << std::endl;
189  }
190 
191  dQdx->put(0.0);
192  dFdx->put(0.0);
193  // Do nothing in the Matrix Free Case.
194  return(true);
195  }
196  else
197  {
198  Report::DevelFatal0().in("HBLoader::loadDAEMatrices") << "This function actually was called in a non matrix free case";
199 
200  if (DEBUG_HB)
201  {
202  Xyce::dout() << std::endl
203  << Xyce::section_divider << std::endl
204  << " HBLoader::loadDAEMatrices: Time dependent with matrix case" << std::endl;
205  }
206 
207  return loadTimeDepDAEMatrices(X,S,dSdt, Store, dQdx,dFdx);
208  }
209 
210 }
211 
212 //-----------------------------------------------------------------------------
213 // Function : HBLoader::loadTimeDepDAEMatrices
214 // Purpose :
215 // Special Notes :
216 // Scope : public
217 // Creator : Todd Coffey, 1414
218 // Creation Date : 09/03/08
219 //-----------------------------------------------------------------------------
220 bool HBLoader::loadTimeDepDAEMatrices( Linear::Vector * X,
221  Linear::Vector * S,
222  Linear::Vector * dSdt,
223  Linear::Vector * Store,
224  Linear::Matrix * dQdx,
225  Linear::Matrix * dFdx)
226 {
227  if (DEBUG_HB)
228  {
229  Xyce::dout() << std::endl
230  << Xyce::section_divider << std::endl
231  << " HBLoader::loadTimeDepDAEMatrices" << std::endl;
232  }
233 
234  //Zero out matrices
235  dQdx->put(0.0);
236  dFdx->put(0.0);
237 
238  Linear::Vector appdSdt( *appNextStaVecPtr_ );
239 
240  Linear::BlockMatrix & bdQdx = *dynamic_cast<Linear::BlockMatrix*>(dQdx);
241  Linear::BlockMatrix & bdFdx = *dynamic_cast<Linear::BlockMatrix*>(dFdx);
242  Linear::BlockVector & bX = *dynamic_cast<Linear::BlockVector*>(X);
243 #ifdef Xyce_FLEXIBLE_DAE_LOADS
244  Linear::BlockVector & bS = *dynamic_cast<Linear::BlockVector*>(S);
245  Linear::BlockVector & bdSdt = *dynamic_cast<Linear::BlockVector*>(dSdt);
246  Linear::BlockVector & bStore = *dynamic_cast<Linear::BlockVector*>(Store);
247 #endif // Xyce_FLEXIBLE_DAE_LOADS
248 
249  int BlockCount = bX.blockCount();
250 
251  for( int i = 0; i < BlockCount; ++i )
252  {
253  if (DEBUG_HB)
254  {
255  Xyce::dout() << "Processing diagonal matrix block " << i << " of " << BlockCount-1 << std::endl;
256  }
257 
258 #ifdef Xyce_FLEXIBLE_DAE_LOADS
259  //Set Time for fast time scale somewhere
260 // state_.fastTime = times_[i];
262 
263  //Update the sources
264  appLoaderPtr_->updateSources();
265 
266  *appVecPtr_ = bX.block(i);
267  *appNextStaVecPtr_ = bS.block(i);
268  appdSdt = bdSdt.block(i);
269  *appNextStoVecPtr_ = bStore.block(i);
270 
271  appLoaderPtr_->loadDAEMatrices( &*appVecPtr_, &*appNextStaVecPtr_, &appdSdt, &*appNextStoVecPtr_, &*appdQdxPtr_, &*appdFdxPtr_);
272 
273  bdQdx.block(i,i).add( *appdQdxPtr_ );
274  bdFdx.block(i,i).add( *appdFdxPtr_ );
275 #else
276  //For now, the matrices are loaded during the loadDAEVectors method
277  //Just copied here
278  bdQdx.block(i,i).add( bmdQdxPtr_->block(i,i) );
279  bdFdx.block(i,i).add( bmdFdxPtr_->block(i,i) );
280 
281 #endif // Xyce_FLEXIBLE_DAE_LOADS
282  }
283 
284  // Fast Time scale discretization terms:
285  // These are d(dQ/dt1)/dx terms, but go into bdFdx.
286  // For this procedure, need to re-use the app matrix, appdQdx.
287  Linear::Matrix & dQdxdt = *appdQdxPtr_;
288 
289  int Start = fastTimeDiscPtr_->Start();
290  int Width = fastTimeDiscPtr_->Width();
291 
292  const std::vector<double> & Coeffs = fastTimeDiscPtr_->Coeffs();
293 
294  for( int i = 0; i < BlockCount; ++i )
295  {
296  if (DEBUG_HB)
297  {
298  Xyce::dout() << "Processing off diagonal matrix blocks on row " << i << " of " << BlockCount-1 << std::endl;
299  }
300 
301  int Loc;
302  int indexT1 = i + Start + periodicTimesOffset_;
303  int indexT2 = indexT1 + Width - 1;
304  double invh2 = 1.0 / (periodicTimes_[indexT2] - periodicTimes_[indexT1]);
305 
306  for( int j = 0; j < Width; ++j )
307  {
308  Loc = i + (j+Start);
309 
310  if( Loc < 0 )
311  {
312  Loc += BlockCount;
313  }
314  else if( Loc > (BlockCount-1) )
315  {
316  Loc -= BlockCount;
317  }
318 
319  dQdxdt.put(0.0);
320  dQdxdt.add( bdQdx.block(Loc,Loc) );
321  dQdxdt.scale( Coeffs[j]*invh2 );
322  bdFdx.block(i,Loc).add( dQdxdt );
323  }
324  }
325 
326  dQdx->fillComplete();
327  dFdx->fillComplete();
328 
329  if (DEBUG_HB)
330  {
331  Xyce::dout() << "HB bX:" << std::endl;
332  bX.printPetraObject(dout());
333  Xyce::dout() << "HB bdQdx:" << std::endl;
334  bdQdx.printPetraObject(dout());
335  Xyce::dout() << "HB bdFdx:" << std::endl;
336  bdFdx.printPetraObject(dout());
337 #ifdef Xyce_FLEXIBLE_DAE_LOADS
338  Xyce::dout() << "HB bS:" << std::endl;
339  bS.printPetraObject(dout());
340  Xyce::dout() << "HB dSdt:" << std::endl;
341  bdSdt.printPetraObject(dout());
342  Xyce::dout() << "HB bStore:" << std::endl;
343  bStore.printPetraObject(dout());
344 #endif // Xyce_FLEXIBLE_DAE_LOADS
345 
346  Xyce::dout() << Xyce::section_divider << std::endl;
347  }
348 
349  return true;
350 }
351 
352 //-----------------------------------------------------------------------------
353 // Function : HBLoader::applyDAEMatrices
354 // Purpose :
355 // Special Notes :
356 // Scope : public
357 // Creator :
358 // Creation Date :
359 //-----------------------------------------------------------------------------
360 bool HBLoader::applyDAEMatrices( Linear::Vector * Xf,
361  Linear::Vector * S,
362  Linear::Vector * dSdt,
363  Linear::Vector * Store,
364  const Linear::Vector & Vf,
365  Linear::Vector * dQdxV,
366  Linear::Vector * dFdxV )
367 {
368  if (!matrixFreeFlag_)
369  {
370  Report::DevelFatal0().in("HBLoader::applyDAEMatrices") << "This function should only be called in the matrix free case";
371  }
372  if (DEBUG_HB)
373  {
374  Xyce::dout() << std::endl
375  << Xyce::section_divider << std::endl
376  << " HBLoader::applyDAEMatrices" << std::endl;
377  }
378 
379  Linear::BlockVector & bXf = *dynamic_cast<Linear::BlockVector*>(Xf);
380  // We have to do something special with Vf because AztecOO (or Belos)
381  // probably used the Epetra_LinearProblem's Epetra_Maps to create the input
382  // vector here. In this case, Vf is just an Linear::Vector and not a
383  // Linear::BlockVector.
384  const Linear::BlockVector bVf(Vf, bXf.blockSize());
385 
386  // Calculate contributions of the portion of the Jacobian relating to the nonlinear devices.
387  permutedIFT(bVf, &*bVtPtr_);
388 
389  Linear::BlockVector * bdQdxV = dynamic_cast<Linear::BlockVector*>(dQdxV);
390  Linear::BlockVector * bdFdxV = dynamic_cast<Linear::BlockVector*>(dFdxV);
391 
392  Teuchos::RCP<Linear::BlockVector> bdQdxVt = hbBuilderPtr_->createTimeDomainBlockVector();
393  Teuchos::RCP<Linear::BlockVector> bdFdxVt = hbBuilderPtr_->createTimeDomainBlockVector();
394 
395  int BlockCount = bVtPtr_->blockCount();
396 
397  for( int i = 0; i < BlockCount; ++i )
398  {
399  if (DEBUG_HB)
400  {
401  Xyce::dout() << "Processing diagonal matrix block " << i << " of " << BlockCount-1 << std::endl;
402  }
403 
404  if (DEBUG_HB)
405  {
406  Xyce::dout() << "bVtPtr block i = " << i << " : " << std::endl;
407  bVtPtr_->block(i).printPetraObject(dout());
408 
409  std::cout << " old way app matrices " << std::endl;
410 
411  Xyce::dout() << "appdQdxPtr_ = " << i << " : " << std::endl;
412  appdQdxPtr_->printPetraObject(dout());
413 
414  Xyce::dout() << "appdFdxPtr_ = " << i << " : " << std::endl;
415  appdFdxPtr_->printPetraObject(dout());
416  }
417 
418  if (DEBUG_HB)
419  {
420  Xyce::dout() << "bdQdxVt block (old style) i = " << i << " : " << std::endl;
421  bdQdxVt->block(i).printPetraObject(dout());
422 
423  Xyce::dout() << "bdFdxVt block (old style) i = " << i << " : " << std::endl;
424  bdFdxVt->block(i).printPetraObject(dout());
425  }
426 
427 
428  // Get the already stored nonlinear time-domain Jacobian matrices
429  if ( vecNLAppdQdxPtr_.size() != 0 )
430  {
431  vecNLAppdQdxPtr_[ i ]->matvec(false, bVtPtr_->block(i), bdQdxVt->block(i));
432  vecNLAppdFdxPtr_[ i ]->matvec(false, bVtPtr_->block(i), bdFdxVt->block(i));
433  }
434 
435  }
436 
437  permutedFFT(*bdQdxVt, bdQdxV);
438  permutedFFT(*bdFdxVt, bdFdxV);
439 
440  int blockCount = bXf.blockCount();
441  int blockSize = bXf.block(0).globalLength();
442 
443 
444  int size_ = freqs_.size();
445  int posFreq = (size_-1)/2;
446  double omega = 2.0 * M_PI * freqs_[posFreq];
447 
448  // Add in contributions from part of Jacobian related to the linear devices.
449  // Act on bVf or Vf, see if we can put it into a MultiVector and
450  // apply linAppdQdxPtr_ and linAppdFdxPtr_ to the transpose.
451  // Linear::MultiVector lindQdxV( *(bVtPtr_->blockPmap()), bXf.blockSize()/2 );
452 
453  Teuchos::RCP<Linear::BlockVector> permlindQdxV = hbBuilderPtr_->createExpandedRealFormTransposeBlockVector();
454  Teuchos::RCP<Linear::BlockVector> permlindFdxV = hbBuilderPtr_->createExpandedRealFormTransposeBlockVector();
455 
456  if (!Teuchos::is_null(linAppdFdxPtr_) )
457  applyLinearMatrices (bVf, *permlindQdxV, *permlindFdxV);
458 
459  if(!Teuchos::is_null(linAppdFdxPtr_))
460  {
461  bdQdxV->update (1.0, *permlindQdxV, 1.0);
462  bdFdxV->update (1.0, *permlindFdxV, 1.0);
463 
464  }
465  Teuchos::RCP<Linear::BlockVector> bQVec = hbBuilderPtr_->createExpandedRealFormTransposeBlockVector();
466 
467  for( int i = 0; i < blockCount; ++i )
468  {
469  // QVec needs to be created here since only one processor owns each block
470  // and we do not know which one it is in parallel.
471  Linear::Vector& QVec = bQVec->block(i);
472  Linear::Vector& freqVec = bdQdxV->block(i);
473 
474  // Only one processor owns each block of the frequency-domain vector
475  if (freqVec.localLength() > 0)
476  {
477  omega = 2.0 * M_PI * freqs_[posFreq];
478 
479  QVec[0] = -freqVec[1]*omega;
480  QVec[1] = freqVec[0]*omega;
481 
482  for (int j=1; j < (blockSize/2+1)/2; ++j)
483  {
484  omega = 2.0 * M_PI * freqs_[posFreq+j];
485  QVec[2*j] = -freqVec[2*j+1]*omega;
486  QVec[2*(blockSize/2-j)] = -freqVec[2*j+1]*omega;
487 
488  QVec[2*j+1] = freqVec[2*j]*omega;
489  QVec[2*(blockSize/2-j)+1] = -freqVec[2*j]*omega;
490  }
491  }
492 
493  bdFdxV->block(i).update(1.0, QVec , 1.0);
494  }
495 
496  if (DEBUG_HB)
497  {
498  Xyce::dout() << "HB bdQdxV:" << std::endl;
499  bdQdxV->printPetraObject(dout());
500  Xyce::dout() << "HB bdFdxV:" << std::endl;
501  bdFdxV->printPetraObject(dout());
502 
503  Xyce::dout() << Xyce::section_divider << std::endl;
504  }
505 
506  return true;
507 }
508 
509 //-----------------------------------------------------------------------------
510 // Function : HBLoader::updateState
511 // Purpose :
512 // Special Notes : ERK. This function needs to be a no-op. The reason
513 // is that the state information needs to be the same
514 // at the time of updateState, loadDAEVectors and
515 // loadDAEMatrices. Thus, they have to all happen inside
516 // of the same "fast time" loop. So, this functionality
517 // has been moved up into loadDAEVectors.
518 //
519 // Note: for similar reasons, loadDAEMatrices is called from
520 // within that function as well.
521 //
522 // Scope : public
523 // Creator : Todd Coffey, 1414
524 // Creation Date : 01/17/07
525 //-----------------------------------------------------------------------------
527  (Linear::Vector * nextSolVectorPtr,
528  Linear::Vector * currSolVectorPtr,
529  Linear::Vector * lastSolVectorPtr,
530  Linear::Vector * nextStaVectorPtr,
531  Linear::Vector * currStaVectorPtr,
532  Linear::Vector * lastStaVectorPtr,
533  Linear::Vector * nextStoVectorPtr,
534  Linear::Vector * currStoVectorPtr,
535  Linear::Vector * lastStoVectorPtr
536  )
537 {
538  bool bsuccess = true;
539 
540  // For HB case, this needs to be a no-op.
541 
542  return bsuccess;
543 }
544 //-----------------------------------------------------------------------------
545 // Function : HBLoader::applyLinearMatrices
546 // Purpose : apply matrices related to linear nodes
547 // Special Notes : This is done in frequency domain
548 // Scope : public
549 //
550 // Creator : Ting Mei and Heidi Thornquist
551 // Creation Date : 04/10/2015
552 //-----------------------------------------------------------------------------
553 bool HBLoader::applyLinearMatrices( const Linear::BlockVector & bVf,
554  Linear::BlockVector & permlindQdxV,
555  Linear::BlockVector & permlindFdxV )
556 {
557 
558  Linear::MultiVector lindQdxV( *(bVtPtr_->blockPmap()), bVf.blockSize()/2 );
559  Linear::MultiVector lindFdxV( *(bVtPtr_->blockPmap()), bVf.blockSize()/2 );
560 
561  int first = bVf.startBlock();
562  int last = bVf.endBlock();
563  Linear::MultiVector permVf( *(bVtPtr_->blockPmap()), bVf.blockSize()/2 );
564 
565  // Now copy over data for the blocks that this processor owns.
566  for (int i=first, lid=0; i<last; i++, lid++)
567  {
568  Linear::Vector & currBlock = bVf.block(i);
569 
570  // Insert zero-th value of the Fourier expansion, only need real value.
571  permVf[0][lid] = currBlock[0];
572  for (int j=1; j<currBlock.localLength()/2; j++)
573  {
574  permVf[j][lid] = currBlock[j+1];
575  }
576  }
577 
578  linAppdQdxPtr_->matvec(false, permVf, lindQdxV);
579  linAppdFdxPtr_->matvec(false, permVf, lindFdxV);
580 
581  int numharms = bVtPtr_->blockCount();
582 
583  // Now copy over data for the blocks that this processor owns.
584  for (int i=first, lid=0; i<last; i++, lid++)
585  {
586  Linear::Vector & currBlock = permlindQdxV.block(i);
587  Linear::Vector & currBlockF = permlindFdxV.block(i);
588 
589  currBlock[0] = lindQdxV[0][lid];
590  currBlock[1] = 0.0;
591 
592  currBlockF[0] = lindFdxV[0][lid];
593  currBlockF[1] = 0.0;
594 
595  for (int j=1; j< (numharms + 1)/2; j++)
596  {
597  currBlock[2*j] = lindQdxV[2*j-1][lid];
598  currBlock[2*(numharms-j)] = lindQdxV[2*j-1][lid];
599 
600  currBlockF[2*j] = lindFdxV[2*j-1][lid];
601  currBlockF[2*(numharms-j)] = lindFdxV[2*j-1][lid];
602 
603  currBlock[2*j+1] = lindQdxV[2*j][lid];
604  currBlock[2*( numharms - j) + 1] = -lindQdxV[2*j][lid];
605 
606  currBlockF[2*j+1] = lindFdxV[2*j][lid];
607  currBlockF[2*( numharms - j) + 1] = -lindFdxV[2*j][lid];
608  }
609  }
610 
611  return true;
612 }
613 
614 
615 //-----------------------------------------------------------------------------
616 // Function : HBLoader::loadDAEVectors
617 // Purpose :
618 // Special Notes :
619 // Scope : public
620 // Creator :
621 // Creation Date :
622 //-----------------------------------------------------------------------------
623 bool HBLoader::loadDAEVectors( Linear::Vector * Xf,
624  Linear::Vector * currX,
625  Linear::Vector * lastX,
626  Linear::Vector * S,
627  Linear::Vector * currS,
628  Linear::Vector * lastS,
629  Linear::Vector * dSdt,
630  Linear::Vector * Store,
631  Linear::Vector * currStore,
632  Linear::Vector * lastStore,
633  Linear::Vector * storeLeadCurrQ,
634  Linear::Vector * nextLeadFVectorPtr,
635  Linear::Vector * currLeadFVectorPtr,
636  Linear::Vector * lastLeadFVectorPtr,
637  Linear::Vector * nextLeadQVectorPtr,
638  Linear::Vector * nextJunctionVVectorPtr,
639  Linear::Vector * currJunctionVVectorPtr,
640  Linear::Vector * lastJunctionVVectorPtr,
641  Linear::Vector * Q,
642  Linear::Vector * F,
643  Linear::Vector * B,
644  Linear::Vector * dFdxdVp,
645  Linear::Vector * dQdxdVp )
646 {
647  if (DEBUG_HB)
648  {
649  Xyce::dout() << std::endl
650  << Xyce::section_divider << std::endl
651  << " HBLoader::loadDAEVectors" << std::endl;
652  }
653 
654  //Zero out vectors
655  appVecPtr_->putScalar(0.0);
656  appNextStaVecPtr_->putScalar(0.0);
657  appCurrStaVecPtr_->putScalar(0.0);
658  appLastStaVecPtr_->putScalar(0.0);
659  Linear::Vector appdSdt( *appNextStaVecPtr_ );
660 
661  appNextStoVecPtr_->putScalar(0.0);
662  appCurrStoVecPtr_->putScalar(0.0);
663  appLastStoVecPtr_->putScalar(0.0);
664  appStoLeadCurrQVecPtr_->putScalar(0.0);
665 
666  Linear::Vector appQ( *appVecPtr_ );
667  Linear::Vector appF( *appVecPtr_ );
668  Linear::Vector appB( *appVecPtr_ );
669 
670  Linear::Vector appdFdxdVp( *appVecPtr_ );
671  Linear::Vector appdQdxdVp( *appVecPtr_ );
672 
673  // This is a temporary load storage vector.
674  Linear::Vector dQdt2( *appVecPtr_ );
675 
676  bXtPtr_->putScalar(0.0);
677 
678  Linear::BlockVector & bXf = *dynamic_cast<Linear::BlockVector*>(Xf);
679 
680  permutedIFT(bXf, &*bXtPtr_);
681 
682  // 12/8/06 tscoffe: Note: "b" at beginning of variable name means Linear::BlockVector
683  Linear::BlockVector & bX = *bXtPtr_;
684  Linear::BlockVector & bS = *dynamic_cast<Linear::BlockVector*>(S);
685  Linear::BlockVector & bcurrS = *dynamic_cast<Linear::BlockVector*>(currS);
686  Linear::BlockVector & blastS = *dynamic_cast<Linear::BlockVector*>(lastS);
687  Linear::BlockVector & bdSdt = *dynamic_cast<Linear::BlockVector*>(dSdt);
688  Linear::BlockVector & bStore = *dynamic_cast<Linear::BlockVector*>(Store);
689  Linear::BlockVector & bcurrStore = *dynamic_cast<Linear::BlockVector*>(currStore);
690  Linear::BlockVector & blastStore = *dynamic_cast<Linear::BlockVector*>(lastStore);
691 
692  Linear::BlockVector & bNextLeadF = *dynamic_cast<Linear::BlockVector*>(nextLeadFVectorPtr);
693  Linear::BlockVector & bCurrLeadF = *dynamic_cast<Linear::BlockVector*>(currLeadFVectorPtr);
694  Linear::BlockVector & bLastLeadF = *dynamic_cast<Linear::BlockVector*>(lastLeadFVectorPtr);
695  Linear::BlockVector & bLeadQ = *dynamic_cast<Linear::BlockVector*>(nextLeadQVectorPtr);
696  Linear::BlockVector & bNextJunctionV = *dynamic_cast<Linear::BlockVector*>(nextJunctionVVectorPtr);
697  Linear::BlockVector & bCurrJunctionV = *dynamic_cast<Linear::BlockVector*>(currJunctionVVectorPtr);
698  Linear::BlockVector & bLastJunctionV = *dynamic_cast<Linear::BlockVector*>(lastJunctionVVectorPtr);
699 
700  Linear::BlockVector * bQ = dynamic_cast<Linear::BlockVector*>(Q);
701  Linear::BlockVector * bF = dynamic_cast<Linear::BlockVector*>(F);
702  Linear::BlockVector * bB = dynamic_cast<Linear::BlockVector*>(B);
703 
704  Linear::BlockVector * bdFdxdVp = dynamic_cast<Linear::BlockVector*>(dFdxdVp);
705  Linear::BlockVector * bdQdxdVp = dynamic_cast<Linear::BlockVector*>(dQdxdVp);
706 
707 
708 // Linear::Vector * storeLeadCurrQ,
709 
710  Linear::BlockVector & bstoreLeadCurrQ = *dynamic_cast<Linear::BlockVector*>(storeLeadCurrQ);
711 
712  Linear::BlockVector bQt(*bXtPtr_);
713  Linear::BlockVector bFt(*bXtPtr_);
714  Linear::BlockVector bBt(*bXtPtr_);
715 
716  Linear::BlockVector bdQdxdVpt(*bXtPtr_);
717  Linear::BlockVector bdFdxdVpt(*bXtPtr_);
718 
719 
720 // Linear::BlockVector bstoreLeadCurrQt(*bXtPt );
721 
722 #ifndef Xyce_FLEXIBLE_DAE_LOADS
723  bmdQdxPtr_->put(0.0);
724  bmdFdxPtr_->put(0.0);
725 #endif
726 
727  int BlockCount = bX.blockCount();
728  int BlockSize = bX.blockSize();
729 
730 #if 0
731  Xyce::dout() << "bX.blockCount = "<< BlockCount <<std::endl;
732  Xyce::dout() << "bX.blockSize = "<< BlockSize <<std::endl;
733 #endif
734 
735  // We are storing the time domain Jacobians, initialize the memory here
736  Teuchos::RCP<Epetra_CrsGraph> linearGraph, nonlinGraph;
737  builder_.getSeparatedGraph(linearGraph, nonlinGraph);
738 
739  if ((int)vecNLAppdQdxPtr_.size()!=BlockCount)
740  {
741  // Ask the time-domain builder to created separated maps and graphs for
742  // linear and nonlinear devices
743  if (Teuchos::is_null(linearGraph) || Teuchos::is_null(nonlinGraph))
744  {
745  builder_.setupSeparatedLSObjects();
746  builder_.getSeparatedGraph(linearGraph, nonlinGraph);
747  }
748 
749  if (nonlinGraph->GlobalMaxNumIndices())
750  {
751  vecNLAppdQdxPtr_.resize( BlockCount );
752  vecNLAppdFdxPtr_.resize( BlockCount );
753  for ( int i=0; i < BlockCount; ++i )
754  {
755  Epetra_CrsMatrix* newdQdx_i = new Epetra_CrsMatrix( Copy, *nonlinGraph );
756  vecNLAppdQdxPtr_[i] = Teuchos::rcp( new Linear::Matrix( newdQdx_i, true ) );
757  Epetra_CrsMatrix* newdFdx_i = new Epetra_CrsMatrix( Copy, *nonlinGraph );
758  vecNLAppdFdxPtr_[i] = Teuchos::rcp( new Linear::Matrix( newdFdx_i, true ) );
759  }
760  }
761  if (linearGraph->GlobalMaxNumIndices())
762  {
763  Epetra_CrsMatrix* newdQdx = new Epetra_CrsMatrix( Copy, *linearGraph );
764  linAppdQdxPtr_ = Teuchos::rcp( new Linear::Matrix( newdQdx, true ) );
765  Epetra_CrsMatrix* newdFdx = new Epetra_CrsMatrix( Copy, *linearGraph );
766  linAppdFdxPtr_ = Teuchos::rcp( new Linear::Matrix( newdFdx, true ) );
767  }
768  }
769 
770  // Now perform implicit application of frequency domain Jacobian.
771  for( int i = 0; i < BlockCount; ++i )
772  {
773  if (DEBUG_HB)
774  {
775  Xyce::dout() << "Processing vectors for block " << i << " of " << BlockCount-1 << std::endl;
776  }
777 
778  //Set Time for fast time scale somewhere
779 // state_.fastTime = times_[i];
781 
782  if (DEBUG_HB)
783  {
784  Xyce::dout() << "Calling updateSources on the appLoader" << std::endl;
785  }
786 
787  //Update the sources
788  appLoaderPtr_->updateSources(); // this is here to handle "fast" sources.
789 
790  *appVecPtr_ = bX.block(i);
791  *appNextStaVecPtr_ = bS.block(i);
792  *appCurrStaVecPtr_ = bcurrS.block(i);
793  *appLastStaVecPtr_ = blastS.block(i);
794  appdSdt = bdSdt.block(i);
795  *appNextStoVecPtr_ = bStore.block(i);
796  *appCurrStoVecPtr_ = bcurrStore.block(i);
797  *appLastStoVecPtr_ = blastStore.block(i);
798  *appStoLeadCurrQVecPtr_ = bstoreLeadCurrQ.block(i);
799  *appNextLeadFVecPtr_ = bNextLeadF.block(i);
800  *appCurrLeadFVecPtr_ = bCurrLeadF.block(i);
801  *appLastLeadFVecPtr_ = bLastLeadF.block(i);
802  *appLeadQVecPtr_ = bLeadQ.block(i);
803  *appNextJunctionVVecPtr_ = bNextJunctionV.block(i);
804  *appCurrJunctionVVecPtr_ = bCurrJunctionV.block(i);
805  *appLastJunctionVVecPtr_ = bLastJunctionV.block(i);
806 
807 //blastStore.block(i); // Need to update to correct block!
808 
809  if (DEBUG_HB)
810  {
811  Xyce::dout() << "Updating State for block " << i << " of " << BlockCount-1 << std::endl;
812  }
813 
814  // Note: This updateState call is here (instead of in the
815  // HBLoader::updateState function) because it has to be called
816  // for the same fast time point.
817  appLoaderPtr_->updateState
818  ( &*appVecPtr_,
819  &*appVecPtr_, // note, this is a placeholder! ERK
820  &*appVecPtr_, // note, this is a placeholder! ERK
823  );
824 
825  bS.block(i) = *appNextStaVecPtr_;
826  bcurrS.block(i) = *appCurrStaVecPtr_;
827  blastS.block(i) = *appLastStaVecPtr_;
828  bStore.block(i) = *appNextStoVecPtr_;
829  bcurrStore.block(i) = *appCurrStoVecPtr_;
830  blastStore.block(i) = *appLastStoVecPtr_;
831 
832  if (DEBUG_HB)
833  {
834  Xyce::dout() << "Calling loadDAEVectors on the appLoader" << std::endl;
835  }
836 
837  // This has to be done because the app loader does NOT zero these vectors out.
838  appQ.putScalar(0.0);
839  appF.putScalar(0.0);
840  appB.putScalar(0.0);
841  appdFdxdVp.putScalar(0.0);
842  appdQdxdVp.putScalar(0.0);
843 
844  appLoaderPtr_->loadDAEVectors
845  ( &*appVecPtr_,
846  &*appVecPtr_, // note, this is a placeholder! ERK
847  &*appVecPtr_, // note, this is a placeholder! ERK
852  &appQ, &appF, &appB,
853  &appdFdxdVp, &appdQdxdVp );
854 
855  bQt.block(i) = appQ;
856  bFt.block(i) = appF;
857  bBt.block(i) = appB;
858 
859  bdQdxdVpt.block(i) = appdQdxdVp;
860  bdFdxdVpt.block(i) = appdFdxdVp;
861 
862  bstoreLeadCurrQ.block(i) = *appStoLeadCurrQVecPtr_;
863 
864  bStore.block(i) = *appNextStoVecPtr_; //lead current get loaded during loadDAEVectors
865 
866  bNextLeadF.block(i) = *appNextLeadFVecPtr_; // lead currents loaded into lead current vector.
867  bLeadQ.block(i) = *appLeadQVecPtr_;
868 
869  // Store the time domain Jacobian for future use.
870  appdQdxPtr_->put(0.0);
871  appdFdxPtr_->put(0.0);
872 
873  // Load dQdx and dFdx into the storage location for this time point.
874  appLoaderPtr_->loadDAEMatrices( &*appVecPtr_, &*appNextStaVecPtr_, &appdSdt, &*appNextStoVecPtr_, &*appdQdxPtr_, &*appdFdxPtr_);
875 
876  // Copy over matrix values into separate storage for linear and nonlinear.
877 
878  // dQdx
879  std::vector< RCP<Linear::Matrix> > outMatrices;
880  if (i==0 && !Teuchos::is_null(linAppdQdxPtr_) )
881  outMatrices.push_back( linAppdQdxPtr_ );
882  if (vecNLAppdQdxPtr_.size() != 0)
883  outMatrices.push_back( vecNLAppdQdxPtr_[i] );
884  Linear::extractValues( *appdQdxPtr_, outMatrices );
885 
886  outMatrices.clear();
887 
888  // dFdx
889  if (i==0 && !Teuchos::is_null(linAppdFdxPtr_))
890  outMatrices.push_back(linAppdFdxPtr_);
891  if (vecNLAppdFdxPtr_.size() != 0)
892  outMatrices.push_back(vecNLAppdFdxPtr_[i]);
893  Linear::extractValues( *appdFdxPtr_, outMatrices );
894  }
895 
896  permutedFFT(bBt, bB);
897 
898  permutedFFT(bQt, bQ);
899 
900  permutedFFT(bFt, bF);
901 
902  permutedFFT(bdQdxdVpt, bdQdxdVp);
903 
904  permutedFFT(bdFdxdVpt, bdFdxdVp);
905 
906  bStoreVecFreqPtr_->putScalar(0.0);
907  bStoreLeadCurrQVecFreqPtr_->putScalar(0.0);
908 
909  permutedFFT(bStore, &*bStoreVecFreqPtr_);
910  permutedFFT(bstoreLeadCurrQ, &*bStoreLeadCurrQVecFreqPtr_);
911 
912  bLeadCurrentVecFreqPtr_->putScalar(0.0);
913  bLeadCurrentQVecFreqPtr_->putScalar(0.0);
914  permutedFFT(bNextLeadF, &*bLeadCurrentVecFreqPtr_);
916 
917  int blockCount = bXf.blockCount();
918  int blockSize = bXf.block(0).globalLength();
919 // double omega = 2.0 * M_PI/ period_;
920 
921  int size_ = freqs_.size();
922  int posFreq = (size_-1)/2;
923  double omega = 2.0 * M_PI * freqs_[posFreq];
924 
925  for( int i = 0; i < blockCount; ++i )
926  {
927  // Create work vectors from the current frequency block vector
928  // NOTE: This needs to be done for each block to make sure that the
929  // map is the same as the bF block.
930  Linear::Vector QVec(bQ->block(i));
931  Linear::Vector freqVec = bQ->block(i);
932 
933  Linear::Vector dQdxdVpVec(bdQdxdVp->block(i));
934  Linear::Vector freqVec1 = bdQdxdVp->block(i);
935 
936  // Only one processor owns each block of the frequency-domain vector
937  if (freqVec.localLength() > 0)
938  {
939  omega = 2.0 * M_PI * freqs_[posFreq];
940 
941  QVec[0] = -freqVec[1]*omega;
942  QVec[1] = freqVec[0]*omega;
943 
944  dQdxdVpVec[0] = -freqVec1[1]*omega;
945  dQdxdVpVec[1] = freqVec1[0]*omega;
946 
947  for (int j=1; j < (blockSize/2+1)/2; ++j)
948  {
949 
950  omega = 2.0 * M_PI * freqs_[posFreq + j];
951 
952  QVec[2*j] = -freqVec[2*j+1]*omega;
953  QVec[2*(blockSize/2-j)] = -freqVec[2*j+1]*omega;
954 
955  QVec[2*j+1] = freqVec[2*j]*omega;
956  QVec[2*(blockSize/2-j)+1] = -freqVec[2*j]*omega;
957 
958  dQdxdVpVec[2*j] = -freqVec1[2*j+1]*omega;
959  dQdxdVpVec[2*(blockSize/2-j)] = -freqVec1[2*j+1]*omega;
960 
961  dQdxdVpVec[2*j+1] = freqVec1[2*j]*omega;
962  dQdxdVpVec[2*(blockSize/2-j)+1] = -freqVec1[2*j]*omega;
963 
964  }
965  }
966 
967  bF->block(i).update(1.0, QVec , 1.0);
968 
969  bdFdxdVp->block(i).update(1.0, dQdxdVpVec, 1.0);
970  }
971 
972  blockCount = bStoreVecFreqPtr_->blockCount();
973  blockSize = bStoreVecFreqPtr_->blockSize();
974 
975  for( int i = 0; i < blockCount; ++i )
976  {
977  // Create work vectors from the current frequency block vector
978  // NOTE: This needs to be done for each block to make sure that the
979  // map is the same as the bF block.
980  Linear::Vector stoLeadCurrdQdtVec(bStoreLeadCurrQVecFreqPtr_->block(i));
981  Linear::Vector stoLeadCurrQVec = bStoreLeadCurrQVecFreqPtr_->block(i);
982 
983  // Only one processor owns each block of the frequency-domain vector
984  if (stoLeadCurrQVec.localLength() > 0)
985  {
986 
987  omega = 2.0 * M_PI * freqs_[posFreq];
988  stoLeadCurrdQdtVec[0] = stoLeadCurrQVec[1]*omega;
989  stoLeadCurrdQdtVec[1] = stoLeadCurrQVec[0]*omega;
990 
991  for (int j=1; j < (blockSize/2+1)/2; ++j)
992  {
993  omega = 2.0 * M_PI * freqs_[posFreq+ j];
994 
995  stoLeadCurrdQdtVec[2*j] = -stoLeadCurrQVec[2*j+1]*omega;
996  stoLeadCurrdQdtVec[2*(blockSize/2-j)] = -stoLeadCurrQVec[2*j+1]*omega;
997 
998  stoLeadCurrdQdtVec[2*j+1] = stoLeadCurrQVec[2*j]*omega;
999  stoLeadCurrdQdtVec[2*(blockSize/2-j)+1] = -stoLeadCurrQVec[2*j]*omega;
1000  }
1001  }
1002 
1003  bStoreVecFreqPtr_->block(i).update(1.0, stoLeadCurrdQdtVec, 1.0);
1004  }
1005 
1006  //
1007  // take care of the lead current vector
1008  blockCount = bLeadCurrentVecFreqPtr_->blockCount();
1009  blockSize = bLeadCurrentVecFreqPtr_->blockSize();
1010 
1011  for( int i = 0; i < blockCount; ++i )
1012  {
1013  // Create work vectors from the current frequency block vector
1014  // NOTE: This needs to be done for each block to make sure that the
1015  // map is the same as the bF block.
1016  Linear::Vector leadCurrdQdtVec(bLeadCurrentQVecFreqPtr_->block(i));
1017  Linear::Vector leadCurrQVec = bLeadCurrentQVecFreqPtr_->block(i);
1018 
1019  // Only one processor owns each block of the frequency-domain vector
1020  if (leadCurrQVec.localLength() > 0)
1021  {
1022 
1023  omega = 2.0 * M_PI * freqs_[posFreq];
1024  leadCurrdQdtVec[0] = leadCurrQVec[1]*omega;
1025  leadCurrdQdtVec[1] = leadCurrQVec[0]*omega;
1026 
1027  for (int j=1; j < (blockSize/2+1)/2; ++j)
1028  {
1029  omega = 2.0 * M_PI * freqs_[posFreq+ j];
1030 
1031  leadCurrdQdtVec[2*j] = -leadCurrQVec[2*j+1]*omega;
1032  leadCurrdQdtVec[2*(blockSize/2-j)] = -leadCurrQVec[2*j+1]*omega;
1033 
1034  leadCurrdQdtVec[2*j+1] = leadCurrQVec[2*j]*omega;
1035  leadCurrdQdtVec[2*(blockSize/2-j)+1] = -leadCurrQVec[2*j]*omega;
1036  }
1037  }
1038 
1039  bLeadCurrentVecFreqPtr_->block(i).update(1.0, leadCurrdQdtVec, 1.0);
1040  }
1041 
1042 
1043 // permutedIFT(*bStoreVecFreqPtr_, &bStore);
1044  if (DEBUG_HB)
1045  {
1046  Xyce::dout() << "HB X Vector" << std::endl;
1047  bX.printPetraObject(std::cout);
1048  // Xyce::dout() << "HB S Vector" << std::endl;
1049  // bS.printPetraObject(Xyce::dout());
1050  // Xyce::dout() << "HB dSdt Vector" << std::endl;
1051  // bdSdt.printPetraObject(Xyce::dout());
1052  Xyce::dout() << "HB Store Vector" << std::endl;
1053  bStore.printPetraObject(std::cout);
1054  Xyce::dout() << "HB Q Vector" << std::endl;
1055  bQ->printPetraObject(std::cout);
1056  Xyce::dout() << "HB F Vector" << std::endl;
1057  bF->printPetraObject(std::cout);
1058  Xyce::dout() << "HB bdFdxdVp Vector" << std::endl;
1059  bdFdxdVp->printPetraObject(std::cout);
1060  Xyce::dout() << "HB bdQdxdVp Vector" << std::endl;
1061  bdQdxdVp->printPetraObject(std::cout);
1062 
1063 #ifndef Xyce_FLEXIBLE_DAE_LOADS
1064  Xyce::dout() << "HB bmdQdx_" << std::endl;
1065  bmdQdxPtr_->printPetraObject(std::cout);
1066  Xyce::dout() << "HB bmdFdx_" << std::endl;
1067  bmdFdxPtr_->printPetraObject(std::cout);
1068 #endif // Xyce_FLEXIBLE_DAE_LOADS
1069  Xyce::dout() << Xyce::section_divider << std::endl;
1070  }
1071 
1072  return true;
1073 }
1074 
1075 //-----------------------------------------------------------------------------
1076 // Function : HBLoader::loadDeviceErrorWeightMask
1077 // Purpose :
1078 // Special Notes :
1079 // Scope : public
1080 // Creator : Eric Keiter
1081 // Creation Date : 11/1/2014
1082 //-----------------------------------------------------------------------------
1083 bool
1085  Linear::Vector * deviceMask) const
1086 {
1087 #if 0
1088  Linear::BlockVector & bDevMask = *dynamic_cast<Linear::BlockVector*>(deviceMask);
1089 
1090 #if 0
1091  Xyce::dout() << "HBLoader::loadDeviceErrorWeightMask. Original (nonblock) deviceMask.size = "
1092  << appVecPtr_->globalLength() <<std::endl;
1093 #endif
1094 
1095  appVecPtr_->putScalar(1.0);
1096  bool returnValue = deviceManager_.loadErrorWeightMask(&*appVecPtr_);
1097 
1098  int blockCount = bDevMask.blockCount();
1099  int blockSize = bDevMask.blockSize();
1100 
1101 #if 0
1102  Xyce::dout() << "bDevMask.blockCount = "<< blockCount <<std::endl;
1103  Xyce::dout() << "bDevMask.blockSize = "<< blockSize <<std::endl;
1104  appVecPtr_->printPetraObject(Xyce::dout());
1105 #endif
1106 
1107  //Teuchos::RCP<N_PDS_ParMap> baseMap = Teuchos::rcp_const_cast<N_PDS_ParMap>( hbBuilderPtr_->getBaseStoreMap() );
1108 
1109  for( int i = 0; i < blockCount; ++i )
1110  {
1111  // See if this variable is owned by the local processor.
1112  // If so, this processor owns the entire j-th block of the vector
1113  //int lid = baseMap->globalToLocalIndex( i );
1114 
1115  Linear::Vector& localVecRef = bDevMask.block(i);
1116 
1117  for (int j=0;j<blockSize;++j)
1118  {
1119  localVecRef[j] = (*appVecPtr_)[i];
1120  }
1121  }
1122 
1123 #if 0
1124  bDevMask.printPetraObject(Xyce::dout());
1125 #endif
1126 
1127  return returnValue;
1128 #else
1129  return true;
1130 #endif
1131 }
1132 
1133 //-----------------------------------------------------------------------------
1134 // Function : HBLoader::permutedFFT
1135 // Purpose :
1136 // Special Notes :
1137 // Scope : public
1138 // Creator : Ting Mei
1139 // Creation Date : 09/05/08
1140 //---------------------------------------------------------------------------
1141 void HBLoader::permutedFFT(const Linear::BlockVector & xt, Linear::BlockVector * xf)
1142 {
1143  // Call the function to compute the permuted FFT from the block system helper functions.
1144  computePermutedDFT( *dftInterface_, xt, xf );
1145 }
1146 
1147 //-----------------------------------------------------------------------------
1148 // Function : HBLoader::permutedIFT
1149 // Purpose :
1150 // Special Notes :
1151 // Scope : public
1152 // Creator : Ting Mei
1153 // Creation Date : 09/05/08
1154 //---------------------------------------------------------------------------
1155 void HBLoader::permutedIFT(const Linear::BlockVector & xf, Linear::BlockVector * xt, int numTimePts_ )
1156 {
1157  // Call the function to compute the permuted IFT from the block system helper functions.
1158  computePermutedIFT( *dftInterface_, xf, xt, numTimePts_);
1159 }
1160 
1161 //-----------------------------------------------------------------------------
1162 // Function : HBLoader::getVoltageLimiterStatus
1163 // Purpose :
1164 // Special Notes :
1165 // Scope : public
1166 // Creator : Eric Keiter
1167 // Creation Date : 6/22/2015
1168 //---------------------------------------------------------------------------
1170 {
1171  return appLoaderPtr_->getVoltageLimiterStatus();
1172 }
1173 
1174 //-----------------------------------------------------------------------------
1175 // Function : HBLoader::setVoltageLimiterStatus
1176 // Purpose :
1177 // Special Notes :
1178 // Scope : public
1179 // Creator : Eric Keiter
1180 // Creation Date : 6/22/2015
1181 //---------------------------------------------------------------------------
1182 void HBLoader::setVoltageLimiterStatus(bool voltageLimterStatus)
1183 {
1184  return appLoaderPtr_->setVoltageLimiterStatus(voltageLimterStatus);
1185 }
1186 
1187 } // namespace Loader
1188 } // namespace Xyce
HBLoader(const Teuchos::RCP< const N_MPDE_Discretization > discPtr, Device::DeviceMgr &device_manager, Linear::Builder &builder)
Teuchos::RCP< Linear::Vector > appLastStaVecPtr_
virtual void setVoltageLimiterStatus(bool voltageLimterStatus)
Teuchos::RCP< Linear::Vector > appCurrStoVecPtr_
Teuchos::RCP< Linear::BlockMatrix > bmdFdxPtr_
Linear::Builder & builder_
std::vector< double > periodicTimes_
bool loadErrorWeightMask(Linear::Vector *deviceMaskPtr)
Pure virtual class to augment a linear system.
void setFastTimes(const std::vector< double > &times)
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)
void setFastTime(double timeVal)
Device::DeviceMgr & deviceManager_
void permutedIFT(const Linear::BlockVector &xf, Linear::BlockVector *xt, int numTimePts_=0)
bool loadDAEMatrices(Linear::Vector *X, Linear::Vector *S, Linear::Vector *dSdt, Linear::Vector *Store, Linear::Matrix *dQdx, Linear::Matrix *dFdx)
Teuchos::RCP< Loader > appLoaderPtr_
Actually a CktLoader.
Teuchos::RCP< const N_MPDE_Discretization > fastTimeDiscPtr_
Teuchos::RCP< Linear::Matrix > appdQdxPtr_
Teuchos::RCP< Linear::Vector > appStoLeadCurrQVecPtr_
void permutedFFT(const Linear::BlockVector &xt, Linear::BlockVector *xf)
Teuchos::RCP< Linear::Vector > appCurrStaVecPtr_
Teuchos::RCP< Linear::HBBuilder > hbBuilderPtr_
Teuchos::RCP< Linear::Matrix > linAppdFdxPtr_
Device::DeviceMgr & deviceManager_
Definition: N_ANP_HB.C:2199
Teuchos::RCP< Linear::Vector > appCurrLeadFVecPtr_
Teuchos::RCP< Linear::BlockVector > bXtPtr_
Teuchos::RCP< Linear::BlockVector > bLeadCurrentVecFreqPtr_
Teuchos::RCP< Linear::BlockVector > bStoreVecFreqPtr_
Teuchos::RCP< Linear::BlockVector > bVtPtr_
Teuchos::RCP< Linear::Vector > appNextJunctionVVecPtr_
Teuchos::RCP< Linear::BlockVector > bStoreLeadCurrQVecFreqPtr_
Teuchos::RCP< Linear::BlockMatrix > bmdQdxPtr_
bool applyLinearMatrices(const Linear::BlockVector &bVf, Linear::BlockVector &permlindQdxV, Linear::BlockVector &permlindFdxV)
void setHBFreqs(const std::vector< double > &freqs)
void registerHBBuilder(Teuchos::RCP< Linear::HBBuilder > hbBuilderPtr)
virtual bool getVoltageLimiterStatus()
Teuchos::RCP< N_UTL_DFTInterfaceDecl< std::vector< double > > > dftInterface_
Teuchos::RCP< Linear::Vector > appLeadQVecPtr_
std::vector< Teuchos::RCP< Linear::Matrix > > vecNLAppdQdxPtr_
Teuchos::RCP< Linear::Vector > appLastLeadFVecPtr_
bool loadTimeDepDAEMatrices(Linear::Vector *X, Linear::Vector *S, Linear::Vector *dSdt, Linear::Vector *Store, Linear::Matrix *dQdx, Linear::Matrix *dFdx)
Teuchos::RCP< Linear::Matrix > appdFdxPtr_
bool applyDAEMatrices(Linear::Vector *X, Linear::Vector *S, Linear::Vector *dSdt, Linear::Vector *Store, const Linear::Vector &V, Linear::Vector *dQdxV, Linear::Vector *dFdxV)
std::vector< double > times_
Linear::Builder & builder_
Definition: N_ANP_HB.C:2200
#define M_PI
Teuchos::RCP< Linear::Vector > appLastStoVecPtr_
std::vector< Teuchos::RCP< Linear::Matrix > > vecNLAppdFdxPtr_
Teuchos::RCP< Linear::BlockVector > bLeadCurrentQVecFreqPtr_
Teuchos::RCP< Linear::Vector > appNextLeadFVecPtr_
Teuchos::RCP< Linear::Vector > appLastJunctionVVecPtr_
std::vector< double > freqs_
Teuchos::RCP< Linear::Vector > appVecPtr_
Teuchos::RCP< Linear::Matrix > linAppdQdxPtr_
bool loadDAEVectors(Linear::Vector *X, Linear::Vector *currX, Linear::Vector *lastX, Linear::Vector *S, Linear::Vector *currS, Linear::Vector *lastS, Linear::Vector *dSdt, Linear::Vector *Store, Linear::Vector *currStore, Linear::Vector *lastStore, Linear::Vector *storeLeadCurrQ, Linear::Vector *nextLeadFVectorPtr, Linear::Vector *currLeadFVectorPtr, Linear::Vector *lastLeadFVectorPtr, Linear::Vector *nextLeadQVectorPtr, Linear::Vector *nextJunctionVVectorPtr, Linear::Vector *currJunctionVVectorPtr, Linear::Vector *lastJunctionVVectorPtr, Linear::Vector *Q, Linear::Vector *F, Linear::Vector *B, Linear::Vector *dFdxdVp, Linear::Vector *dQdxdVp)
Teuchos::RCP< Linear::Vector > appCurrJunctionVVecPtr_
Teuchos::RCP< Linear::Vector > appNextStoVecPtr_
Teuchos::RCP< Linear::Vector > appNextStaVecPtr_
bool loadDeviceErrorWeightMask(Linear::Vector *deviceMask) const