Xyce  6.1
N_DEV_MutIndNonLin2.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_DEV_MutIndNonLin2.C,v $
27 //
28 // Purpose :
29 //
30 // Special Notes :
31 //
32 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
33 //
34 // Creation Date : 03/21/2005
35 //
36 // Revision Information:
37 // ---------------------
38 //
39 // Revision Number: $Revision: 1.63.2.1 $
40 //
41 // Revision Date : $Date: 2015/04/02 18:20:11 $
42 //
43 // Current Owner : $Author: tvrusso $
44 //-------------------------------------------------------------------------
45 
46 #include <Xyce_config.h>
47 
48 // #define Xyce_NO_MUTIND_MASK 1
49 // ---------- Standard Includes ----------
50 #include <fstream>
51 #include <algorithm>
52 #include <vector>
53 #include <set>
54 
55 // ---------- Xyce Includes ----------
56 #include <N_DEV_Const.h>
57 #include <N_DEV_DeviceOptions.h>
58 #include <N_DEV_DeviceMaster.h>
59 #include <N_DEV_ExternData.h>
60 #include <N_DEV_MatrixLoadData.h>
61 #include <N_DEV_MutIndNonLin2.h>
62 #include <N_DEV_SolverState.h>
63 #include <N_DEV_Message.h>
64 #include <N_ERH_ErrorMgr.h>
65 
66 #include <N_DEV_MutIndNonLin.h>
67 
68 #include <N_LAS_Vector.h>
69 #include <N_LAS_Matrix.h>
70 #include <N_UTL_FeatureTest.h>
71 #include <N_UTL_Math.h>
72 
73 using Teuchos::rcp;
74 
75 namespace Xyce {
76 namespace Device {
77 
78 
79 namespace MutIndNonLin2 {
80 
81 
83 {
84  p.addPar ("COUP_VAL",1.0,&MutIndNonLin2::Instance::mutualCup)
86  .setUnit(U_NONE)
87  .setCategory(CAT_NONE)
88  .setDescription("Coupling coefficient");
89 
90  p.addPar ("NONLINEARCOUPLING",0.0,&MutIndNonLin2::Instance::nonlinFlag)
92  .setUnit(U_NONE)
93  .setCategory(CAT_NONE)
94  .setDescription("Nonlinear coupling flag");
95 
96  p.addPar ("COUPLEDMutIndNonLin",std::vector<std::string>(),&MutIndNonLin2::Instance::inductorNames)
97  .setUnit(U_NONE)
98  .setCategory(CAT_NONE)
99  .setDescription("" );
100 
101  p.addPar ("COUPLEDINDUCTANCE",std::vector<double>(),&MutIndNonLin2::Instance::inductorInductances)
102  .setUnit(U_NONE)
103  .setCategory(CAT_NONE)
104  .setDescription("");
105 
106  p.addPar ("NODE1",std::vector<std::string>(),&MutIndNonLin2::Instance::inductorsNode1)
107  .setUnit(U_NONE)
108  .setCategory(CAT_NONE)
109  .setDescription("");
110 
111  p.addPar ("NODE2",std::vector<std::string>(),&MutIndNonLin2::Instance::inductorsNode2)
112  .setUnit(U_NONE)
113  .setCategory(CAT_NONE)
114  .setDescription("");
115 
116  p.addPar ("COUPLING",std::vector<double>(),&MutIndNonLin2::Instance::couplingCoefficient)
117  .setUnit(U_NONE)
118  .setCategory(CAT_NONE)
119  .setDescription("Coupling coefficient");
120 
121  p.addPar ("COUPLEDINDUCTOR",std::vector<std::string>(),&MutIndNonLin2::Instance::couplingInductor)
122  .setUnit(U_NONE)
123  .setCategory(CAT_NONE)
124  .setDescription("");
125 }
126 
128 {
129  p.addPar ("A",1000.0,&MutIndNonLin2::Model::A)
130  .setUnit(U_AMPMM1)
131  .setCategory(CAT_MATERIAL)
132  .setDescription("Thermal energy parameter");
133 
134  p.addPar ("AREA",0.1,&MutIndNonLin2::Model::Area)
135  .setUnit(U_CM2)
136  .setCategory(CAT_GEOMETRY)
137  .setDescription("Mean magnetic cross-sectional area");
138 
139  p.addPar ("ALPHA",5.0e-5,&MutIndNonLin2::Model::Alpha)
140  .setUnit(U_NONE)
141  .setCategory(CAT_GEOMETRY)
142  .setDescription("Domain coupling parameter");
143 
144  p.addPar ("BETAH",0.0001,&MutIndNonLin2::Model::BetaH)
145  .setUnit(U_NONE)
146  .setCategory(CAT_NONE)
147  .setDescription("Modeling constant");
148 
149  p.addPar ("BETAM",3.125e-5,&MutIndNonLin2::Model::BetaM)
150  .setUnit(U_NONE)
151  .setCategory(CAT_NONE)
152  .setDescription("Modeling constant");
153 
154  p.addPar ("C",0.2,&MutIndNonLin2::Model::C)
155  .setUnit(U_NONE)
156  .setCategory(CAT_MATERIAL)
157  .setDescription("Domain flesing parameter");
158 
159  p.addPar ("DELV",0.1,&MutIndNonLin2::Model::DeltaV)
160  .setUnit(U_VOLT)
161  .setCategory(CAT_NONE)
162  .setDescription("Smoothing coefficient for voltage difference over first inductor");
163 
164  p.addPar ("GAP",0.0,&MutIndNonLin2::Model::Gap)
165  .setUnit(U_CM)
166  .setCategory(CAT_GEOMETRY)
167  .setDescription("Effective air gap");
168 
169  p.addPar ("K",500.0,&MutIndNonLin2::Model::Kirr)
170  .setUnit(U_AMPMM1)
171  .setCategory(CAT_MATERIAL)
172  .setDescription("Domain anisotropy parameter");
173 
174  p.addPar ("KIRR",500.0,&MutIndNonLin2::Model::Kirr)
175  .setUnit(U_AMPMM1)
176  .setCategory(CAT_MATERIAL)
177  .setDescription("Domain anisotropy parameter");
178 
179  p.addPar ("MS",1.0e+6,&MutIndNonLin2::Model::Ms)
180  .setUnit(U_AMPMM1)
181  .setCategory(CAT_MATERIAL)
182  .setDescription("Saturation magnetization");
183 
185  .setUnit(U_NONE)
186  .setCategory(CAT_NONE)
187  .setDescription("for pspice compatibility -- ignored");
188 
190  .setUnit(U_NONE)
191  .setCategory(CAT_NONE)
192  .setDescription("for pspice compatibility -- ignored");
193 
194  p.addPar ("PATH",1.0,&MutIndNonLin2::Model::Path)
195  .setUnit(U_CM)
196  .setCategory(CAT_GEOMETRY)
197  .setDescription("Total mean magnetic path");
198 
199  p.addPar ("VINF",1.0,&MutIndNonLin2::Model::Vinf)
200  .setUnit(U_VOLT)
201  .setCategory(CAT_NONE)
202  .setDescription("Smoothing coefficient for voltage difference over first inductor");
203 
204  p.addPar ("TNOM",27.0,&MutIndNonLin2::Model::tnom)
205  .setUnit(U_DEGC)
206  .setCategory(CAT_MATERIAL)
207  .setDescription("Reference temperature");
208 
210  .setUnit(U_NONE)
211  .setCategory(CAT_MATERIAL)
212  .setDescription("First order temperature coeff.");
213 
215  .setUnit(U_NONE)
216  .setCategory(CAT_MATERIAL)
217  .setDescription("Second order temperature coeff.");
218 
219  p.addPar ("PZEROTOL",0.1,&MutIndNonLin2::Model::pZeroTol)
220  .setUnit(U_NONE)
221  .setCategory(CAT_NONE)
222  .setDescription("Tolerance for nonlinear zero crossing");
223 
224  p.addPar ("MVARSCALING",1.0,&MutIndNonLin2::Model::mVarScaling)
225  .setUnit(U_NONE)
226  .setCategory(CAT_NONE)
227  .setDescription("M-variable scaling.");
228 
229  p.addPar ("RVARSCALING",1.0,&MutIndNonLin2::Model::rVarScaling)
230  .setUnit(U_NONE)
231  .setCategory(CAT_NONE)
232  .setDescription("R-variable scaling");
233 
234  p.addPar ("MEQNSCALING",1.0,&MutIndNonLin2::Model::mEqScaling)
235  .setUnit(U_NONE)
236  .setCategory(CAT_NONE)
237  .setDescription("M-equation scaling");
238 
239  p.addPar ("REQNSCALING",1.0,&MutIndNonLin2::Model::rEqScaling)
240  .setUnit(U_NONE)
241  .setCategory(CAT_NONE)
242  .setDescription("R-equation scaling");
243 
244  p.addPar ("OUTPUTSTATEVARS",0,&MutIndNonLin2::Model::outputStateVars)
245  .setUnit(U_NONE)
246  .setCategory(CAT_NONE)
247  .setDescription("Flag to save state variables" );
248 
249  p.addPar ("INCLUDEDELTAM",0,&MutIndNonLin2::Model::includeDeltaM)
251  .setUnit(U_NONE)
252  .setCategory(CAT_NONE)
253  .setDescription("Flag to make M calculation implicit" );
254 
255  p.addPar ("USERKINTEGRATION",0,&MutIndNonLin2::Model::useRKIntegration)
257  .setUnit(U_NONE)
258  .setCategory(CAT_NONE)
259  .setDescription("Flag to use 4th order Runge-Kutta integration for dM/dH" );
260 
261  p.addPar ("USESTATEDERIV",0,&MutIndNonLin2::Model::useStateDeriv)
262  .setUnit(U_NONE)
263  .setCategory(CAT_NONE)
264  .setDescription("Flag to use state vector for derivatives" );
265 
266  p.addPar ("VOLTAGELIMITERFLAG",0,&MutIndNonLin2::Model::voltageLimiterFlag)
267  .setUnit(U_NONE)
268  .setCategory(CAT_NONE)
269  .setDescription("Flag to use voltage limiting on Mag and R internal variables" );
270 
271  p.addPar ("MAGLIMITTHRES",0.1,&MutIndNonLin2::Model::magLimitThres)
272  .setUnit(U_NONE)
273  .setCategory(CAT_NONE)
274  .setDescription("Threshold over which newton interation changes in Mag are limited." );
275 
276  p.addPar ("RLIMITTHRES",0.1,&MutIndNonLin2::Model::rLimitThres)
277  .setUnit(U_NONE)
278  .setCategory(CAT_NONE)
279  .setDescription("Threshold over which newton interation changes in R are limited." );
280 
281  p.addPar ("FACTORMS",0,&MutIndNonLin2::Model::factorMS)
282  .setUnit(U_NONE)
283  .setCategory(CAT_NONE)
284  .setDescription("Flag to save state variables" );
285 }
286 
287 
288 //-----------------------------------------------------------------------------
289 // Function : Instance::Instance
290 // Purpose : constructor
291 // Special Notes :
292 // Scope : public
293 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
294 // Creation Date : 03/21/2005
295 //-----------------------------------------------------------------------------
297  const Configuration & configuration,
298  const InstanceBlock & IB,
299  Model & Iiter,
300  const FactoryBlock & factory_block)
301  : DeviceInstance(IB, configuration.getInstanceParameters(), factory_block),
302  model_(Iiter),
303  temp(getDeviceOptions().temp.getImmutableValue<double>()),
304  P(0.0),
305  dP_dM(0.0),
306  dP_dBranchCurrentSum(0.0),
307  dP_dV1Pos(0.0),
308  dP_dV1Neg(0.0),
309  branchCurrentSum(0.0),
310  mEquFval(0.0),
311  MagVar(0.0),
312  oldBranchCurrentSum(0.0),
313  MagVarUpdate(0.0),
314  lastMagUpdate(0.0),
315  PPreviousStep(0.0),
316  includeDeltaM(false),
317  useRKIntegration(false),
318  outputStateVarsFlag( false )
319 {
320  if (DEBUG_DEVICE)
321  {
322  Xyce::dout() << "In Instance constructor" << std::endl;
323  }
324 
325  // for a simple case of 2 leads, we have 3 internal vars (I_branch, H, M)
326  // and one state var (I_branch)
327  numExtVars = 2;
328  numIntVars = 3;
329  numStateVars = 0;
330  tempGiven = false;
331 
332  const int ibev = IB.numExtVars;
333  const int ibiv = IB.numIntVars;
334 
335  // Set params to constant default values:
336  setDefaultParams ();
337 
338  // Set params according to instance line and constant defaults from metadata:
339  setParams (IB.params);
340 
341  // these two vectors are used for 4th order RK
342  // the vectors are one shorted than one would think you would need because
343  // current step is held by local variables (thus, this just data from steps
344  // n-1, n-2 and n-3
345  branchCurrentSumHistory.resize(3);
346  PFunctionHistory.resize(3);
347 
348  // if the model card askes for the delta M calculation to be implicit, then
349  // change the includeDeltaM flag
351  {
352  includeDeltaM = true;
353  }
354 
356  {
357  useRKIntegration = true;
358  }
359 
360  // now load the instance data vector
361  for( int i=0; i<inductorNames.size(); ++i )
362  {
363  InductorInstanceData * inductorData = new InductorInstanceData();
364  inductorData->name = inductorNames[i];
365  inductorData->L = inductorInductances[i];
366  inductorData->baseL = inductorInductances[i];
367  inductorData->ICGiven = false;
368  inductorData->inductorCurrentOffsets.resize( inductorNames.size() );
369 
370  instanceData.push_back( inductorData );
371  }
372  numInductors = instanceData.size();
373 
374  inductorCurrents.resize( numInductors );
375  inductanceVals.resize( numInductors );
376  LOI.resize( numInductors );
377  LO.resize( numInductors );
378  for( int i=0; i<numInductors; ++i)
379  {
380  LO[i].resize( numInductors );
381  }
382  dHe_dI.resize(numInductors);
383  dManp_dI.resize(numInductors);
384  ddelM_dI.resize(numInductors);
385  dMirrp_dI.resize(numInductors);
386  dP_dI.resize( numInductors );
387 
389 
390  // Calculate any parameters specified as expressions:
392 
393  // calculate dependent (ie computed) params and check for errors:
394  processParams ();
395 
396  // if the user has requested output of the internal vars H and M
397  // then open a file for that output.
398  if( model_.outputStateVars > 0 )
399  {
400  outputStateVarsFlag = true;
401  std::string filename( "Inductor_" );
402  filename.append( getName().getEncodedName() );
403  filename.append( ".dat" );
404  // convert any embeded ':' or '%' characters to '_'
405  replace( filename.begin(), filename.end(), '%', '_' );
406  replace( filename.begin(), filename.end(), ':', '_' );
407 
408  outputFileStreamPtr = rcp( new std::ofstream() );
409  outputFileStreamPtr->open( filename.c_str() );
410  if( !(*outputFileStreamPtr) )
411  {
412  UserError(*this) << "Could not open file " << filename << " for output of state variables";
413  }
414  else {
415  (*outputFileStreamPtr).setf(std::ios::scientific, std::ios::floatfield );
416  (*outputFileStreamPtr).width(20);
417  (*outputFileStreamPtr).precision(12);
418  }
419  }
420 
421  // update internal/external/state variable counts
422  numExtVars = 2*numInductors; // 2 nodes Vin, Vout per inductor
423  numIntVars = numInductors; // branch current per inductor
424  // if we're including the deltaM and deltaHapp variables, then there will be two more internal vars
425  if(includeDeltaM)
426  {
427  numIntVars+=1;
428  }
429 
430  // allocate space for InductorOffsets
431  deltaMEquInductorOffsets.resize(numInductors);
432 
433  // set up the jacobian stamp
434  // for an individual inductor with the two interal variables would be
435  //
436  // V1 V2 Ib
437  // kcl1 1
438  // kcl2 -1
439  // branch 1 -1 L/dt
440  //
441  // for a collection of these, the internal variable, branch equations,
442  // must be at the end of a given stamp row as well as the internal
443  // vars for M and R in this non-linear version.
444  //
445  // So for N inductors the samp is:
446  //
447  // V1 V2 V3 V4 ... V2N I1 I2 ... IN M
448  // kcl1 1
449  // kcl2 -1
450  // kcl3 1
451  // kcl4 -1
452  // branch1 1 -1 L/dt c ... c
453  // branch2 1 -1 c L/dt ... c
454  // delta M x x ... x x
455  //
456  // where "c" is an induced current change and "x" are
457  // values which must be computed.
458 
459  jacStamp.resize( numExtVars + numIntVars);
460 
461  for( int i=0; i< numInductors; ++i )
462  {
463  //
464  // allocate space
465  //
466  // kcl V+ node
467  jacStamp[2*i].resize(1);
468  // kcl V- node
469  jacStamp[2*i+1].resize(1);
470  if( i == 0 )
471  {
472  jacStamp[2*numInductors].resize(numInductors + 2);
473  }
474  else
475  {
476  jacStamp[2*numInductors + i].resize(numInductors + 2);
477  }
478 
479  //
480  // fill in dependency
481  //
482  // kcl V+ node
483  jacStamp[2*i ][0] = 2*numInductors + i;
484  // kcl V- node
485  jacStamp[2*i+1][0] = 2*numInductors + i;
486 
487  if( i==0 )
488  {
489  jacStamp[2*numInductors ][0] = 0;
490  jacStamp[2*numInductors ][1] = 1;
491  for( int j=0; j<numInductors; ++j )
492  {
493  jacStamp[2*numInductors][j+2] = 2*numInductors + j;
494  }
495  }
496  else
497  {
498  jacStamp[2*numInductors + i][0] = 2*i;
499  jacStamp[2*numInductors + i][1] = 2*i + 1;
500  for( int j=0; j<numInductors; ++j )
501  {
502  jacStamp[2*numInductors + i][j+2] = 2*numInductors + j;
503  }
504  }
505  }
506 
507  if(includeDeltaM)
508  {
509  // now the deltaM equation
510  jacStamp[ 3*numInductors ].resize(numInductors + 1);
511  // deltaM offsets to I_1 ... I_n and deltaM
512  for(int i=0; i<numInductors; ++i)
513  {
514  jacStamp[ 3*numInductors ][i]=2*numInductors+i;
515  }
517 
518  }
519  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS))
520  {
521  Xyce::dout() << "Instance::Instance----------" << std::endl;
522  Xyce::dout() << "numExtVars = " << numExtVars << ", " << ibev << std::endl
523  << "numIntVars = " << numIntVars << ", " << ibiv << std::endl
524  << "numStateVars = " << numStateVars << std::endl
525  << "numInductors = " << numInductors << std::endl
526  << "jacStamp = " << std::endl;
527  for( int i = 0; i<jacStamp.size(); ++i )
528  {
529  Xyce::dout() << "jacStamp[ " << i << " ] = { ";
530  for( int j=0; j<jacStamp[i].size(); ++j)
531  {
532  Xyce::dout() << jacStamp[i][j];
533  if( j != ( jacStamp[i].size() -1 ) )
534  {
535  Xyce::dout() << ", ";
536  }
537  }
538  Xyce::dout() << " }" << std::endl;
539  }
540  }
541 }
542 
543 
544 //-----------------------------------------------------------------------------
545 // Function : Instance::~Instance
546 // Purpose : destructor
547 // Special Notes :
548 // Scope : public
549 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
550 // Creation Date : 03/21/2005
551 //-----------------------------------------------------------------------------
553 {
554  // Close output file if we opened one
555  if( outputStateVarsFlag && outputFileStreamPtr.get() && outputFileStreamPtr->is_open() )
556  {
557  outputFileStreamPtr->close();
558  }
559 
560  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
561  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
562  for ( ; currentInductor != endInductor ; ++currentInductor)
563  {
564  if (*currentInductor != NULL)
565  {
566  delete *currentInductor;
567  *currentInductor = NULL;
568  }
569  }
570  instanceData.clear();
571 }
572 
573 //-----------------------------------------------------------------------------
574 // Function : Instance::registerLIDs
575 // Purpose :
576 // Special Notes :
577 // Scope : public
578 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
579 // Creation Date : 03/21/2005
580 //-----------------------------------------------------------------------------
581 void Instance::registerLIDs(const std::vector<int> & intLIDVecRef,
582  const std::vector<int> & extLIDVecRef)
583 {
584  AssertLIDs(intLIDVecRef.size() == numIntVars);
585  AssertLIDs(extLIDVecRef.size() == numExtVars);
586 
587  // copy over the global ID lists.
588  intLIDVec = intLIDVecRef;
589  extLIDVec = extLIDVecRef;
590 
591  // Now use these lists to obtain the indices into the
592  // linear algebra entities. This assumes an order.
593  // For the matrix indices, first do the rows.
594  // get the current values of the inductances and currentOffsets
595  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
596  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
597  int i = 0;
598  int j = 0;
599  while( currentInductor != endInductor )
600  {
601  (*currentInductor)->li_Pos = extLIDVec[ i++ ];
602  (*currentInductor)->li_Neg = extLIDVec[ i++ ];
603  (*currentInductor)->li_Branch = intLIDVec[ j++ ];
604  currentInductor++;
605  }
606 
607  if(includeDeltaM)
608  {
609  // now get deltaHapp and deltaM
610  //li_deltaHappVar = intLIDVec[ j++ ];
611  li_deltaMagVar = intLIDVec[ j++ ];
612  }
613 
614  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS))
615  {
616  Xyce::dout() << "Instance::registerLIDs------------------------" << std::endl;
617  currentInductor = instanceData.begin();
618  i=0;
619  while( currentInductor != endInductor )
620  {
621  Xyce::dout() << "Inductor [ " << i++ << " ] "
622  << " li_Pos = " << (*currentInductor)->li_Pos
623  << " li_Neg = " << (*currentInductor)->li_Neg
624  << " li_Branch = " << (*currentInductor)->li_Branch << std::endl;
625  currentInductor++;
626  }
627  if(includeDeltaM)
628  {
629  Xyce::dout() << " li_deltaMagVar = " << li_deltaMagVar << std::endl;
630  }
631  }
632 }
633 
634 //-----------------------------------------------------------------------------
635 // Function : Instance::loadNodeSymbols
636 // Purpose :
637 // Special Notes :
638 // Scope : public
639 // Creator : Eric R. Keiter, SNL, Parallel Computational Sciences
640 // Creation Date : 05/13/05
641 //-----------------------------------------------------------------------------
642 void Instance::loadNodeSymbols(Util::SymbolTable &symbol_table) const
643 {
644  for (std::vector<InductorInstanceData *>::const_iterator it = instanceData.begin(), end = instanceData.end(); it != end; ++it)
645  addInternalNode(symbol_table, (*it)->li_Branch, getName(), (*it)->name + "_branch");
646 }
647 
648 //-----------------------------------------------------------------------------
649 // Function : Instance::registerStateLIDs
650 // Purpose :
651 // Special Notes :
652 // Scope : public
653 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
654 // Creation Date : 03/21/2005
655 //-----------------------------------------------------------------------------
656 void Instance::registerStateLIDs( const std::vector<int> & staLIDVecRef )
657 {
658  AssertLIDs(staLIDVecRef.size() == numStateVars);
659 
660  // copy over the global ID lists.
661  staLIDVec = staLIDVecRef;
662 
663  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS))
664  {
665  Xyce::dout() << "Instance::registerStateLIDs-------------------" << std::endl;
666  }
667 }
668 
669 //-----------------------------------------------------------------------------
670 // Function : Instance::jacobianStamp
671 // Purpose :
672 // Special Notes :
673 // Scope : public
674 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
675 // Creation Date : 03/21/2005
676 //-----------------------------------------------------------------------------
677 const std::vector< std::vector<int> > & Instance::jacobianStamp() const
678 {
679  return jacStamp;
680 }
681 
682 //-----------------------------------------------------------------------------
683 // Function : Instance::registerJacLIDs
684 // Purpose :
685 // Special Notes :
686 // Scope : public
687 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
688 // Creation Date : 03/21/2005
689 //-----------------------------------------------------------------------------
690 void Instance::registerJacLIDs( const std::vector< std::vector<int> > & jacLIDVec )
691 {
692  DeviceInstance::registerJacLIDs( jacLIDVec );
693  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS))
694  {
695  Xyce::dout() << "Instance::registerJacLIDs ----------------------------" << std::endl;
696 
697  Xyce::dout() << "jacLIDVec = " << std::endl;
698  for( int i = 0; i<jacStamp.size(); ++i )
699  {
700  Xyce::dout() << "jacLIDVec[ " << i << " ] = { ";
701  for( int j=0; j<jacLIDVec[i].size(); ++j)
702  {
703  Xyce::dout() << jacLIDVec[i][j];
704  if( j != ( jacLIDVec[i].size() -1 ) )
705  {
706  Xyce::dout() << ", ";
707  }
708  }
709  Xyce::dout() << " }" << std::endl;
710  }
711  }
712 
713  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
714  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
715  int i = 0;
716  while( currentInductor != endInductor )
717  {
718  (*currentInductor)->APosEquBraVarOffset = jacLIDVec[ 2*i ][ 0 ];
719  (*currentInductor)->ANegEquBraVarOffset = jacLIDVec[ 2*i + 1 ][ 0 ];
720  (*currentInductor)->vPosOffset = jacLIDVec[ 2*numInductors + i ][ 0 ];
721  (*currentInductor)->vNegOffset = jacLIDVec[ 2*numInductors + i ][ 1 ];
722 
723  (*currentInductor)->ABraEquPosNodeOffset = jacLIDVec[ 2*numInductors + i ][ 0 ];
724  (*currentInductor)->ABraEquNegNodeOffset = jacLIDVec[ 2*numInductors + i ][ 1 ];
725  for( int j=0; j<numInductors; ++j )
726  {
727  if( i == j )
728  {
729  (*currentInductor)->ABraEquBraVarOffset = jacLIDVec[ 2*numInductors + i ][ j + 2 ];
730  }
731  (*currentInductor)->inductorCurrentOffsets[ j ] = jacLIDVec[ 2*numInductors + i ][ j + 2 ];
732  }
733 
734  currentInductor++;
735  i++;
736  }
737 
738  if(includeDeltaM)
739  {
740  // now get the deltaM
741  for( int i=0; i<numInductors; i++)
742  {
743  deltaMEquInductorOffsets[i] = jacLIDVec[ 3*numInductors ][ i ];
744  }
745  deltaMEquDeltaMOffset = jacLIDVec[ 3*numInductors ][ numInductors ];
746  }
747 
748  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS))
749  {
750  currentInductor = instanceData.begin();
751  i=0;
752  while( currentInductor != endInductor )
753  {
754  Xyce::dout() << "Inductor [ " << i << " ] " << (*currentInductor)->name << std::endl
755  << " APosEquBraVarOffset = " << (*currentInductor)->APosEquBraVarOffset << std::endl
756  << " ANegEquBraVarOffset = " << (*currentInductor)->ANegEquBraVarOffset << std::endl
757  << " vPosOffset = " << (*currentInductor)->vPosOffset << std::endl
758  << " vNegOffset = " << (*currentInductor)->vNegOffset << std::endl
759  << " ABraEquPosNodeOffset = " << (*currentInductor)->ABraEquPosNodeOffset << std::endl
760  << " ABraEquNegNodeOffset = " << (*currentInductor)->ABraEquNegNodeOffset << std::endl
761  << " ABraEquBraVarOffset = " << (*currentInductor)->ABraEquBraVarOffset << std::endl
762  << " magOffset = " << (*currentInductor)->magOffset << std::endl;
763  Xyce::dout() << "\tInductor branch offsets = { ";
764  for( int j=0; j<numInductors ; ++j )
765  {
766  Xyce::dout() << (*currentInductor)->inductorCurrentOffsets[ j ] << ", ";
767  }
768  Xyce::dout() << "} " << std::endl;
769  i++;
770  currentInductor++;
771  }
772  Xyce::dout() << std::endl;
773 
774  if(includeDeltaM)
775  {
776  Xyce::dout() << "deltaMEquInductorOffsets = ";
777  for( int i=0; i<numInductors; i++ )
778  {
779  Xyce::dout() << deltaMEquInductorOffsets[i] << " ";
780  }
781  Xyce::dout() //<< "deltaMEquDeltaHappOffset = " << deltaMEquDeltaHappOffset
782  << " deltaMEquDeltaMOffset = " << deltaMEquDeltaMOffset << std::endl;
783  }
784  }
785 }
786 
787 
788 //-----------------------------------------------------------------------------
789 // Function : Instance::processParams
790 // Purpose :
791 // Special Notes :
792 // Scope : public
793 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
794 // Creation Date : 03/21/2005
795 //-----------------------------------------------------------------------------
797 {
798  // now set the temperature related stuff.
799  if (tempGiven)
800  {
802  }
803 
804  return true;
805 }
806 
807 //-----------------------------------------------------------------------------
808 // Function : Instance::updateTemperature
809 // Purpose :
810 // Special Notes :
811 // Scope : public
812 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
813 // Creation Date : 03/21/2005
814 //-----------------------------------------------------------------------------
815 bool Instance::updateTemperature ( const double & temp)
816 {
817  bool bsuccess = true;
818 
819  // current temp difference from reference temp.
820  double difference = temp - model_.tnom;
821 
822  std::vector< InductorInstanceData* >::iterator currentData = instanceData.begin();
823  while( currentData != instanceData.end() )
824  {
825  double factor = 1.0 + (model_.tempCoeff1)*difference +
826  (model_.tempCoeff2)*difference*difference;
827  (*currentData)->L = ((*currentData)->baseL)*factor;
828  currentData++;
829  }
830 
831  // now that the inductances have changed we need to update the matrix.
833 
834  return bsuccess;
835 }
836 
837 //-----------------------------------------------------------------------------
838 // Function : Instance::updateIntermediateVars
839 // Purpose : updates a set of common variables used by RHS and jacobian
840 // Special Notes :
841 // Scope : public
842 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
843 // Creation Date : 03/21/2005
844 //-----------------------------------------------------------------------------
846 {
847  bool bsuccess = true;
848  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
849  {
850  Xyce::dout() << "Instance::updateIntermediateVars" << std::endl;
851  }
852  Linear::Vector & solVector = *(extData.nextSolVectorPtr);
853  Linear::Vector & lastSolVector = *(extData.lastSolVectorPtr);
854 
855  // some parameters in the model class that we will use often
856  const double A = model_.A;
857  const double Alpha = model_.Alpha;
858  const double Area = model_.Area;
859  const double BetaH = model_.BetaH;
860  const double BetaM = model_.BetaM;
861  const double C = model_.C;
862  const double DeltaV = model_.DeltaV;
863  const double Gap = model_.Gap;
864  const double Ms = model_.Ms;
865  const double Kirr = model_.Kirr;
866  const double Path = model_.Path;
867  const double Vinf = model_.Vinf;
868 
869  double latestMag; //NoMag = solVector[ li_MagVar ];
870 
871  //sum of currents through the inductors
872  branchCurrentSum = 0.0;
873  for(int i=0; i<numInductors; i++)
874  {
875  branchCurrentSum += solVector[instanceData[i]->li_Branch] * inductanceVals[ i ];
876  }
877 
878  latestMag = MagVar + MagVarUpdate;
879 
880  // used in voltage drop over first inductor
881  double V1Pos = solVector[(instanceData[0])->li_Pos];
882  double V1Neg = solVector[(instanceData[0])->li_Neg];
883 
884  double qV = (DeltaV / Vinf) * (V1Pos - V1Neg);
885 
886  double tanh_qV = 0.0;
887  if (fabs(qV) < CONSTTANH_THRESH)
888  {
889  tanh_qV = tanh(qV);
890  }
891  else if (qV < 0.0)
892  {
893  tanh_qV = -1.0;
894  }
895  else
896  {
897  tanh_qV = 1.0;
898  }
899 
900  double Happ = branchCurrentSum / Path;
901 
902 #ifdef MS_FACTORING2
903  double H = Happ - (Gap / Path) * latestMag * Ms;
904  double He = H + Alpha * latestMag * Ms;
905 #else
906  double H = Happ - (Gap / Path) * latestMag;
907  double He = H + Alpha * latestMag;
908 #endif
909 
910  double Heo = BetaH*A;
911 
912  // terms that come up frequently
913  double gap_path = Gap / Path;
914  double He2 = He*He;
915  double Heo2 = Heo*Heo;
916  double sq_Heo2He2 = sqrt(Heo2 + He2);
917 
918  double delM0 = model_.BetaM * Ms;
919  double Man = Ms * He / ( A + sq_Heo2He2 );
920 #ifdef MS_FACTORING2
921  double delM = Man - latestMag*Ms;
922 #else
923  double delM = Man - latestMag;
924 #endif
925 
926  double delM2 = delM*delM;
927  double delM02 = delM0*delM0;
928  double sq_delM02delM2 = sqrt( delM02 + delM2 );
929 
930  double Pold = P;
931 
932  #ifdef MS_FACTORING2
933  double Mirrp = (delM * tanh_qV + sq_delM02delM2 ) / (2*( Kirr- Alpha * sq_delM02delM2));
934  double Manp = Ms * (A + Heo2/sq_Heo2He2) / pow(A + sq_Heo2He2, 2.0);
935  P = ( C * Manp + (1 - C) * Mirrp) / ((1 + (gap_path - Alpha) * C * Manp + gap_path * (1-C) * Mirrp)*Ms);
936  #else
937  double Mirrp = (delM * tanh_qV + sq_delM02delM2 ) / (2*( Kirr- Alpha * sq_delM02delM2));
938  double Manp = Ms*(A + Heo2/sq_Heo2He2) / pow(A + sq_Heo2He2, 2.0);
939  P = ( C * Manp + (1 - C) * Mirrp) / (1 + (gap_path - Alpha) * C * Manp + gap_path * (1-C) * Mirrp);
940 #endif
941 
942 
943  // at this point we have P so now we can update mag.
944  //
945  // The problem is that if deltaM is too big, then we need to shrink
946  // the time step. One way to control this is to set a max time
947  // step. But what we really need to do is calculate deltaM and
948  // then if it's over some fraction of Ms then turn on the limiting
949  // flag (or bail on the step but I think turning on limiting is
950  // safer and if we hit maxItter with it on then we'll get that step
951  // rejected.
952 
953  if( useRKIntegration )
954  {
955  // use 4th order runga-kutta to estimate MagVarUpdate
957  MagVarUpdate = stepLen * ( PFunctionHistory[0] +
958  2*PFunctionHistory[1] +
959  2*PFunctionHistory[2] +
960  P) / 6;
961  }
962  else
963  {
964  // forward euler method
965  MagVarUpdate = P * (branchCurrentSum - oldBranchCurrentSum) / model_.Path;
966 
967  // trap
968  double MagVarUpdateWithTrap = 0.5 * (P + PPreviousStep) * (branchCurrentSum - oldBranchCurrentSum) / model_.Path;
969  origFlag = true;
970  if( fabs( MagVarUpdate ) > 0.25 * Ms )
971  {
972  // step was too big, so
973  // turn on limiting
974  origFlag = false;
975  }
976  }
977 
978  latestMag = MagVar + MagVarUpdate;
979 
980  if(includeDeltaM)
981  {
982  // in this case we're scaling MagVarUpdate by Ms because it's being solved
983  // with the full system and big changes in
984  MagVarUpdate /=model_.Ms;
985  }
986 
987  double dP_Denom = 1.0 + (gap_path - Alpha)*C*Manp + gap_path * (1.0-C) * Mirrp;
988 
989  // now get dP_dI for each inductor
990  for( int i=0; i<numInductors; i++)
991  {
992 
993  dHe_dI[ i ] = inductanceVals[ i ] / Path;
994 
995  dManp_dI[i] = ( -Ms * He / (pow(A + sq_Heo2He2, 2.0)*sq_Heo2He2)) *
996  ( (Heo2 / (Heo2 + He2)) + (2.0*(A + Heo2 / sq_Heo2He2)/(A+sq_Heo2He2)) ) * dHe_dI[i];
997 
998  ddelM_dI[i] = (Ms / (A + sq_Heo2He2)) * (1.0 - He2/((A + sq_Heo2He2)*sq_Heo2He2)) * dHe_dI[i];
999 
1000  dMirrp_dI[i] = (1.0/(2.0*(Kirr - Alpha*sq_delM02delM2))) *
1001  (tanh_qV + delM/sq_delM02delM2 +
1002  (2.0*Alpha*delM*(delM*tanh_qV +
1003  sq_delM02delM2)/(2.0*(Kirr-Alpha*sq_delM02delM2)*sq_delM02delM2))) * ddelM_dI[i];
1004 
1005  dP_dI[i] = (1.0/dP_Denom) * (C * dManp_dI[i] + (1.0-C) * dMirrp_dI[i]) -
1006  ( (C*Manp + (1.0-C)*Mirrp)/pow(dP_Denom,2.0) ) *
1007  ( (gap_path - Alpha)*C*dManp_dI[i] + gap_path*(1.0-C)*dMirrp_dI[i] );
1008 
1009  }
1010 
1011  return true;
1012 }
1013 
1014 //-----------------------------------------------------------------------------
1015 // Function : Instance::updateInductanceMatrix()
1016 // Purpose : A matrix of inductances is used often enough that it
1017 // calculated and stored as a member variable here
1018 // If and inductance ever changes say from updating
1019 // the temperature or a parameter udpate, then this
1020 // routine must be called again.
1021 // Special Notes :
1022 // Scope : public
1023 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1024 // Creation Date : 03/21/2005
1025 //-----------------------------------------------------------------------------
1027 {
1028  std::vector< InductorInstanceData* >::iterator
1029  currentInductor = instanceData.begin();
1030  std::vector< InductorInstanceData* >::iterator
1031  endInductor = instanceData.end();
1032 
1033  // collec the inductances
1034  int i=0;
1035  while( currentInductor != endInductor )
1036  {
1037  inductanceVals[ i ] = (*currentInductor)->L;
1038  i++;
1039  currentInductor++;
1040  }
1041 
1042  double Area = model_.Area;
1043  double Path = model_.Path;
1044 
1045  // compute the inductance matrix
1046  for( i=0; i<numInductors; ++i)
1047  {
1048  for( int j=0; j<numInductors; ++j)
1049  {
1050  // 4.0e-7 * M_PI is a magnetic constant, the permeability of free space [Henries/m]
1051  LO[i][j] = mutualCup * 4.0e-7 * M_PI * (Area / Path) * inductanceVals[i] * inductanceVals[j];
1052  }
1053  }
1054 
1055 }
1056 
1057 //-----------------------------------------------------------------------------
1058 // Function : Instance::updatePrimaryState
1059 // Purpose :
1060 // Special Notes :
1061 // Scope : public
1062 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1063 // Creation Date : 03/21/2005
1064 //-----------------------------------------------------------------------------
1066 {
1067  bool bsuccess = true;
1068  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1069  {
1070  Xyce::dout() << "Instance::updatePrimaryState---------------" << std::endl
1071  << "\tname = " << getName() << std::endl;
1072  }
1073  // udate dependent parameters
1075 
1076 #if 0
1077  // don't need to do this as we're not using the state vector
1078  Linear::Vector & solVector = *(extData.nextSolVectorPtr);
1079  Linear::Vector & staVector = *(extData.nextStaVectorPtr);
1080 
1081  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1082  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1083  int i = 0;
1084  while( currentInductor != endInductor )
1085  {
1086  double current = solVector[ ( (*currentInductor)->li_Branch) ];
1087  if( (getSolverState().dcopFlag) && ((*currentInductor)->ICGiven) )
1088  {
1089  current = (*currentInductor)->IC;
1090  }
1091  // place this value for the charge in the state vector.
1092  staVector[((*currentInductor)->li_currentState)] = current;
1093  currentInductor++;
1094  i++;
1095  }
1096 #endif
1097 
1098  return bsuccess;
1099 }
1100 
1101 //-----------------------------------------------------------------------------
1102 // Function : Instance::updateSecondaryState
1103 // Purpose :
1104 // Special Notes :
1105 // Scope : public
1106 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1107 // Creation Date : 03/21/2005
1108 //-----------------------------------------------------------------------------
1110 {
1111  bool bsuccess = true;
1112  return bsuccess;
1113 }
1114 
1115 //-----------------------------------------------------------------------------
1116 // Function : Instance::acceptStep
1117 // Purpose : This function updates the value of MagVar
1118 //
1119 // Scope : public
1120 // Creator : Rich Schiek, SNL, Electrical Systems Modeling
1121 // Creation Date : 01/25/2011
1122 //-----------------------------------------------------------------------------
1124 {
1125  if (!getSolverState().dcopFlag)
1126  {
1127  if(includeDeltaM)
1128  {
1130  }
1131  else
1132  {
1133  MagVar += MagVarUpdate;
1134  }
1136  PPreviousStep = P;
1137  if( fabs(MagVar) > 2*model_.Ms )
1138  {
1139  MagVar = 0.0;
1140  }
1141 
1142  if( useRKIntegration )
1143  {
1144  // fill in history for RK integration of dM/dH
1145  for(int i=0; i<2; i++)
1146  {
1149  }
1152  }
1154  }
1155 }
1156 
1157 //-----------------------------------------------------------------------------
1158 // Function : Instance::loadDAEQVector
1159 //
1160 // Purpose : Loads the Q-vector contributions for a single
1161 // instance.
1162 //
1163 // Special Notes : The "Q" vector is part of a standard DAE formalism in
1164 // which the system of equations is represented as:
1165 //
1166 // f(x) = dQ(x)/dt + F(x) - B(t) = 0
1167 //
1168 // Scope : public
1169 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1170 // Creation Date : 03/21/2005
1171 //-----------------------------------------------------------------------------
1173 {
1174  bool bsuccess = true;
1175 
1176  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1177  {
1178  Xyce::dout() << "Instance::loadDAEQVector------------------------" << std::endl
1179  << "\tname = " << getName() << std::endl;
1180  }
1181 
1182  Linear::Vector * daeQVecPtr = extData.daeQVectorPtr;
1183  Linear::Vector & solVector = *(extData.nextSolVectorPtr);
1184 
1185  // update LOI -- the following product
1186  // I = column vector of currents
1187  // L = row vector of inductances
1188  // LO = matrix = mutualCup * sqrt( L' * L )
1189  // LOI = column vector = mutualCup * sqrt( L' * L ) * I
1190  // LOI[1] = mutualCup * sqrt(L[1]*L[1])*I[1]) +
1191  // mutualCup * sqrt(L[1]*L[2])*I[2]) + ...
1192  // mutualCup * sqrt(L[1]*L[n])*I[n])
1193 
1194  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1195  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1196  int i = 0;
1197  while( currentInductor != endInductor )
1198  {
1199  if( (getSolverState().dcopFlag) && (*currentInductor)->ICGiven == true )
1200  {
1201  inductorCurrents[ i ] = (*currentInductor)->IC;
1202  }
1203  else
1204  {
1205  inductorCurrents[ i ] = solVector[ (*currentInductor)->li_Branch ];
1206  }
1207  i++;
1208  currentInductor++;
1209  }
1210 
1211  for( i = 0; i < numInductors; ++i )
1212  {
1213  LOI[ i ] = 0;
1214  for( int j = 0; j < numInductors; ++j )
1215  {
1216  LOI[i] += LO[i][j] * inductorCurrents[j];
1217  }
1218  }
1219 
1220  // loop over each inductor and load it's Q vector components
1221  // and each inductor's contribution to the R equ.
1222  currentInductor = instanceData.begin();
1223  endInductor = instanceData.end();
1224  i = 0;
1225  while( currentInductor != endInductor )
1226  {
1227 
1228  (*daeQVecPtr)[((*currentInductor)->li_Branch)] += LOI[ i ];
1229  i++;
1230  currentInductor++;
1231  }
1232 
1233  return bsuccess;
1234 }
1235 
1236 
1237 //-----------------------------------------------------------------------------
1238 // Function : Instance::loadDAEFVector
1239 //
1240 // Purpose : Loads the F-vector contributions for a single
1241 // instance.
1242 //
1243 // Special Notes : See the special notes for loadDAEFVector.
1244 //
1245 // Same as loadRHS, but without the capacitor
1246 // currents.
1247 //
1248 // Scope : public
1249 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1250 // Creation Date : 03/21/2005
1251 //-----------------------------------------------------------------------------
1253 {
1254  bool bsuccess=true;
1255 
1256  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1257  {
1258  Xyce::dout() << "Instance::loadDAEFVector------------------------" << std::endl
1259  << "\tname = " << getName() << std::endl;
1260  }
1261 
1262  Linear::Vector * daeFVecPtr = extData.daeFVectorPtr;
1263  Linear::Vector & solVector = *(extData.nextSolVectorPtr);
1264 
1265  double Gap = model_.Gap;
1266  double Path = model_.Path;
1267 
1268 
1269  // used in scaling the branch equation;
1270  double mid = 1.0 + (1.0 - (Gap / Path))*P;
1271 
1272  // loop over each inductor and load it's F vector components
1273  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1274  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1275  int i=0;
1276  while( currentInductor != endInductor )
1277  {
1278  double current = solVector[(*currentInductor)->li_Branch];
1279  double vNodePos = solVector[(*currentInductor)->li_Pos];
1280  double vNodeNeg = solVector[(*currentInductor)->li_Neg];
1281 
1282 
1283  (*daeFVecPtr)[((*currentInductor)->li_Pos)] += current;
1284 
1285  (*daeFVecPtr)[((*currentInductor)->li_Neg)] += -current;
1286 
1287  (*daeFVecPtr)[((*currentInductor)->li_Branch)] += -((vNodePos - vNodeNeg)/mid);
1288 
1289  currentInductor++;
1290  i++;
1291  }
1292 
1293  if(includeDeltaM)
1294  {
1295  // the deltaHapp equation
1296  //(*daeFVecPtr)[li_deltaHappVar] += solVector[li_deltaHappVar];
1297  //(*daeFVecPtr)[li_deltaHappVar] -= HappVarUpdate;
1298 
1299  // the deltaM equation
1300  (*daeFVecPtr)[li_deltaMagVar] += solVector[li_deltaMagVar];
1301  (*daeFVecPtr)[li_deltaMagVar] -= MagVarUpdate;
1302  }
1303 
1304  return bsuccess;
1305 }
1306 
1307 //-----------------------------------------------------------------------------
1308 // Function : Instance::loadDAEdQdx
1309 //
1310 // Purpose : Loads the Q-vector contributions for a single
1311 // instance.
1312 // Scope : public
1313 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1314 // Creation Date : 03/21/2005
1315 //-----------------------------------------------------------------------------
1317 {
1318  bool bsuccess = true;
1319  int i;
1320 
1321  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1322  {
1323  Xyce::dout() << "Instance::loadDAEdQdx-----------------------" << std::endl
1324  << "\tname = " << getName() << std::endl;
1325  }
1326  Linear::Matrix * dQdxMatPtr = extData.dQdxMatrixPtr;
1327 
1328  // loop over each inductor and load it's Q vector components
1329  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1330  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1331  i = 0;
1332  while( currentInductor != endInductor )
1333  {
1334  for( int j=0; j<numInductors; ++j )
1335  {
1336  (*dQdxMatPtr)[((*currentInductor)->li_Branch)]
1337  [(*currentInductor)->inductorCurrentOffsets[j]] += LO[i][j];
1338  }
1339  i++;
1340  currentInductor++;
1341  }
1342 
1343  return bsuccess;
1344 }
1345 
1346 
1347 
1348 //-----------------------------------------------------------------------------
1349 // Function : Instance::loadDAEdFdx ()
1350 //
1351 // Purpose : Loads the F-vector contributions for a single
1352 // instance.
1353 //
1354 // Special Notes :
1355 //
1356 // Scope : public
1357 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1358 // Creation Date : 03/21/2005
1359 //-----------------------------------------------------------------------------
1361 {
1362  bool bsuccess = true;
1363 
1364  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1365  {
1366  Xyce::dout() << "Instance::loadDAEdFdx----------------------" << std::endl
1367  << "\tname = " << getName() << std::endl;
1368  }
1369 
1370  Linear::Vector & solVector = *(extData.nextSolVectorPtr);
1371  Linear::Vector & lastSolVector = *(extData.lastSolVectorPtr);
1372  Linear::Matrix * dFdxMatPtr = extData.dFdxMatrixPtr;
1373 
1374  // pull these parameters up from the model class to make it easier
1375  // to view the equations.
1376  const double Gap = model_.Gap;
1377  const double Path = model_.Path;
1378 
1379 
1380  // loop over each inductor and load it's dFdx components
1381 #ifdef MS_FACTORING2
1382  double mid = 1.0 + (1.0 - (Gap/Path))*P*model_.Ms;
1383 #else
1384  double mid = 1.0 + (1.0 - (Gap/Path))*P;
1385 #endif
1386 
1387  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1388  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1389  while( currentInductor != endInductor )
1390  {
1391  // do the normal work for an inductor
1392  (*dFdxMatPtr)[((*currentInductor)->li_Pos)] [((*currentInductor)->APosEquBraVarOffset)] += 1.0;
1393  (*dFdxMatPtr)[((*currentInductor)->li_Neg)] [((*currentInductor)->ANegEquBraVarOffset)] += -1.0;
1394  (*dFdxMatPtr)[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquPosNodeOffset)] += -1.0/mid;
1395  (*dFdxMatPtr)[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquNegNodeOffset)] += 1.0/mid;
1396 
1397  double delV = solVector[(*currentInductor)->li_Pos] - solVector[(*currentInductor)->li_Neg];
1398 
1399  for( int j = 0; j<numInductors; ++j )
1400  {
1401 
1402  (*dFdxMatPtr)[((*currentInductor)->li_Branch)][(*currentInductor)->inductorCurrentOffsets[j]] +=
1403  delV * (1.0 - (Gap/Path)) * dP_dI[j]/(mid*mid);
1404 
1405  }
1406 
1407  currentInductor++;
1408  }
1409 
1410  if(includeDeltaM)
1411  {
1412  (*dFdxMatPtr)[li_deltaMagVar][deltaMEquDeltaMOffset] = 1.0;
1413 
1414  for( int i=0; i<numInductors; i++ )
1415  {
1416  (*dFdxMatPtr)[li_deltaMagVar][deltaMEquInductorOffsets[i]] =
1417  -((inductanceVals[i]*(solVector[ instanceData[i]->li_Branch ] - lastSolVector[ instanceData[i]->li_Branch ] ) * dP_dI[i])
1418  + P*inductanceVals[i])/(model_.Path*model_.Ms );
1419  }
1420  }
1421 
1422  return bsuccess;
1423 }
1424 
1425 //-----------------------------------------------------------------------------
1426 // Function : Instance::outputPlotFiles
1427 // Purpose : If requested by the use in the model statement,
1428 // this routine outputs values of the internal
1429 // state variables M, H and R to a file
1430 // named "Inductor_name.dat". File is opened
1431 // and closed in the contructor and destructor.
1432 // Special Notes :
1433 // Scope : public
1434 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1435 // Creation Date : 03/21/2005
1436 //-----------------------------------------------------------------------------
1438 {
1439  bool bsuccess = true;
1440  if( outputStateVarsFlag && outputFileStreamPtr.get() && (*outputFileStreamPtr) )
1441  {
1442  Linear::Vector & solVector = *(extData.nextSolVectorPtr);
1443 
1444 #ifdef MS_FACTORING2
1445  double latestMag = MagVar*model_.Ms;
1446 #else
1447  double latestMag = MagVar;
1448 #endif
1449 
1450  if( includeDeltaM )
1451  {
1452  (*outputFileStreamPtr)
1453  << getSolverState().currTime << " "
1454  << latestMag
1455  << std::endl;
1456  }
1457  else
1458  {
1459  (*outputFileStreamPtr)
1460  << getSolverState().currTime << " "
1461  << latestMag
1462  << std::endl;
1463  }
1464 
1465  }
1466  return bsuccess;
1467 }
1468 
1469 
1470 //-----------------------------------------------------------------------------
1471 // Function : Instance::setIC
1472 // Purpose :
1473 // Special Notes :
1474 // Scope : public
1475 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1476 // Creation Date : 03/21/2005
1477 //-----------------------------------------------------------------------------
1479 {
1480  int i_bra_sol;
1481  int i_f_state;
1482 
1483  bool bsuccess = true;
1484  return bsuccess;
1485 }
1486 
1487 //-----------------------------------------------------------------------------
1488 // Function : Instance::varTypes
1489 // Purpose :
1490 // Special Notes :
1491 // Scope : public
1492 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1493 // Creation Date : 03/21/2005
1494 //-----------------------------------------------------------------------------
1495 void Instance::varTypes( std::vector<char> & varTypeVec )
1496 {
1497  varTypeVec.resize(numInductors);
1498  for(int i=0; i<numInductors; i++)
1499  {
1500  varTypeVec[i] = 'I';
1501  }
1502 }
1503 
1504 
1505 
1506 //-----------------------------------------------------------------------------
1507 // Function : Instance::Pcalc
1508 // Purpose :
1509  // this is a templated function for a complicated term P(M,I_1... I_n) that relates
1510  // the magnetic saturation of the mutual indcutor to the individual currents
1511  // through the inductors. We'll need dP_dM and this tempated function automates
1512  // that calculation via Sacado
1513 // Special Notes :
1514 // Scope : public
1515 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1516 // Creation Date : 10/13/2011
1517 //-----------------------------------------------------------------------------
1518 template <typename ScalarT>
1520  const ScalarT & Mag, const ScalarT & CurrentSum,
1521  const ScalarT & Vpos, const ScalarT & Vneg)
1522  // independent variable M
1523  // independent variable Sum(n_i * I_i) windings * current through each inductor
1524  // independent variable Vpos, Vneg -- voltage drop over first inductor
1525 {
1526  // some parameters in the model class that we will use often
1527  const double A = model_.A;
1528  const double Alpha = model_.Alpha;
1529  const double Area = model_.Area;
1530  const double BetaH = model_.BetaH;
1531  const double BetaM = model_.BetaM;
1532  const double C = model_.C;
1533  const double DeltaV = model_.DeltaV;
1534  const double Gap = model_.Gap;
1535  const double Ms = model_.Ms;
1536  const double Kirr = model_.Kirr;
1537  const double Path = model_.Path;
1538  const double Vinf = model_.Vinf;
1539 
1540  ScalarT qV = (DeltaV / Vinf) * (Vpos - Vneg);
1541 
1542  ScalarT tanh_qV = 0.0;
1543  if (fabs(qV) < CONSTTANH_THRESH)
1544  {
1545  ScalarT tanh_qV = tanh(qV);
1546  }
1547  else if (qV < 0.0)
1548  {
1549  tanh_qV = -1.0;
1550  }
1551  else
1552  {
1553  tanh_qV = 1.0;
1554  }
1555 
1556  ScalarT Happ = CurrentSum / Path;
1557 
1558 #ifdef MS_FACTORING2
1559  ScalarT H = Happ - (Gap / Path) * Mag * Ms;
1560  ScalarT He = H + Alpha * Mag * Ms;
1561 #else
1562  ScalarT H = Happ - (Gap / Path) * Mag;
1563  ScalarT He = H + Alpha * Mag;
1564 #endif
1565 
1566  ScalarT Heo = BetaH*A;
1567 
1568  // terms that come up frequently
1569  ScalarT gap_path = Gap / Path;
1570  ScalarT He2 = He*He;
1571  ScalarT Heo2 = Heo*Heo;
1572  ScalarT sq_Heo2He2 = sqrt(Heo2 + He2);
1573 
1574  ScalarT delM0 = model_.BetaM * Ms;
1575  ScalarT Man = Ms * He / ( A + sq_Heo2He2 );
1576 #ifdef MS_FACTORING2
1577  ScalarT delM = Man - Mag*Ms;
1578 #else
1579  ScalarT delM = Man - Mag;
1580 #endif
1581 
1582  ScalarT delM2 = delM*delM;
1583  ScalarT delM02 = delM0*delM0;
1584  ScalarT sq_delM02delM2 = sqrt( delM02 + delM2 );
1585 
1586 #ifdef MS_FACTORING2
1587  ScalarT Mirrp = (delM * tanh_qV + sq_delM02delM2 ) / (2*( Kirr- Alpha * sq_delM02delM2));
1588  ScalarT Manp = Ms * (A + Heo2/sq_Heo2He2) / pow(A + sq_Heo2He2, 2.0);
1589  ScalarT Pval = ( C * Manp + (1 - C) * Mirrp) / ((1 + (gap_path - Alpha) * C * Manp + gap_path * (1-C) * Mirrp)*Ms);
1590 #else
1591  ScalarT Mirrp = (delM * tanh_qV + sq_delM02delM2 ) / (2*( Kirr- Alpha * sq_delM02delM2));
1592  ScalarT Manp = Ms*(A + Heo2/sq_Heo2He2) / pow(A + sq_Heo2He2, 2.0);
1593  ScalarT Pval = ( C * Manp + (1 - C) * Mirrp) / (1 + (gap_path - Alpha) * C * Manp + gap_path * (1-C) * Mirrp);
1594 #endif
1595 
1596  return Pval;
1597 } // end of Pcalc() function
1598 
1599 
1600 //-----------------------------------------------------------------------------
1601 // Function : Model::processParams
1602 // Purpose :
1603 // Special Notes :
1604 // Scope : public
1605 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1606 // Creation Date : 03/21/2005
1607 //-----------------------------------------------------------------------------
1609 {
1610  return true;
1611 }
1612 
1613 //----------------------------------------------------------------------------
1614 // Function : Model::processInstanceParams
1615 // Purpose :
1616 // Special Notes :
1617 // Scope : public
1618 // Creator : Dave Shirely, PSSI
1619 // Creation Date : 03/23/06
1620 //----------------------------------------------------------------------------
1622 {
1623  std::vector<Instance*>::iterator iter;
1624  std::vector<Instance*>::iterator first = instanceContainer.begin();
1625  std::vector<Instance*>::iterator last = instanceContainer.end();
1626 
1627  for (iter=first; iter!=last; ++iter)
1628  {
1629  (*iter)->processParams();
1630  }
1631 
1632  return true;
1633 }
1634 
1635 //-----------------------------------------------------------------------------
1636 // Function : Model::Model
1637 // Purpose : block constructor
1638 // Special Notes :
1639 // Scope : public
1640 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1641 // Creation Date : 03/21/2005
1642 //-----------------------------------------------------------------------------
1644  const Configuration & configuration,
1645  const ModelBlock & MB,
1646  const FactoryBlock & factory_block)
1647  : DeviceModel(MB, configuration.getModelParameters(), factory_block),
1648  A(0.0),
1649  Alpha(0.0),
1650  Area(0.0),
1651  BetaH(0.0),
1652  BetaM(0.0),
1653  C(0.0),
1654  DeltaV(0.0),
1655  Gap(0.0),
1656  Kirr(0.0),
1657  Ms(0.0),
1658  Path(0.0),
1659  Vinf(0.0),
1660  tempCoeff1(0.0),
1661  tempCoeff2(0.0),
1662  outputStateVars(0),
1663  useRKIntegration(0),
1664  includeDeltaM(0),
1665  tnom(getDeviceOptions().tnom)
1666 {
1667  setLevel(2);
1668 
1669 
1670  // Set params to constant default values:
1671  setDefaultParams ();
1672 
1673  // Set params according to .model line and constant defaults from metadata:
1674  setModParams (MB.params);
1675 
1676  // scale gap, path and area from cm and cm^2 to m and m^2
1677  Gap *= 1.0e-2;
1678  Path *= 1.0e-2;
1679  Area *= 1.0e-4;
1680 
1681  // Set any non-constant parameter defaults:
1682 
1683  if (!given("TNOM"))
1684  {
1686  }
1687 
1688  // Calculate any parameters specified as expressions:
1690 
1691  // calculate dependent (ie computed) params and check for errors:
1692  processParams ();
1693 }
1694 
1695 //-----------------------------------------------------------------------------
1696 // Function : Model::~Model
1697 // Purpose : destructor
1698 // Special Notes :
1699 // Scope : public
1700 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1701 // Creation Date : 03/21/2005
1702 //-----------------------------------------------------------------------------
1704 {
1705  std::vector<Instance*>::iterator iter;
1706  std::vector<Instance*>::iterator first = instanceContainer.begin();
1707  std::vector<Instance*>::iterator last = instanceContainer.end();
1708 
1709  for (iter=first; iter!=last; ++iter)
1710  {
1711  delete (*iter);
1712  }
1713 }
1714 
1715 // additional Declarations
1716 
1717 //-----------------------------------------------------------------------------
1718 // Function : Model::printOutInstances
1719 // Purpose : debugging tool.
1720 // Special Notes :
1721 // Scope : public
1722 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1723 // Creation Date : 03/21/2005
1724 //-----------------------------------------------------------------------------
1725 std::ostream &Model::printOutInstances(std::ostream &os) const
1726 {
1727  std::vector<Instance*>::const_iterator iter;
1728  std::vector<Instance*>::const_iterator first = instanceContainer.begin();
1729  std::vector<Instance*>::const_iterator last = instanceContainer.end();
1730 
1731  int i, isize;
1732  isize = instanceContainer.size();
1733 
1734  os << std::endl;
1735  os << "Number of MutIndNonLin instances: " << isize << std::endl;
1736  os << " name=\t\tmodelName\tParameters" << std::endl;
1737  for (i=0, iter=first; iter!=last; ++iter, ++i)
1738  {
1739  os << " " << i << ": " << (*iter)->getName() << "\t";
1740  os << getName();
1741  os << std::endl;
1742  }
1743 
1744  os << std::endl;
1745 
1746  return os;
1747 }
1748 
1749 //-----------------------------------------------------------------------------
1750 // Function : Model::forEachInstance
1751 // Purpose :
1752 // Special Notes :
1753 // Scope : public
1754 // Creator : David Baur
1755 // Creation Date : 2/4/2014
1756 //-----------------------------------------------------------------------------
1757 /// Apply a device instance "op" to all instances associated with this
1758 /// model
1759 ///
1760 /// @param[in] op Operator to apply to all instances.
1761 ///
1762 ///
1763 void Model::forEachInstance(DeviceInstanceOp &op) const /* override */
1764 {
1765  for (std::vector<Instance *>::const_iterator it = instanceContainer.begin(); it != instanceContainer.end(); ++it)
1766  op(*it);
1767 }
1768 
1769 
1770 Device *Traits::factory(const Configuration &configuration, const FactoryBlock &factory_block)
1771 {
1772 
1773  return new DeviceMaster<Traits>(configuration, factory_block, factory_block.solverState_, factory_block.deviceOptions_);
1774 }
1775 
1777 {
1779  .registerDevice("min", 2)
1780  .registerModelType("min", 2);
1781 }
1782 
1783 } // namespace MutIndNonLin2
1784 } // namespace Device
1785 } // namespace Xyce
const InstanceName & getName() const
static void loadModelParameters(ParametricData< Model > &model_parameters)
void registerLIDs(const std::vector< int > &intLIDVecRef, const std::vector< int > &extLIDVecRef)
Teuchos::RCP< std::ofstream > outputFileStreamPtr
const SolverState & solverState_
Descriptor & addPar(const char *parName, T default_value, T U::*varPtr)
Adds the parameter description to the parameter map.
Definition: N_DEV_Pars.h:1429
std::vector< std::string > couplingInductor
Linear::Vector * nextSolVectorPtr
bool given(const std::string &parameter_name) const
void varTypes(std::vector< char > &varTypeVec)
Linear::Vector * daeQVectorPtr
Pure virtual class to augment a linear system.
void registerStateLIDs(const std::vector< int > &staLIDVecRef)
void addInternalNode(Util::SymbolTable &symbol_table, int index, const InstanceName &instance_name, const std::string &lead_name)
virtual void forEachInstance(DeviceInstanceOp &op) const
Apply a device instance "op" to all instances associated with this model.
void loadNodeSymbols(Util::SymbolTable &symbol_table) const
Populates and returns the store name map.
#define AssertLIDs(cmp)
virtual void registerJacLIDs(const JacobianStamp &jacLIDVec)
bool updateTemperature(const double &temp_tmp)
std::vector< Instance * > instanceContainer
std::vector< Param > params
Parameters from the line.
Linear::Vector * lastSolVectorPtr
std::vector< InductorInstanceData * > instanceData
void setParams(const std::vector< Param > &params)
const std::string & getName() const
The FactoryBlock contains parameters needed by the device, instance and model creation functions...
const DeviceOptions & getDeviceOptions() const
static Device * factory(const Configuration &configuration, const FactoryBlock &factory_block)
const std::vector< std::vector< int > > & jacobianStamp() const
const DeviceOptions & deviceOptions_
Linear::Vector * nextStaVectorPtr
static Config< T > & addConfiguration()
Adds the device to the Xyce device configuration.
Linear::Matrix * dFdxMatrixPtr
static void loadInstanceParameters(ParametricData< Instance > &instance_parameters)
The Device class is an interface for device implementations.
Definition: N_DEV_Device.h:101
virtual std::ostream & printOutInstances(std::ostream &os) const
bool processInstanceParams()
processInstanceParams
std::vector< std::vector< int > > jacStamp
std::vector< std::string > inductorsNode1
Class Configuration contains device configuration data.
std::vector< std::string > inductorNames
#define CONSTTANH_THRESH
Definition: N_DEV_Const.h:114
const SolverState & getSolverState() const
void registerJacLIDs(const std::vector< std::vector< int > > &jacLIDVec)
#define M_PI
std::vector< int > inductorCurrentOffsets
Instance(const Configuration &configuration, const InstanceBlock &IB, Model &Iiter, const FactoryBlock &factory_block)
Linear::Vector * daeFVectorPtr
ScalarT Pcalc(const ScalarT &Mag, const ScalarT &CurrentSum, const ScalarT &Vpos, const ScalarT &Vneg)
std::vector< std::vector< double > > LO
ModelBlock represents a .MODEL line from the netlist.
Manages parameter binding for class C.
Definition: N_DEV_Pars.h:214
std::vector< std::string > inductorsNode2
InstanceBlock represent a device instance line from the netlist.
std::vector< Param > params
Linear::Matrix * dQdxMatrixPtr
void setModParams(const std::vector< Param > &params)