Xyce  6.1
N_DEV_MutIndNonLin.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_MutIndNonLin.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.174 $
40 //
41 // Revision Date : $Date: 2015/04/24 20:25:45 $
42 //
43 // Current Owner : $Author: dgbaur $
44 //-------------------------------------------------------------------------
45 #include <Xyce_config.h>
46 
47 // ---------- Standard Includes ----------
48 #include <fstream>
49 #include <algorithm>
50 #include <vector>
51 #include <set>
52 
53 // ---------- Xyce Includes ----------
54 #include <N_DEV_DeviceOptions.h>
55 #include <N_DEV_DeviceMaster.h>
56 #include <N_DEV_ExternData.h>
57 #include <N_DEV_MatrixLoadData.h>
58 #include <N_DEV_MutIndNonLin.h>
59 #include <N_DEV_SolverState.h>
60 #include <N_DEV_Message.h>
61 #include <N_ERH_ErrorMgr.h>
62 #include <N_UTL_FeatureTest.h>
63 #include <N_UTL_Math.h>
64 
65 //This contains important constants like permitivity of free space
66 #include <N_DEV_Const.h>
67 
68 #include <N_LAS_Vector.h>
69 #include <N_LAS_Matrix.h>
70 
71 using Teuchos::rcp;
72 
73 namespace Xyce {
74 namespace Device {
75 
76 namespace MutIndNonLin {
77 
79 {
80  p.addPar("COUP_VAL",1.0,&MutIndNonLin::Instance::mutualCup)
82  .setUnit(U_NONE)
83  .setCategory(CAT_NONE)
84  .setDescription("Coupling coefficient");
85 
86  p.addPar("NONLINEARCOUPLING",0.0,&MutIndNonLin::Instance::nonlinFlag)
88  .setUnit(U_NONE)
89  .setCategory(CAT_NONE)
90  .setDescription("Nonlinear coupling flag");
91 
92  p.addPar("COUPLEDMutIndNonLin",std::vector<std::string>(),&MutIndNonLin::Instance::inductorNames)
93  .setUnit(U_NONE)
94  .setCategory(CAT_NONE)
95  .setDescription("");
96 
97  p.addPar("COUPLEDINDUCTANCE",std::vector<double>(),&MutIndNonLin::Instance::inductorInductances)
98  .setUnit(U_NONE)
99  .setCategory(CAT_NONE)
100  .setDescription("");
101 
102  p.addPar("NODE1",std::vector<std::string>(),&MutIndNonLin::Instance::inductorsNode1)
103  .setUnit(U_NONE)
104  .setCategory(CAT_NONE)
105  .setDescription("");
106 
107  p.addPar("NODE2",std::vector<std::string>(),&MutIndNonLin::Instance::inductorsNode2)
108  .setUnit(U_NONE)
109  .setCategory(CAT_NONE)
110  .setDescription("");
111 
112  p.addPar("COUPLING",std::vector<double>(),&MutIndNonLin::Instance::couplingCoefficient)
113  .setUnit(U_NONE)
114  .setCategory(CAT_NONE)
115  .setDescription("Coupling coefficient");
116 
117  p.addPar("COUPLEDINDUCTOR",std::vector<std::string>(),&MutIndNonLin::Instance::couplingInductor)
118  .setUnit(U_NONE)
119  .setCategory(CAT_NONE)
120  .setDescription("");
121 }
122 
124 {
125  p.addPar("A",1000.0,&MutIndNonLin::Model::A)
126  .setUnit(U_AMPMM1)
127  .setCategory(CAT_MATERIAL)
128  .setDescription("Thermal energy parameter");
129 
130  p.addPar("AREA",0.1,&MutIndNonLin::Model::Area)
131  .setUnit(U_CM2)
132  .setCategory(CAT_GEOMETRY)
133  .setDescription("Mean magnetic cross-sectional area");
134 
135  p.addPar("ALPHA",5.0e-5,&MutIndNonLin::Model::Alpha)
136  .setUnit(U_NONE)
137  .setCategory(CAT_GEOMETRY)
138  .setDescription("Domain coupling parameter");
139 
140  p.addPar("BETAH",0.0001,&MutIndNonLin::Model::BetaH)
141  .setUnit(U_NONE)
142  .setCategory(CAT_NONE)
143  .setDescription("Modeling constant");
144 
145  p.addPar("BETAM",3.125e-5,&MutIndNonLin::Model::BetaM)
146  .setUnit(U_NONE)
147  .setCategory(CAT_NONE)
148  .setDescription("Modeling constant");
149 
150  p.addPar("C",0.2,&MutIndNonLin::Model::C)
151  .setUnit(U_NONE)
152  .setCategory(CAT_MATERIAL)
153  .setDescription("Domain flexing parameter");
154 
155  p.addPar("CLIM",0.005,&MutIndNonLin::Model::CLim)
156  .setUnit(U_NONE)
157  .setCategory(CAT_MATERIAL)
158  .setDescription("Value below which domain flexing parameter will be treated as zero.");
159 
160  p.addPar("DELVSCALING",1.0e3,&MutIndNonLin::Model::DeltaVScaling)
161  .setUnit(U_VOLT)
162  .setCategory(CAT_NONE)
163  .setDescription("Smoothing coefficient for voltage difference over first inductor");
164 
165  p.addPar("CONSTDELVSCALING",false,&MutIndNonLin::Model::UseConstantDeltaVScaling)
166  .setUnit(U_VOLT)
167  .setCategory(CAT_NONE)
168  .setDescription("Use constant scaling factor to smooth voltage difference over first inductor");
169 
170  p.addPar("INCLUDEMEQU",true,&MutIndNonLin::Model::includeMEquation)
172  .setUnit(U_NONE)
173  .setCategory(CAT_NONE)
174  .setDescription("Flag to include the magnetics in the solution.");
175 
176  p.addPar("GAP",0.0,&MutIndNonLin::Model::Gap)
177  .setUnit(U_CM)
178  .setCategory(CAT_GEOMETRY)
179  .setDescription("Effective air gap");
180 
181  p.addPar("K",500.0,&MutIndNonLin::Model::Kirr)
182  .setUnit(U_AMPMM1)
183  .setCategory(CAT_MATERIAL)
184  .setDescription("Domain anisotropy parameter");
185 
186  p.addPar("KIRR",500.0,&MutIndNonLin::Model::Kirr)
187  .setUnit(U_AMPMM1)
188  .setCategory(CAT_MATERIAL)
189  .setDescription("Domain anisotropy parameter");
190 
191  p.addPar("MS",1.0e+6,&MutIndNonLin::Model::Ms)
192  .setUnit(U_AMPMM1)
193  .setCategory(CAT_MATERIAL)
194  .setDescription("Saturation magnetization");
195 
197  .setUnit(U_NONE)
198  .setCategory(CAT_NONE)
199  .setDescription("for pspice compatibility -- ignored");
200 
202  .setUnit(U_NONE)
203  .setCategory(CAT_NONE)
204  .setDescription("for pspice compatibility -- ignored");
205 
206  p.addPar("PATH",1.0,&MutIndNonLin::Model::Path)
207  .setUnit(U_CM)
208  .setCategory(CAT_GEOMETRY)
209  .setDescription("Total mean magnetic path");
210 
211  p.addPar("TNOM",27.0,&MutIndNonLin::Model::tnom)
212  .setUnit(U_DEGC)
213  .setCategory(CAT_MATERIAL)
214  .setDescription("Reference temperature");
215 
217  .setUnit(U_NONE)
218  .setCategory(CAT_MATERIAL)
219  .setDescription("First order temperature coeff.");
220 
222  .setUnit(U_NONE)
223  .setCategory(CAT_MATERIAL)
224  .setDescription("Second order temperature coeff.");
225 
226  p.addPar("PZEROTOL",0.1,&MutIndNonLin::Model::pZeroTol)
227  .setUnit(U_NONE)
228  .setCategory(CAT_NONE)
229  .setDescription("Tolerance for nonlinear zero crossing");
230 
231  p.addPar("MVARSCALING",1.0,&MutIndNonLin::Model::mVarScaling)
232  .setUnit(U_NONE)
233  .setCategory(CAT_NONE)
234  .setDescription("M-variable scaling.");
235 
236  p.addPar("RVARSCALING",1.0,&MutIndNonLin::Model::rVarScaling)
237  .setUnit(U_NONE)
238  .setCategory(CAT_NONE)
239  .setDescription("R-variable scaling");
240 
241  p.addPar("MEQNSCALING",1.0,&MutIndNonLin::Model::mEqScaling)
242  .setUnit(U_NONE)
243  .setCategory(CAT_NONE)
244  .setDescription("M-equation scaling");
245 
246  p.addPar("REQNSCALING",1.0,&MutIndNonLin::Model::rEqScaling)
247  .setUnit(U_NONE)
248  .setCategory(CAT_NONE)
249  .setDescription("R-equation scaling");
250 
251  p.addPar("OUTPUTSTATEVARS",0,&MutIndNonLin::Model::outputStateVars)
252  .setUnit(U_NONE)
253  .setCategory(CAT_NONE)
254  .setDescription("Flag to save state variables");
255 
256  p.addPar("FACTORMS",0,&MutIndNonLin::Model::factorMS)
257  .setUnit(U_NONE)
258  .setCategory(CAT_NONE)
259  .setDescription("Flag to save state variables");
260 
261  p.addPar("BHSIUNITS",0,&MutIndNonLin::Model::BHSiUnits)
262  .setUnit(U_NONE)
263  .setCategory(CAT_NONE)
264  .setDescription("Flag to report B and H in SI units");
265 }
266 
267 
268 // Class Instance
269 
270 //-----------------------------------------------------------------------------
271 // Function : Instance::Instance
272 // Purpose : constructor
273 // Special Notes :
274 // Scope : public
275 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
276 // Creation Date : 03/21/2005
277 //-----------------------------------------------------------------------------
279  const Configuration & configuration,
280  const InstanceBlock & IB,
281  Model & Iiter,
282  const FactoryBlock & factory_block)
283  : DeviceInstance(IB, configuration.getInstanceParameters(), factory_block),
284  model_(Iiter),
285  temp(getDeviceOptions().temp.getImmutableValue<double>()),
286  outputStateVarsFlag( false ),
287  maxVoltageDrop(1.0e-10)
288 {
289  scalingRHS = 1.0;
290  numExtVars = 2;
291  numIntVars = 3;
292 
293  numStateVars = 2;
294  setNumStoreVars(3);
295 
296  tempGiven = false;
297 
298  const int ibev = IB.numExtVars;
299  const int ibiv = IB.numIntVars;
300 
301  // Set params to constant default values:
302  setDefaultParams ();
303 
304  // Set params according to instance line and constant defaults from metadata:
305  setParams (IB.params);
306 
307  // now load the instance data vector
308  for( int i=0; i<inductorNames.size(); ++i )
309  {
310  InductorInstanceData * inductorData = new InductorInstanceData();
311  inductorData->name = inductorNames[i];
312  inductorData->L = inductorInductances[i];
313  inductorData->baseL = inductorInductances[i];
314  inductorData->ICGiven = false;
315  inductorData->inductorCurrentOffsets.resize( inductorNames.size() );
316 
317  instanceData.push_back( inductorData );
318  }
319  numInductors = instanceData.size();
320 
321  // set up the device connectivity map
322  // each simple inductor in this mutual inductor
323  // is maked as a connection (given a common, non-zero
324  // value in devConMap)
325  devConMap.resize(2*numInductors);
326  for(int i=0; i<numInductors; i++)
327  {
328  devConMap[i] = devConMap[i+1] = (i+1);
329  }
330 
331  mEquInductorOffsets.resize( numInductors );
332  rEquInductorOffsets.resize( numInductors );
333  inductorCurrents.resize( numInductors );
334  inductanceVals.resize( numInductors );
335  LOI.resize( numInductors );
336  LO.resize( numInductors );
337  for( int i=0; i<numInductors; ++i)
338  {
339  LO[i].resize( numInductors );
340  }
341 
342  // set up the device connectivity map
343  // each simple inductor in this mutual inductor
344  // is maked as a connection (given a common, non-zero
345  // value in devConMap)
346  devConMap.resize(2*numInductors);
347  for(int i=0; i<numInductors; i++)
348  {
349  devConMap[i] = devConMap[i+1] = (i+1);
350  }
351 
353 
354  // Calculate any parameters specified as expressions:
356 
357  // calculate dependent (ie computed) params and check for errors:
358  processParams ();
359 
360  // if the user has requested output of the state variables M, H and R
361  // then open a file for that output.
362  if( model_.outputStateVars > 0 )
363  {
364  outputStateVarsFlag = true;
365  std::string filename( "Inductor_" );
366  filename.append( getName().getEncodedName() );
367  filename.append( ".dat" );
368  // convert any embeded ':' or '%' characters to '_'
369  replace( filename.begin(), filename.end(), '%', '_' );
370  replace( filename.begin(), filename.end(), ':', '_' );
371 
372  outputFileStreamPtr = rcp( new std::ofstream() );
373  outputFileStreamPtr->open( filename.c_str() );
374  if( !(*outputFileStreamPtr) )
375  {
376  UserError(*this) << "Could not open file " << filename << " for output of state variables";
377  }
378  else {
379  (*outputFileStreamPtr).setf(std::ios::scientific, std::ios::floatfield );
380  (*outputFileStreamPtr).width(20);
381  (*outputFileStreamPtr).precision(12);
382  }
383  }
384 
385  // size some vectors needed in loadDAEdFdx
386  dHe_dI.resize( numInductors );
387  dManp_dI.resize( numInductors );
388  ddelM_dI.resize( numInductors );
389  dMirrp_dI.resize( numInductors );
390  dP_dI.resize( numInductors );
391 
392  // Check if the magnetic moment equation should be dropped form the
393  // system of equations because domain flexing is essentially zero
394  if( model_.C <= model_.CLim )
395  {
396  model_.includeMEquation=false;
397  }
398 
399  // update internal/external/state variable counts
401  int numExtraEquations = 1; // Allways use the R equation
403  {
404  numExtraEquations++; // Also add the M equation
405  }
406  numIntVars = numInductors + numExtraEquations;
407  numStateVars = 5; // extra state variables for M and dM/dt and R, H & B
408  //numStateVars += 2*numInductors; individual inductors no longer need state / store space
409 
410  // set up the jacobian stamp
411  // for an individual inductor with the two interal variables would be
412  //
413  // V1 V2 Ib
414  // kcl1 1
415  // kcl2 -1
416  // branch 1 -1 L/dt
417  //
418  // for a collection of these, the internal variable, branch equations,
419  // must be at the end of a given stamp row as well as the internal
420  // vars for M and R in this non-linear version.
421  //
422  // So for N inductors the samp is:
423  //
424  // V1 V2 V3 V4 ... V2N I1 I2 ... IN M R
425  // kcl1 1
426  // kcl2 -1
427  // kcl3 1
428  // kcl4 -1
429  // branch1 1 -1 L/dt c ... c x
430  // branch2 x x 1 -1 c L/dt ... c x
431  // M equ x x x x ... x x x
432  // R equ x x ... x x
433  //
434  // where "c" is an induced current change and "x" are
435  // values which must be computed.
436 
437  jacStamp.resize( 2 * numInductors + numIntVars);
438 
439  for( int i=0; i< numInductors; ++i )
440  {
441  //
442  // allocate space
443  //
444  // kcl V+ node
445  jacStamp[2*i].resize(1);
446  // kcl V- node
447  jacStamp[2*i+1].resize(1);
448  // branch node -- every branch needs to access the first
449  // inductor's V+ and V-, so they all contribute there
450  if( i == 0 )
451  {
453  {
454  jacStamp[2*numInductors].resize(numInductors + 3);
455  }
456  else
457  {
458  jacStamp[2*numInductors].resize(numInductors + 2);
459  }
460  }
461  else
462  {
464  {
465  jacStamp[2*numInductors + i].resize(numInductors + 5);
466  }
467  else
468  {
469  jacStamp[2*numInductors + i].resize(numInductors + 4);
470  }
471  }
472 
473  //
474  // fill in dependency
475  //
476  // kcl V+ node
477  jacStamp[2*i ][0] = 2*numInductors + i;
478  // kcl V- node
479  jacStamp[2*i+1][0] = 2*numInductors + i;
480 
481  if( i==0 )
482  {
483  jacStamp[2*numInductors ][0] = 0;
484  jacStamp[2*numInductors ][1] = 1;
485  for( int j=0; j<numInductors; ++j )
486  {
487  jacStamp[2*numInductors][j+2] = 2*numInductors + j;
488  }
490  {
491  jacStamp[2*numInductors][numInductors+2] = 3*numInductors;
492  }
493  }
494  else
495  {
496  jacStamp[2*numInductors + i][0] = 0;
497  jacStamp[2*numInductors + i][1] = 1;
498  jacStamp[2*numInductors + i][2] = 2*i;
499  jacStamp[2*numInductors + i][3] = 2*i + 1;
500  for( int j=0; j<numInductors; ++j )
501  {
502  jacStamp[2*numInductors + i][j+4] = 2*numInductors + j;
503  }
504  // only add this term if M equation is in use
506  {
507  jacStamp[2*numInductors + i][numInductors+4] = 3*numInductors;
508  }
509  }
510  }
511 
512  int rOffset=0;
514  {
515  // now the M equation
516  jacStamp[ 3*numInductors ].resize(numInductors + 4);
517  // M offsets for V+ and V- on first inductor
518  jacStamp[ 3*numInductors ][0] = 0;
519  jacStamp[ 3*numInductors ][1] = 1;
520  // M offsets to each inductor's branch equ.
521  for(int i=0; i<numInductors; ++i)
522  {
523  jacStamp[ 3*numInductors ][i+2]=2*numInductors+i;
524  }
525  // M offsets to M and R
526  jacStamp[ 3*numInductors ][numInductors + 2] = 3*numInductors;
527  jacStamp[ 3*numInductors ][numInductors + 3] = 3*numInductors + 1;
528  rOffset=1;
529  }
530 
531  // and the R equations
532  jacStamp[ 3*numInductors + rOffset].resize(numInductors + 1);
533 
534  // R offsets to each inductor's branch equ.
535  for(int i=0; i<numInductors; ++i)
536  {
537  jacStamp[ 3*numInductors + rOffset ][i]=2*numInductors+i;
538  }
539 
540  // R offsets to R
541  jacStamp[ 3*numInductors + rOffset][numInductors ] = 3*numInductors + rOffset;
542 
543  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS))
544  {
545  Xyce::dout() << "Instance::Instance----------" << std::endl;
546  Xyce::dout() << "numExtVars = " << numExtVars << ", " << ibev << std::endl
547  << "numIntVars = " << numIntVars << ", " << ibiv << std::endl
548  << "numStateVars = " << numStateVars << std::endl
549  << "numInductors = " << numInductors << std::endl
550  << "jacStamp = " << std::endl;
551  for( int i = 0; i<jacStamp.size(); ++i )
552  {
553  Xyce::dout() << "jacStamp[ " << i << " ] = { ";
554  for( int j=0; j<jacStamp[i].size(); ++j)
555  {
556  Xyce::dout() << jacStamp[i][j];
557  if( j != ( jacStamp[i].size() -1 ) )
558  {
559  Xyce::dout() << ", ";
560  }
561  }
562  Xyce::dout() << " }" << std::endl;
563  }
564  }
565 }
566 
567 //-----------------------------------------------------------------------------
568 // Function : Instance::~Instance
569 // Purpose : destructor
570 // Special Notes :
571 // Scope : public
572 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
573 // Creation Date : 03/21/2005
574 //-----------------------------------------------------------------------------
576 {
577  // Close output file if we opened one
578  if( outputStateVarsFlag && outputFileStreamPtr.get() && outputFileStreamPtr->is_open() )
579  {
580  outputFileStreamPtr->close();
581  }
582 
583  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
584  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
585  for ( ; currentInductor != endInductor ; ++currentInductor)
586  {
587  if (*currentInductor != NULL)
588  {
589  delete *currentInductor;
590  *currentInductor = NULL;
591  }
592  }
593  instanceData.clear();
594 }
595 //-----------------------------------------------------------------------------
596 // Function : Instance::processParams
597 // Purpose :
598 // Special Notes :
599 // Scope : public
600 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
601 // Creation Date : 03/21/2005
602 //-----------------------------------------------------------------------------
604 {
605  // now set the temperature related stuff.
606  if (tempGiven)
607  {
609  }
610 
611  return true;
612 }
613 
614 
615 //-----------------------------------------------------------------------------
616 // Function : Instance::registerLIDs
617 // Purpose :
618 // Special Notes :
619 // Scope : public
620 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
621 // Creation Date : 03/21/2005
622 //-----------------------------------------------------------------------------
623 void Instance::registerLIDs(const std::vector<int> & intLIDVecRef,
624  const std::vector<int> & extLIDVecRef)
625 {
626  AssertLIDs(intLIDVecRef.size() == numIntVars);
627  AssertLIDs(extLIDVecRef.size() == numExtVars);
628 
629  // copy over the global ID lists.
630  intLIDVec = intLIDVecRef;
631  extLIDVec = extLIDVecRef;
632 
633  // Now use these lists to obtain the indices into the
634  // linear algebra entities. This assumes an order.
635  // For the matrix indices, first do the rows.
636  // get the current values of the inductances and currentOffsets
637  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
638  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
639  int i = 0;
640  int j = 0;
641  while( currentInductor != endInductor )
642  {
643  (*currentInductor)->li_Pos = extLIDVec[ i++ ];
644  (*currentInductor)->li_Neg = extLIDVec[ i++ ];
645  (*currentInductor)->li_Branch = intLIDVec[ j++ ];
646  currentInductor++;
647  }
648 
649  // now get the M and R local id's
651  {
652  li_MagVar = intLIDVec[ j++ ];
653  }
654  li_RVar = intLIDVec[ j++ ];
655 
656  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS))
657  {
658  Xyce::dout() << "Instance::registerLIDs------------------------" << std::endl;
659  currentInductor = instanceData.begin();
660  i=0;
661  while( currentInductor != endInductor )
662  {
663  Xyce::dout() << "Inductor [ " << i++ << " ] "
664  << " li_Pos = " << (*currentInductor)->li_Pos
665  << " li_Neg = " << (*currentInductor)->li_Neg
666  << " li_Branch = " << (*currentInductor)->li_Branch << std::endl;
667  currentInductor++;
668  }
669  Xyce::dout() << " li_MagVar = " << li_MagVar << std::endl
670  << " li_RVar = " << li_RVar << std::endl;
671  }
672 }
673 
674 //-----------------------------------------------------------------------------
675 // Function : Instance::loadNodeSymbols
676 // Purpose :
677 // Special Notes :
678 // Scope : public
679 // Creator : Eric R. Keiter, SNL, Parallel Computational Sciences
680 // Creation Date : 05/13/05
681 //-----------------------------------------------------------------------------
682 void Instance::loadNodeSymbols(Util::SymbolTable &symbol_table) const
683 {
684  for (std::vector<InductorInstanceData *>::const_iterator it = instanceData.begin(), end = instanceData.end(); it != end; ++it ) {
685  std::string branchInductorName = getName().getSubcircuitName() + ":" + (*it)->name;
686  if( getName().getSubcircuitName() == "" )
687  branchInductorName = (*it)->name;
688  InstanceName bInductorIName = InstanceName( branchInductorName );
689  std::string encodedName = spiceInternalName( bInductorIName, "branch");
690  addInternalNode(symbol_table, (*it)->li_Branch, getName(), (*it)->name + "_branch");
691  addInternalNode(symbol_table, (*it)->li_Branch, encodedName);
692  }
693 
695  {
696  addInternalNode(symbol_table, li_MagVar, getName(), "m");
697  addInternalNode(symbol_table, li_MagVar, getName().getEncodedName() + "_m");
698  }
699  addInternalNode(symbol_table, li_RVar, getName(), "r");
700  addInternalNode(symbol_table, li_RVar, getName().getEncodedName() + "_r");
701 
702  addStateNode(symbol_table, li_MagVarState, getName(), "ms");
703  addStateNode(symbol_table, li_MagVarDerivState, getName(), "dmdt");
704 
705  addStoreNode(symbol_table, li_RVarStore, getName().getEncodedName() + "_r");
706  addStoreNode(symbol_table, li_BVarStore, getName().getEncodedName() + "_b");
707  addStoreNode(symbol_table, li_HVarStore, getName().getEncodedName() + "_h");
708 }
709 
710 //-----------------------------------------------------------------------------
711 // Function : Instance::registerStateLIDs
712 // Purpose :
713 // Special Notes :
714 // Scope : public
715 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
716 // Creation Date : 03/21/2005
717 //-----------------------------------------------------------------------------
718 void Instance::registerStateLIDs( const std::vector<int> & staLIDVecRef )
719 {
720  AssertLIDs(staLIDVecRef.size() == numStateVars);
721 
722  // copy over the global ID lists.
723  staLIDVec = staLIDVecRef;
724  int i = 0;
725 
726  li_MagVarState = staLIDVec[i++];
727  li_MagVarDerivState = staLIDVec[i++];
728 
729  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS))
730  {
731  Xyce::dout() << "Instance::registerStateLIDs-------------------" << std::endl;
732 
733  Xyce::dout() << "li_MagVarState = " << li_MagVarState << std::endl
734  << "li_MagVarDerivState = " << li_MagVarDerivState << std::endl
735  ;
736  }
737 }
738 
739 //-----------------------------------------------------------------------------
740 // Function : Instance::registerStoreLIDs
741 // Purpose :
742 // Special Notes :
743 // Scope : public
744 // Creator : Richard Schiek, SNL
745 // Creation Date : 8/17/2012
746 //-----------------------------------------------------------------------------
747 void Instance::registerStoreLIDs(const std::vector<int> & stoLIDVecRef )
748 {
749  AssertLIDs(stoLIDVecRef.size() == getNumStoreVars());
750 
751  // copy over the global ID lists.
752  stoLIDVec = stoLIDVecRef;
753 
754  li_RVarStore = stoLIDVec[0];
755  li_BVarStore = stoLIDVec[1];
756  li_HVarStore = stoLIDVec[2];
757 }
758 
759 //-----------------------------------------------------------------------------
760 // Function : Instance::jacobianStamp
761 // Purpose :
762 // Special Notes :
763 // Scope : public
764 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
765 // Creation Date : 03/21/2005
766 //-----------------------------------------------------------------------------
767 const std::vector< std::vector<int> > & Instance::jacobianStamp() const
768 {
769  return jacStamp;
770 }
771 
772 //-----------------------------------------------------------------------------
773 // Function : Instance::registerJacLIDs
774 // Purpose :
775 // Special Notes :
776 // Scope : public
777 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
778 // Creation Date : 03/21/2005
779 //-----------------------------------------------------------------------------
780 void Instance::registerJacLIDs( const std::vector< std::vector<int> > & jacLIDVec )
781 {
782  DeviceInstance::registerJacLIDs( jacLIDVec );
783  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS))
784  {
785  Xyce::dout() << "Instance::registerJacLIDs ----------------------------" << std::endl;
786 
787  Xyce::dout() << "jacLIDVec = " << std::endl;
788  for( int i = 0; i<jacStamp.size(); ++i )
789  {
790  Xyce::dout() << "jacLIDVec[ " << i << " ] = { ";
791  for( int j=0; j<jacLIDVec[i].size(); ++j)
792  {
793  Xyce::dout() << jacLIDVec[i][j];
794  if( j != ( jacLIDVec[i].size() -1 ) )
795  {
796  Xyce::dout() << ", ";
797  }
798  }
799  Xyce::dout() << " }" << std::endl;
800  }
801  }
802 
803  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
804  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
805  // int numInductors = instanceData.size(); // don't need this as it's defined at class level
806  int i = 0;
807  while( currentInductor != endInductor )
808  {
809  (*currentInductor)->APosEquBraVarOffset = jacLIDVec[ 2*i ][ 0 ];
810  (*currentInductor)->ANegEquBraVarOffset = jacLIDVec[ 2*i + 1 ][ 0 ];
811  (*currentInductor)->vPosOffset = jacLIDVec[ 2*numInductors + i ][ 0 ];
812  (*currentInductor)->vNegOffset = jacLIDVec[ 2*numInductors + i ][ 1 ];
813  int extraOffset = 2;
814  if( i == 0)
815  {
816  extraOffset = 0;
817  }
818  (*currentInductor)->ABraEquPosNodeOffset = jacLIDVec[ 2*numInductors + i ][ 0 + extraOffset ];
819  (*currentInductor)->ABraEquNegNodeOffset = jacLIDVec[ 2*numInductors + i ][ 1 + extraOffset ];
820  for( int j=0; j<numInductors; ++j )
821  {
822  if( i == j )
823  {
824  (*currentInductor)->ABraEquBraVarOffset = jacLIDVec[ 2*numInductors + i ][ j + 2 + extraOffset ];
825  }
826  (*currentInductor)->inductorCurrentOffsets[ j ] = jacLIDVec[ 2*numInductors + i ][ j + 2 + extraOffset ];
827  }
829  {
830  (*currentInductor)->magOffset = jacLIDVec[ 2*numInductors + i ][ numInductors + 2 + extraOffset ];
831  }
832  currentInductor++;
833  i++;
834  }
835 
836  int rOffset=0;
838  {
839  // now get the M equation offsets
840  mEquVPosOffset = jacLIDVec[ 3*numInductors ][0];
841  mEquVNegOffset = jacLIDVec[ 3*numInductors ][1];
842  for( i=0; i<numInductors; ++i )
843  {
844  mEquInductorOffsets[i] = jacLIDVec[ 3*numInductors ][ i + 2];
845  }
846  mEquMOffset = jacLIDVec[ 3*numInductors ][ numInductors + 2 ];
847  mEquROffset = jacLIDVec[ 3*numInductors ][ numInductors + 3 ];
848 
849  rOffset=1;
850  }
851 
852  // now get the R equation offsets
853  for( i=0; i<numInductors; ++i )
854  {
855  rEquInductorOffsets[i] = jacLIDVec[ 3*numInductors + rOffset ][ i ];
856  }
857  rEquROffset = jacLIDVec[ 3*numInductors + rOffset ][ numInductors ];
858 
859  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS))
860  {
861  currentInductor = instanceData.begin();
862  i=0;
863  while( currentInductor != endInductor )
864  {
865  Xyce::dout() << "Inductor [ " << i << " ] " << (*currentInductor)->name << std::endl
866  << " APosEquBraVarOffset = " << (*currentInductor)->APosEquBraVarOffset << std::endl
867  << " ANegEquBraVarOffset = " << (*currentInductor)->ANegEquBraVarOffset << std::endl
868  << " vPosOffset = " << (*currentInductor)->vPosOffset << std::endl
869  << " vNegOffset = " << (*currentInductor)->vNegOffset << std::endl
870  << " ABraEquPosNodeOffset = " << (*currentInductor)->ABraEquPosNodeOffset << std::endl
871  << " ABraEquNegNodeOffset = " << (*currentInductor)->ABraEquNegNodeOffset << std::endl
872  << " ABraEquBraVarOffset = " << (*currentInductor)->ABraEquBraVarOffset << std::endl
873  << " magOffset = " << (*currentInductor)->magOffset << std::endl;
874  Xyce::dout() << "\tInductor branch offsets = { ";
875  for( int j=0; j<numInductors ; ++j )
876  {
877  Xyce::dout() << (*currentInductor)->inductorCurrentOffsets[ j ] << ", ";
878  }
879  Xyce::dout() << "} " << std::endl;
880  i++;
881  currentInductor++;
882  }
883 
884  Xyce::dout() << "mEquVPosOffset = " << mEquVPosOffset << "\tmEquVNegOffset = " << mEquVNegOffset << std::endl;
885  Xyce::dout() << "mEquInductorOffsets = ";
886  for(i=0;i<numInductors; ++i)
887  {
888  Xyce::dout() << mEquInductorOffsets[i] << ", ";
889  }
890  Xyce::dout() << std::endl
891  << "mEquMOffset = " << mEquMOffset << "\tmEquROffset = " << mEquROffset << std::endl;
892 
893  Xyce::dout() << "rEquInductorOffsets = ";
894  for(i=0;i<numInductors; ++i)
895  {
896  Xyce::dout() << rEquInductorOffsets[i] << ", ";
897  }
898  Xyce::dout() << std::endl
899  << "rEquROffset = " << rEquROffset << std::endl;
900  }
901 }
902 
903 //-----------------------------------------------------------------------------
904 // Function : Instance::updateTemperature
905 // Purpose :
906 // Special Notes :
907 // Scope : public
908 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
909 // Creation Date : 03/21/2005
910 //-----------------------------------------------------------------------------
911 bool Instance::updateTemperature ( const double & temp)
912 {
913  bool bsuccess = true;
914 
915  // current temp difference from reference temp.
916  double difference = temp - model_.tnom;
917 
918  std::vector< InductorInstanceData* >::iterator currentData = instanceData.begin();
919  while( currentData != instanceData.end() )
920  {
921  double factor = 1.0 + (model_.tempCoeff1)*difference +
922  (model_.tempCoeff2)*difference*difference;
923  (*currentData)->L = ((*currentData)->baseL)*factor;
924  currentData++;
925  }
926 
927  // now that the inductances have changed we need to update the matrix.
929 
930  return bsuccess;
931 }
932 
933 //-----------------------------------------------------------------------------
934 // Function : Instance::updateIntermediateVars
935 // Purpose : updates a set of common variables used by RHS and jacobian
936 // Special Notes :
937 // Scope : public
938 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
939 // Creation Date : 03/21/2005
940 //-----------------------------------------------------------------------------
942 {
943  bool bsuccess = true;
944 
945  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
946  {
947  Xyce::dout() << "Instance::updateIntermediateVars " << std::endl;
948  }
949 
950  Linear::Vector & solVector = *(extData.nextSolVectorPtr);
951  Linear::Vector & staVector = *(extData.nextStaVectorPtr);
952 
953  // some parameters in the model class that we will use often
954  const double A = model_.A;
955  const double Alpha = model_.Alpha;
956  const double Area = model_.Area;
957  const double BetaH = model_.BetaH;
958  const double BetaM = model_.BetaM;
959  const double C = model_.C;
960  const double DeltaVScaling = model_.DeltaVScaling;
961  const double Gap = model_.Gap;
962  const double Ms = model_.Ms;
963  const double Kirr = model_.Kirr;
964  const double Path = model_.Path;
965 
966  const double mVarScaling = model_.mVarScaling;
967  const double rVarScaling = model_.rVarScaling;
968  const double mEqScaling = model_.mEqScaling;
969  const double rEqScaling = model_.rEqScaling;
970 
971  // calculate the voltage drop over the first inductor
972  // as this is needed later
973  double Vpos = solVector[(instanceData[0])->li_Pos];
974  double Vneg = solVector[(instanceData[0])->li_Neg];
975 
976  // voltage drop over first inductor.
977  double voltageDrop= Vpos - Vneg;
978 
979  // only update maxVoltageDrop when system has converged or we may
980  // get wildly wrong values.
981  Linear::Vector & lastSolVector = *(extData.currSolVectorPtr);
982  double lastVoltageDrop = lastSolVector[(instanceData[0])->li_Pos] - lastSolVector[(instanceData[0])->li_Neg];
983  if ( (getSolverState().newtonIter == 0) && (fabs(lastVoltageDrop) > maxVoltageDrop) )
984  {
985  maxVoltageDrop=fabs(lastVoltageDrop);
986  }
987 
988  // approximate the sgn( voltageDrop ) with
989  // tanh ( scalefactor * voltageDrop / maxVoltageDrop )
991  {
992  qV = DeltaVScaling * voltageDrop;
993  }
994  else
995  {
996  qV = DeltaVScaling * voltageDrop / maxVoltageDrop;
997  }
998 
999  double tanh_qV = 0.0;
1000 
1001  if ( (fabs(qV) < CONSTTANH_THRESH) )
1002  {
1003  tanh_qV = tanh(qV);
1004  }
1005  else if (qV < 0.0)
1006  {
1007  tanh_qV = -1.0;
1008  }
1009  else
1010  {
1011  tanh_qV = 1.0;
1012  }
1013 
1014  Happ = 0.0;
1015  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1016  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1017  int il=0;
1018  while( currentInductor != endInductor )
1019  {
1020  Happ += solVector[(*currentInductor)->li_Branch] * inductanceVals[ il ];
1021  il++;
1022  currentInductor++;
1023  }
1024  Happ /= Path;
1025 
1026  double latestMag=0.0;
1027  if( model_.includeMEquation )
1028  {
1029  latestMag = mVarScaling * solVector[ li_MagVar ];
1030  }
1031  else
1032  {
1033  latestMag = mVarScaling * staVector[ li_MagVarState ];
1034  }
1035 
1036  if( model_.factorMS )
1037  {
1038  latestMag *= Ms;
1039  }
1040 
1041  double H = Happ - (Gap / Path) * latestMag;
1042 
1043  He = H + Alpha * latestMag;
1044 
1045  Heo = BetaH*A;
1046 
1047  // terms that come up frequently
1048  const double gap_path = Gap / Path;
1049  const double He2 = He*He;
1050  const double Heo2 = Heo*Heo;
1051  const double sq_Heo2He2 = sqrt(Heo2 + He2);
1052 
1053  delM0 = model_.BetaM * Ms;
1054  double Man = Ms * He / ( A + sq_Heo2He2 );
1055  delM = Man - latestMag;
1056 
1057  // terms that come up frequently
1058  const double delM2 = delM*delM;
1059  const double delM02 = delM0*delM0;
1060  const double sq_delM02delM2 = sqrt( delM02 + delM2 );
1061 
1062  if( model_.factorMS )
1063  {
1064  Mirrp = (delM * tanh_qV + sq_delM02delM2 ) / (2*( Kirr- Alpha * sq_delM02delM2));
1065  Manp = Ms * (A + Heo2/sq_Heo2He2) / pow(A + sq_Heo2He2, 2.0);
1066  P = ( C * Manp + (1 - C) * Mirrp) / ((1 + (gap_path - Alpha) * C * Manp + gap_path * (1-C) * Mirrp)*Ms);
1067  }
1068  else
1069  {
1070  Mirrp = (delM * tanh_qV + sq_delM02delM2 ) / (2*( Kirr- Alpha * sq_delM02delM2));
1071  Manp = Ms * (A + Heo2/sq_Heo2He2) / pow(A + sq_Heo2He2, 2.0);
1072  P = ( C * Manp + (1 - C) * Mirrp) / (1 + (gap_path - Alpha) * C * Manp + gap_path * (1-C) * Mirrp);
1073  }
1074 
1075  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1076  {
1077  Xyce::dout() << "\tA = " << A << std::endl
1078  << "\tArea = " << Area << std::endl
1079  << "\tPath = " << Path << std::endl
1080  << "\tGap = " << Gap << std::endl
1081  << "\tC = " << C << std::endl
1082  << "\tVpos = " << Vpos << std::endl
1083  << "\tVneg = " << Vneg << std::endl
1084  << "\tvoltageDrop = " << voltageDrop << std::endl
1085  << "\tqV = " << qV << std::endl
1086  << "\tdelM0 = " << delM0 << std::endl
1087  << "\tHapp = " << Happ
1088  << "\tlatestMag = " << latestMag
1089  << "\tlatestR = " << rVarScaling * solVector[ li_RVar ] << std::endl
1090  << "\tHe = " << He << std::endl
1091  << "\tH = " << H << std::endl
1092  << "\tHeo = " << Heo << std::endl
1093  << "\tMan = " << Man << std::endl
1094  << "\tdelM = " << delM << std::endl
1095  << "\tMirrp = " << Mirrp << std::endl
1096  << "\tManp = " << Manp << std::endl
1097  << "\tP = " << P << std::endl
1098  << "\tgetSolverState().newtonIter = " << getSolverState().newtonIter << std::endl
1099  << std::endl;
1100  }
1101 
1102  // now calculate important derivative quantities
1103 
1104  double dHe_dM = ((Alpha - gap_path) * mVarScaling);
1105 
1106  double dManp_dM = ( -Ms * He / (pow(A + sq_Heo2He2, 2.0)*sq_Heo2He2)) *
1107  ( (Heo2 / (Heo2 + He2)) + (2.0*(A + Heo2 / sq_Heo2He2)/(A+sq_Heo2He2)) ) * dHe_dM;
1108 
1109  double ddelM_dM = ( dHe_dM*Ms/(A + sq_Heo2He2) ) * (1.0 - He2 / ((A + sq_Heo2He2)*sq_Heo2He2)) - mVarScaling;
1110 
1111  double dMirrp_dM = (1.0/(2.0*(Kirr - Alpha*sq_delM02delM2))) *
1112  (tanh_qV + delM/sq_delM02delM2 +
1113  (2.0*Alpha*delM*(delM*tanh_qV + sq_delM02delM2)
1114  /(2.0*(Kirr-Alpha*sq_delM02delM2)*sq_delM02delM2))) * ddelM_dM;
1115 
1116  double dP_Denom=0.0;
1117  if( model_.factorMS )
1118  {
1119  dP_Denom = 1.0 + (gap_path - Alpha)*C*Manp + gap_path * (1.0-C) * Mirrp;
1120 
1121  dP_dM = (1.0/dP_Denom) * (C * dManp_dM + (1.0-C) * dMirrp_dM) -
1122  ( (C*Manp + (1.0-C)*Mirrp)/pow(dP_Denom,2.0) ) *
1123  ( (gap_path - Alpha)*C*dManp_dM + gap_path*(1.0-C)*dMirrp_dM );
1124  dP_dM /= Ms;
1125  }
1126  else
1127  {
1128  dP_Denom = 1.0 + (gap_path - Alpha)*C*Manp + gap_path * (1.0-C) * Mirrp;
1129 
1130  dP_dM = (1.0/dP_Denom) * (C * dManp_dM + (1.0-C) * dMirrp_dM) -
1131  ( (C*Manp + (1.0-C)*Mirrp)/pow(dP_Denom,2.0) ) *
1132  ( (gap_path - Alpha)*C*dManp_dM + gap_path*(1.0-C)*dMirrp_dM );
1133  }
1134 
1135  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1136  {
1137  Xyce::dout() << "\tA = " << A << std::endl
1138  << "\tAlpha = " << Alpha << std::endl
1139  << "\tC = " << C << std::endl
1140  << "\tGap = " << Gap << std::endl
1141  << "\tMs = " << Ms << std::endl
1142  << "\tKirr = " << Kirr << std::endl
1143  << "\tPath = " << Path << std::endl
1144  << "\tHe2 = " << He2 << std::endl
1145  << "\tHeo2 = " << Heo2 << std::endl
1146  << "\tdelM2 = " << delM2 << std::endl
1147  << "\tdelM02 = " << delM02 << std::endl
1148  << "\tdHe_dM = " << dHe_dM << std::endl
1149  << "\tdManp_dM = " << dManp_dM << std::endl
1150  << "\tddelM_dM = " << ddelM_dM << std::endl
1151  << "\tdMirrp_dM = " << dMirrp_dM << std::endl
1152  << "\tdP_dM = " << dP_dM << std::endl
1153  << "\tdenom 1+(1-lg/lt)P = " << (1+(1-Gap/Path)*P) << std::endl;
1154  }
1155 
1156  // % Now find (dP/dI_i): (this is nearly identical to dP/dM)
1157  currentInductor = instanceData.begin();
1158 
1159  for( int i=0; i<numInductors; ++i )
1160  {
1161 
1162  dHe_dI[ i ] = inductanceVals[ i ] / Path;
1163  dManp_dI[i] = ( -Ms * He / (pow(A + sq_Heo2He2, 2.0)*sq_Heo2He2)) *
1164  ( (Heo2 / (Heo2 + He2)) + (2.0*(A + Heo2 / sq_Heo2He2)/(A+sq_Heo2He2)) ) * dHe_dI[i];
1165  ddelM_dI[i] = (Ms / (A + sq_Heo2He2)) * (1.0 - He2/((A + sq_Heo2He2)*sq_Heo2He2)) * dHe_dI[i];
1166  dMirrp_dI[i] = (1.0/(2.0*(Kirr - Alpha*sq_delM02delM2))) *
1167  (tanh_qV + delM/sq_delM02delM2 +
1168  (2.0*Alpha*delM*(delM*tanh_qV +
1169  sq_delM02delM2)/(2.0*(Kirr-Alpha*sq_delM02delM2)*sq_delM02delM2))) * ddelM_dI[i];
1170  dP_dI[i] = (1.0/dP_Denom) * (C * dManp_dI[i] + (1.0-C) * dMirrp_dI[i]) -
1171  ( (C*Manp + (1.0-C)*Mirrp)/pow(dP_Denom,2.0) ) *
1172  ( (gap_path - Alpha)*C*dManp_dI[i] + gap_path*(1.0-C)*dMirrp_dI[i] );
1173 
1174  if( model_.factorMS )
1175  {
1176  dP_dI[i] /= Ms;
1177  }
1178 
1179  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1180  {
1181  Xyce::dout() << "\tdHe_dI[ " << i << " ] =" << dHe_dI[ i ] << std::endl
1182  << "\tdManp_dI[ " << i << " ] = " << dManp_dI[i] << std::endl
1183  << "\tddelM_dI[ " << i << " ] = " << ddelM_dI[i] << std::endl
1184  << "\tMirrp_dI[ " << i << " ] = " << dMirrp_dI[i] << std::endl
1185  << "\tdP_dI[ " << i << " ] = " << dP_dI[i] << std::endl;
1186  }
1187  currentInductor++;
1188  }
1189 
1190  // Now find (dP/dV_1):
1191  double dMirrp_dVp = (delM * DeltaVScaling * (1.0-pow(tanh_qV,2.0))) /
1192  (2.0 * (Kirr - Alpha * sq_delM02delM2));
1193  double dMirrp_dVn = -dMirrp_dVp;
1194 
1195  dP_dVp = (1.0/dP_Denom) * ((1.0-C) * dMirrp_dVp) -
1196  ( (C*Manp + (1.0-C)*Mirrp)/pow(dP_Denom,2.0) ) * ( gap_path*(1.0-C)*dMirrp_dVp );
1197 
1198  if( model_.factorMS )
1199  {
1200  dP_dVp /= Ms;
1201  }
1202 
1203  dP_dVn = -dP_dVp;
1204 
1205  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1206  {
1207  Xyce::dout() << "\tdMirrp_dVp = " << dMirrp_dVp << std::endl
1208  << "\tdMirrp_dVn = " << dMirrp_dVn << std::endl
1209  << "\tdP_dVp = " << dP_dVp << std::endl
1210  << "\tdP_dVn = " << dP_dVn << std::endl;
1211  }
1212 
1213  return true;
1214 }
1215 
1216 //-----------------------------------------------------------------------------
1217 // Function : Instance::updateInductanceMatrix()
1218 // Purpose : A matrix of inductances is used often enough that it
1219 // calculated and stored as a member variable here
1220 // If and inductance ever changes say from updating
1221 // the temperature or a parameter udpate, then this
1222 // routine must be called again.
1223 // Special Notes :
1224 // Scope : public
1225 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1226 // Creation Date : 03/21/2005
1227 //-----------------------------------------------------------------------------
1229 {
1230  std::vector< InductorInstanceData* >::iterator
1231  currentInductor = instanceData.begin();
1232  std::vector< InductorInstanceData* >::iterator
1233  endInductor = instanceData.end();
1234 
1235  // collec the inductances
1236  int i=0;
1237  while( currentInductor != endInductor )
1238  {
1239  inductanceVals[ i ] = (*currentInductor)->L;
1240  i++;
1241  currentInductor++;
1242  }
1243 
1244  double Area = model_.Area;
1245  double Path = model_.Path;
1246 
1247  // compute the inductance matrix
1248  for( i=0; i<numInductors; ++i)
1249  {
1250  for( int j=0; j<numInductors; ++j)
1251  {
1252  // 4.0e-7 * M_PI is a magnetic constant, the permeability of free space [Henries/m]
1253  LO[i][j] = mutualCup * 4.0e-7 * M_PI * (Area / Path) * inductanceVals[i] * inductanceVals[j];
1254  }
1255  }
1256 
1257 }
1258 
1259 
1260 //-----------------------------------------------------------------------------
1261 // Function : Instance::updatePrimaryState
1262 // Purpose :
1263 // Special Notes :
1264 // Scope : public
1265 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1266 // Creation Date : 03/21/2005
1267 //-----------------------------------------------------------------------------
1269 {
1270  bool bsuccess = true;
1271  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1272  {
1273  Xyce::dout() << "Instance::updatePrimaryState---------------" << std::endl
1274  << "\tname = " << getName() << std::endl;
1275  }
1276 
1278 
1279  Linear::Vector & solVector = *(extData.nextSolVectorPtr);
1280  Linear::Vector & staVector = *(extData.nextStaVectorPtr);
1281  Linear::Vector & stoVector = *(extData.nextStoVectorPtr);
1282  double mVarScaling = model_.mVarScaling;
1283  double rVarScaling = model_.rVarScaling;
1284 
1285  // place current values of mag, H and R in state vector
1286  // must unscale them as the rest of the class assumes
1287  // that these aren't scaled yet.
1288  double latestMag = 0.0;
1289  if( model_.includeMEquation )
1290  {
1291  staVector[ li_MagVarState ] = solVector[ li_MagVar ];
1292  latestMag = mVarScaling * solVector[ li_MagVar ];
1293  }
1294  else
1295  {
1296  latestMag = mVarScaling * staVector[ li_MagVarState ];
1297  }
1298  stoVector[ li_RVarStore ] = solVector[ li_RVar ];
1299 
1300  if( model_.factorMS )
1301  {
1302  latestMag *= model_.Ms;
1303  }
1304 
1305  // B and H are quantities that we can calculate from M and R. Store them in the state vector
1306  // for output if the user requests it.
1307  stoVector[ li_HVarStore ] = model_.HCgsFactor * (Happ - (model_.Gap / model_.Path) * latestMag);
1308  stoVector[ li_BVarStore ] = model_.BCgsFactor * (4.0e-7 * M_PI * (stoVector[ li_HVarStore ] + latestMag));
1309 
1310  return bsuccess;
1311 }
1312 
1313 //-----------------------------------------------------------------------------
1314 // Function : Instance::updateSecondaryState
1315 // Purpose :
1316 // Special Notes :
1317 // Scope : public
1318 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1319 // Creation Date : 03/21/2005
1320 //-----------------------------------------------------------------------------
1322 {
1323  bool bsuccess = true;
1324  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1325  {
1326  Xyce::dout() << "Instance::updateSecondaryState-------------" << std::endl
1327  << "\tname = " << getName() << std::endl;
1328  }
1329 
1330  Linear::Vector & staVector = *(extData.nextStaVectorPtr);
1331  Linear::Vector & staDerivVec = *(extData.nextStaDerivVectorPtr);
1332 
1333  // copy derivitive of Mag from result vector into state vector
1334  staVector[ li_MagVarDerivState ] = staDerivVec[ li_MagVarState ];
1335 
1336  return bsuccess;
1337 }
1338 
1339 //-----------------------------------------------------------------------------
1340 // Function : Instance::loadDAEQVector
1341 //
1342 // Purpose : Loads the Q-vector contributions for a single
1343 // instance.
1344 //
1345 // Special Notes : The "Q" vector is part of a standard DAE formalism in
1346 // which the system of equations is represented as:
1347 //
1348 // f(x) = dQ(x)/dt + F(x) - B(t) = 0
1349 //
1350 // Scope : public
1351 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1352 // Creation Date : 03/21/2005
1353 //-----------------------------------------------------------------------------
1355 {
1356  bool bsuccess = true;
1357  double mVarScaling = model_.mVarScaling;
1358  double rVarScaling = model_.rVarScaling;
1359  double mEqScaling = model_.mEqScaling;
1360  double rEqScaling = model_.rEqScaling;
1361 
1362  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1363  {
1364  Xyce::dout() << "Instance::loadDAEQVector------------------------" << std::endl
1365  << "\tname = " << getName() << std::endl;
1366  }
1367 
1368  Linear::Vector & staVector = *(extData.nextStaVectorPtr);
1369  Linear::Vector & solVector = *(extData.nextSolVectorPtr);
1370  double * qVec = extData.daeQVectorRawPtr;
1371 
1372  // update LOI -- the following product
1373  // I = column vector of currents
1374  // L = row vector of inductances
1375  // LO = matrix = mutualCup * sqrt( L' * L )
1376  // LOI = column vector = mutualCup * sqrt( L' * L ) * I
1377  // LOI[1] = mutualCup * sqrt(L[1]*L[1])*I[1]) +
1378  // mutualCup * sqrt(L[1]*L[2])*I[2]) + ...
1379  // mutualCup * sqrt(L[1]*L[n])*I[n])
1380 
1381  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1382  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1383  int i = 0;
1384  while( currentInductor != endInductor )
1385  {
1386  if( (getSolverState().dcopFlag) && (*currentInductor)->ICGiven == true )
1387  {
1388  inductorCurrents[ i ] = (*currentInductor)->IC;
1389  }
1390  else
1391  {
1392  inductorCurrents[ i ] = solVector[ (*currentInductor)->li_Branch ];
1393  }
1394  i++;
1395  currentInductor++;
1396  }
1397 
1398  for( i = 0; i < numInductors; ++i )
1399  {
1400  LOI[ i ] = 0;
1401  for( int j = 0; j < numInductors; ++j )
1402  {
1403  LOI[i] += LO[i][j] * inductorCurrents[j];
1404  }
1405  }
1406 
1407  // loop over each inductor and load it's Q vector components
1408  // and each inductor's contribution to the R equ.
1409  currentInductor = instanceData.begin();
1410  endInductor = instanceData.end();
1411  i = 0;
1412  while( currentInductor != endInductor )
1413  {
1414 
1415  qVec[((*currentInductor)->li_Branch)] += LOI[ i ];
1416 
1417  double current = inductorCurrents[ i ];
1418  double windings = (*currentInductor)->L;
1419 
1420  qVec[ li_RVar ] += rEqScaling * current * windings;
1421  i++;
1422  currentInductor++;
1423  }
1424 
1425  // load M terms if needed
1426  if( model_.includeMEquation )
1427  {
1428  double latestMag = mVarScaling * staVector[ li_MagVarState ];
1429 
1430  // M equation
1431  if(!getSolverState().dcopFlag)
1432  {
1433  qVec[ li_MagVar ] += mEqScaling * latestMag;
1434  }
1435  }
1436  return bsuccess;
1437 }
1438 
1439 
1440 //-----------------------------------------------------------------------------
1441 // Function : Instance::loadDAEFVector
1442 //
1443 // Purpose : Loads the F-vector contributions for a single
1444 // instance.
1445 //
1446 // Special Notes :
1447 //
1448 // Scope : public
1449 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1450 // Creation Date : 03/21/2005
1451 //-----------------------------------------------------------------------------
1453 {
1454  bool bsuccess=true;
1455  double mVarScaling = model_.mVarScaling;
1456  double rVarScaling = model_.rVarScaling;
1457  double mEqScaling = model_.mEqScaling;
1458  double rEqScaling = model_.rEqScaling;
1459 
1460  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1461  {
1462  Xyce::dout() << "Instance::loadDAEFVector------------------------" << std::endl
1463  << "\tname = " << getName() << std::endl;
1464  }
1465 
1466  Linear::Vector & staVector = *(extData.nextStaVectorPtr);
1467  Linear::Vector & solVector = *(extData.nextSolVectorPtr);
1468  Linear::Vector & staDerivVec = *(extData.nextStaDerivVectorPtr);
1469  Linear::Vector & stoVector = *(extData.nextStoVectorPtr);
1470 
1471  double * fVec = extData.daeFVectorRawPtr;
1472 
1473  double latestR = rVarScaling * stoVector[ li_RVarStore ];
1474 
1475  if(getSolverState().dcopFlag)
1476  {
1477  //enforce R = 0 in dc op
1478  latestR = 0.0;
1479  }
1480  // load M terms if needed
1481  if( model_.includeMEquation )
1482  {
1483  // for the M equation
1484  fVec[li_MagVar] -= mEqScaling * P * latestR / (model_.Path);
1485 
1486  // if |P| is near zero, then the M equation becomes dM/dt = 0, or M is
1487  // constant. In this case we'll add a diagonal element for M so that
1488  // sole dM/dt element in dQ/dX doesn't cause a time step too small error
1489  // Since P is normally very large, we'll test for |P| <= 1.0.
1490 
1491  if( fabs( P ) <= model_.pZeroTol )
1492  {
1493  fVec[li_MagVar] -= mVarScaling * staVector[ li_MagVarState ];
1494  }
1495  }
1496 
1497  // for the R equation
1498  fVec[li_RVar] -= rEqScaling * rVarScaling * stoVector[ li_RVarStore ];
1499 
1500  // used in scaling the branch equation;
1501  double mid=1.0;
1502  if( model_.factorMS )
1503  {
1504  mid = 1.0 + (1.0 - ((model_.Gap) / (model_.Path)))*P*(model_.Ms);
1505  }
1506  else
1507  {
1508  mid = 1.0 + (1.0 - ((model_.Gap) / (model_.Path)))*P;
1509  }
1510 
1511  // loop over each inductor and load it's F vector components
1512  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1513  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1514  int i=0;
1515  while( currentInductor != endInductor )
1516  {
1517  double current = solVector[(*currentInductor)->li_Branch];
1518  double vNodePos = solVector[(*currentInductor)->li_Pos];
1519  double vNodeNeg = solVector[(*currentInductor)->li_Neg];
1520 
1521 
1522  fVec[((*currentInductor)->li_Pos)] += scalingRHS * current;
1523 
1524  fVec[((*currentInductor)->li_Neg)] += -scalingRHS * current;
1525 
1526  fVec[((*currentInductor)->li_Branch)] += -((vNodePos - vNodeNeg)/mid);
1527  double windings = (*currentInductor)->L;
1528 
1529  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1530  {
1531  Xyce::dout() << " Inductor = " << (*currentInductor)->name
1532  << " li_Pos = " << (*currentInductor)->li_Pos
1533  << " li_Neg = " << (*currentInductor)->li_Neg
1534  << " li_Branch = " << (*currentInductor)->li_Branch
1535  << "\tPos/Neg current*windings = " << scalingRHS*current*windings
1536  << "\tBranch = " << ((vNodePos - vNodeNeg)/mid)
1537  << std::endl;
1538  }
1539  currentInductor++;
1540  i++;
1541  }
1542 
1543  return bsuccess;
1544 }
1545 
1546 //-----------------------------------------------------------------------------
1547 // Function : Instance::loadDAEdQdx
1548 //
1549 // Purpose : Loads the Q-vector contributions for a single instance.
1550 //
1551 // Scope : public
1552 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1553 // Creation Date : 03/21/2005
1554 //-----------------------------------------------------------------------------
1556 {
1557  bool bsuccess = true;
1558 
1559  double mVarScaling = model_.mVarScaling;
1560  double rVarScaling = model_.rVarScaling;
1561  double mEqScaling = model_.mEqScaling;
1562  double rEqScaling = model_.rEqScaling;
1563 
1564  int i;
1565 
1566  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1567  {
1568  Xyce::dout() << "Instance::loadDAEdQdx-----------------------" << std::endl
1569  << "\tname = " << getName() << std::endl;
1570  }
1571 
1572  Linear::Matrix * dQdxMatPtr = extData.dQdxMatrixPtr;
1573  // update the M equation if it's needed
1574  if( model_.includeMEquation )
1575  {
1576  // update M equation
1577  if(!getSolverState().dcopFlag)
1578  {
1579  (*dQdxMatPtr)[li_MagVar][mEquMOffset] += mEqScaling * mVarScaling;
1580  }
1581  }
1582 
1583  // update the R equation
1584  for( i = 0; i< numInductors; i++ )
1585  {
1586  (*dQdxMatPtr)[li_RVar][rEquInductorOffsets[i] ] += rEqScaling * inductanceVals[i];
1587  }
1588 
1589  // loop over each inductor and load it's Q vector components
1590  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1591  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1592  i = 0;
1593  while( currentInductor != endInductor )
1594  {
1595  for( int j=0; j<numInductors; ++j )
1596  {
1597  (*dQdxMatPtr)[((*currentInductor)->li_Branch)]
1598  [(*currentInductor)->inductorCurrentOffsets[j]] += LO[i][j];
1599  }
1600  i++;
1601  currentInductor++;
1602  }
1603 
1604  return bsuccess;
1605 }
1606 
1607 
1608 
1609 //-----------------------------------------------------------------------------
1610 // Function : Instance::loadDAEdFdx ()
1611 //
1612 // Purpose : Loads the F-vector contributions for a single instance.
1613 //
1614 // Special Notes :
1615 //
1616 // Scope : public
1617 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1618 // Creation Date : 03/21/2005
1619 //-----------------------------------------------------------------------------
1621 {
1622  bool bsuccess = true;
1623  double mVarScaling = model_.mVarScaling;
1624  double rVarScaling = model_.rVarScaling;
1625  double mEqScaling = model_.mEqScaling;
1626  double rEqScaling = model_.rEqScaling;
1627 
1628  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1629  {
1630  Xyce::dout() << "Instance::loadDAEdFdx----------------------" << std::endl
1631  << "\tname = " << getName() << std::endl;
1632  }
1633 
1634  Linear::Vector & staVector = *(extData.nextStaVectorPtr);
1635  Linear::Vector & solVector = *(extData.nextSolVectorPtr);
1636  Linear::Vector & staDerivVec = *(extData.nextStaDerivVectorPtr);
1637  Linear::Vector & stoVector = *(extData.nextStoVectorPtr);
1638  Linear::Matrix * dFdxMatPtr = extData.dFdxMatrixPtr;
1639 
1640  // udate dependent parameters
1641  //updateIntermediateVars();
1642 
1643  // pull these parameters up from the model class to make it easier
1644  // to view the equations.
1645  const double Gap = model_.Gap;
1646  const double Path = model_.Path;
1647 
1648  // terms that come up frequently
1649  double latestR = rVarScaling * stoVector[ li_RVarStore ];
1650 
1651  // inlcude M equation if it's part of the system
1652  if( model_.includeMEquation )
1653  {
1654  // terms for the M equation
1655  if(!getSolverState().dcopFlag)
1656  {
1657  (*dFdxMatPtr)[ li_MagVar ][ mEquMOffset ] -= mEqScaling * dP_dM * latestR / Path; // d/dM
1658  (*dFdxMatPtr)[ li_MagVar ][ mEquROffset ] -= mEqScaling * P * rVarScaling / Path; // d/dR
1659  (*dFdxMatPtr)[ li_MagVar ][ mEquVPosOffset ] -= mEqScaling * dP_dVp * latestR / Path; // d/dV_+
1660  (*dFdxMatPtr)[ li_MagVar ][ mEquVNegOffset ] -= mEqScaling * dP_dVn * latestR / Path; // d/dV_-
1661  for( int i = 0; i<numInductors; ++i)
1662  {
1663  (*dFdxMatPtr)[ li_MagVar ][mEquInductorOffsets[i] ] -=
1664  mEqScaling * dP_dI[i] * latestR / Path; // d/dI_i;
1665  }
1666  }
1667  else
1668  {
1669  // the above load for the M equation is basically zero in the dc op. We
1670  // need something on the diagonal for M to make the matrix non-singular
1671  (*dFdxMatPtr)[ li_MagVar ][ mEquMOffset ] += getSolverState().pdt_;
1672  }
1673 
1674  // if |P| is near zero, then the M equation becomes dM/dt = 0, or M is
1675  // constant. In this case we'll add a diagonal element for M so that
1676  // sole dM/dt element in dQ/dX doesn't cause a time step too small error
1677  // Since P is normally very large, we'll test for |P| <= 1.0.
1678  if( fabs( P ) <= model_.pZeroTol )
1679  {
1680  (*dFdxMatPtr)[ li_MagVar ][ mEquMOffset ] += 1.0;
1681  }
1682  }
1683 
1684  // update the R equation
1685 
1686  (*dFdxMatPtr)[ li_RVar ][rEquROffset] -= rEqScaling * rVarScaling;
1687 
1688  // loop over each inductor and load it's dFdx components
1689  double mid=0.0;
1690  if( model_.factorMS )
1691  {
1692  mid = 1.0 + (1.0 - ((model_.Gap) / (model_.Path)))*P*(model_.Ms);
1693  }
1694  else
1695  {
1696  mid = 1.0 + (1.0 - ((model_.Gap) / (model_.Path)))*P;
1697  }
1698 
1699  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1700  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1701  while( currentInductor != endInductor )
1702  {
1703  // do the normal work for an inductor
1704 
1705  (*dFdxMatPtr)[((*currentInductor)->li_Pos)] [((*currentInductor)->APosEquBraVarOffset)] += scalingRHS;
1706 
1707  (*dFdxMatPtr)[((*currentInductor)->li_Neg)] [((*currentInductor)->ANegEquBraVarOffset)] += -scalingRHS;
1708 
1709  (*dFdxMatPtr)[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquPosNodeOffset)] += -1.0/mid;
1710 
1711  (*dFdxMatPtr)[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquNegNodeOffset)] += 1.0/mid;
1712 
1713  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1714  {
1715  Xyce::dout()
1716  << "(*currentInductor)->li_Pos = " << (*currentInductor)->li_Pos << std::endl
1717  << "(*currentInductor)->li_Neg = " << (*currentInductor)->li_Neg << std::endl
1718  << "(*currentInductor)->li_Branch = " << (*currentInductor)->li_Branch << std::endl
1719  << "(*currentInductor)->APosEquBraVarOffset = " << (*currentInductor)->APosEquBraVarOffset << std::endl
1720  << "(*currentInductor)->ANegEquBraVarOffset = " << (*currentInductor)->ANegEquBraVarOffset << std::endl
1721  << "(*currentInductor)->ABraEquPosNodeOffset = " << (*currentInductor)->ABraEquPosNodeOffset << std::endl
1722  << "(*currentInductor)->ABraEquNegNodeOffset = " << (*currentInductor)->ABraEquNegNodeOffset << std::endl
1723  << "(*dFdxMatPtr)["<<((*currentInductor)->li_Pos)<<"] ["<<((*currentInductor)->APosEquBraVarOffset)<<"] = " << scalingRHS << std::endl
1724  << "(*dFdxMatPtr)["<<((*currentInductor)->li_Neg)<<"] ["<<((*currentInductor)->ANegEquBraVarOffset)<<"] = " << -scalingRHS << std::endl
1725  << "(*dFdxMatPtr)["<<((*currentInductor)->li_Branch)<<"]["<<((*currentInductor)->ABraEquPosNodeOffset)<<"] = " << -1/mid << std::endl
1726  << "(*dFdxMatPtr)["<<((*currentInductor)->li_Branch)<<"]["<<((*currentInductor)->ABraEquNegNodeOffset)<<"] = " << 1/mid << std::endl;
1727  }
1728 
1729  double delV = solVector[(*currentInductor)->li_Pos] - solVector[(*currentInductor)->li_Neg];
1730 
1731  for( int j = 0; j<numInductors; ++j )
1732  {
1733 
1734  (*dFdxMatPtr)[((*currentInductor)->li_Branch)][(*currentInductor)->inductorCurrentOffsets[j]] +=
1735  delV * (1.0 - (Gap/Path)) * dP_dI[j]/(mid*mid);
1736 
1737  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1738  {
1739  Xyce::dout() << "(*dFdxMatPtr)[((*currentInductor)->li_Branch)][(*currentInductor)->inductorCurrentOffsets[j]] = " << delV * (1 - (Gap/Path)) * dP_dI[j]/(mid*mid) << std::endl;
1740  }
1741  }
1742  if( model_.includeMEquation )
1743  {
1744  (*dFdxMatPtr)[(*currentInductor)->li_Branch][(*currentInductor)->magOffset] += delV * (1.0 - (Gap/Path)) * dP_dM/(mid*mid);
1745  }
1746  (*dFdxMatPtr)[(*currentInductor)->li_Branch][(*currentInductor)->vPosOffset] += delV * (1.0 - (Gap/Path)) * dP_dVp/(mid*mid);
1747 
1748  (*dFdxMatPtr)[(*currentInductor)->li_Branch][(*currentInductor)->vNegOffset] += delV * (1.0 - (Gap/Path)) * dP_dVn/(mid*mid);
1749 
1750  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
1751  {
1752 // Xyce::dout() << "(*dFdxMatPtr)[(*currentInductor)->li_Branch][(*currentInductor)->magOffset] = " << delV * (1 - (Gap/Path)) * dP_dM/(mid*mid) << std::endl
1753 // << "(*dFdxMatPtr)[(*currentInductor)->li_Branch][(*currentInductor)->vPosOffset] = " << delV * (1 - (Gap/Path)) * dP_dVp/(mid*mid) << std::endl
1754 // << "(*dFdxMatPtr)[(*currentInductor)->li_Branch][(*currentInductor)->vNegOffset] = " << delV * (1 - (Gap/Path)) * dP_dVn/(mid*mid) << std::endl;
1755  }
1756  currentInductor++;
1757  }
1758 
1759  return bsuccess;
1760 }
1761 
1762 //-----------------------------------------------------------------------------
1763 // Function : Instance::outputPlotFiles
1764 // Purpose : If requested by the use in the model statement,
1765 // this routine outputs values of the internal
1766 // state variables M, H and R to a file
1767 // named "Inductor_name.dat". File is opened
1768 // and closed in the contructor and destructor.
1769 // Special Notes :
1770 // Scope : public
1771 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1772 // Creation Date : 03/21/2005
1773 //-----------------------------------------------------------------------------
1774 bool Instance::outputPlotFiles(bool force_final_output)
1775 {
1776  bool bsuccess = true;
1777  if( outputStateVarsFlag && outputFileStreamPtr.get() && (*outputFileStreamPtr) )
1778  {
1779  Linear::Vector & solVector = *(extData.nextSolVectorPtr);
1780  Linear::Vector & staVector = *(extData.nextStaVectorPtr);
1781  Linear::Vector & stoVector = *(extData.nextStoVectorPtr);
1782  double mVarScaling = model_.mVarScaling;
1783  double rVarScaling = model_.rVarScaling;
1784 
1785  double latestMag = mVarScaling * staVector[ li_MagVarState ];
1786  if( model_.factorMS )
1787  {
1788  latestMag *= model_.Ms;
1789  }
1790  double latestR = rVarScaling * stoVector[ li_RVarStore ];
1791  (*outputFileStreamPtr)
1792  << getSolverState().currTime_ << " "
1793  << latestMag << "\t "
1794  << latestR << "\t "
1795  << staVector[ li_BVarStore ] << "\t "
1796  << staVector[ li_HVarStore ]
1797  << std::endl;
1798 
1799  }
1800  return bsuccess;
1801 }
1802 
1803 
1804 //-----------------------------------------------------------------------------
1805 // Function : Instance::setIC
1806 // Purpose :
1807 // Special Notes :
1808 // Scope : public
1809 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1810 // Creation Date : 03/21/2005
1811 //-----------------------------------------------------------------------------
1813 {
1814  int i_bra_sol;
1815  int i_f_state;
1816 
1817  bool bsuccess = true;
1818  return bsuccess;
1819 }
1820 
1821 //-----------------------------------------------------------------------------
1822 // Function : Instance::varTypes
1823 // Purpose :
1824 // Special Notes :
1825 // Scope : public
1826 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1827 // Creation Date : 03/21/2005
1828 //-----------------------------------------------------------------------------
1829 void Instance::varTypes( std::vector<char> & varTypeVec )
1830 {
1831  varTypeVec.resize(numInductors+2);
1832  for(int i=0; i<numInductors; i++)
1833  {
1834  varTypeVec[i] = 'I';
1835  }
1836  // I don't know what should be used for non I,V vars.
1837  varTypeVec[numInductors] = 'M';
1838  varTypeVec[numInductors+1] = 'R';
1839 }
1840 
1841 //-----------------------------------------------------------------------------
1842 // Function : Model::processParams
1843 // Purpose :
1844 // Special Notes :
1845 // Scope : public
1846 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1847 // Creation Date : 03/21/2005
1848 //-----------------------------------------------------------------------------
1850 {
1851  return true;
1852 }
1853 
1854 //----------------------------------------------------------------------------
1855 // Function : Model::processInstanceParams
1856 // Purpose :
1857 // Special Notes :
1858 // Scope : public
1859 // Creator : Dave Shirely, PSSI
1860 // Creation Date : 03/23/06
1861 //----------------------------------------------------------------------------
1863 {
1864  std::vector<Instance*>::iterator iter;
1865  std::vector<Instance*>::iterator first = instanceContainer.begin();
1866  std::vector<Instance*>::iterator last = instanceContainer.end();
1867 
1868  for (iter=first; iter!=last; ++iter)
1869  {
1870  (*iter)->processParams();
1871  }
1872 
1873  return true;
1874 }
1875 
1876 //-----------------------------------------------------------------------------
1877 // Function : Model::Model
1878 // Purpose : block constructor
1879 // Special Notes :
1880 // Scope : public
1881 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1882 // Creation Date : 03/21/2005
1883 //-----------------------------------------------------------------------------
1885  const Configuration & configuration,
1886  const ModelBlock & MB,
1887  const FactoryBlock & factory_block)
1888  : DeviceModel(MB, configuration.getModelParameters(), factory_block),
1889  A(0.0),
1890  Alpha(0.0),
1891  Area(0.0),
1892  BetaH(0.0),
1893  BetaM(0.0),
1894  C(0.0),
1895  DeltaVScaling(0.0),
1896  Gap(0.0),
1897  Kirr(0.0),
1898  Ms(0.0),
1899  Path(0.0),
1900  tempCoeff1(0.0),
1901  tempCoeff2(0.0),
1902  outputStateVars(0),
1903  factorMS(0),
1904  BCgsFactor( 10000.0 ),
1905  HCgsFactor( 0.012566370614359 ), // 4 pi / 1000
1906  UseConstantDeltaVScaling( false ),
1907  tnom(getDeviceOptions().tnom)
1908 {
1909  setLevel(1);
1910 
1911 
1912  // Set params to constant default values:
1913  setDefaultParams ();
1914 
1915  // Set params according to .model line and constant defaults from metadata:
1916  setModParams (MB.params);
1917 
1918  // scale gap, path and area from cm and cm^2 to m and m^2
1919  Gap *= 1.0e-2;
1920  Path *= 1.0e-2;
1921  Area *= 1.0e-4;
1922 
1923  if( BHSiUnits != 0 )
1924  {
1925  // user requested SI units over the default of CGS units. Change
1926  // conversion factor to unity.
1927  BCgsFactor=1.0;
1928  HCgsFactor=1.0;
1929  }
1930 
1931  // Set any non-constant parameter defaults:
1932  // when Ms factoring is off, scaling of M/R is still needed.
1933  if( factorMS == 0 )
1934  {
1935  mVarScaling=1.0e3;
1936  rVarScaling=1.0e3;
1937  mEqScaling=1.0e-3;
1938  rEqScaling=1.0e-3;
1939  }
1940 
1941 
1942  if (!given("TNOM"))
1943  {
1945  }
1946 
1947  // Calculate any parameters specified as expressions:
1949 
1950  // calculate dependent (ie computed) params and check for errors:
1951  processParams ();
1952 }
1953 
1954 //-----------------------------------------------------------------------------
1955 // Function : Model::~Model
1956 // Purpose : destructor
1957 // Special Notes :
1958 // Scope : public
1959 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1960 // Creation Date : 03/21/2005
1961 //-----------------------------------------------------------------------------
1963 {
1964  std::vector<Instance*>::iterator iter;
1965  std::vector<Instance*>::iterator first = instanceContainer.begin();
1966  std::vector<Instance*>::iterator last = instanceContainer.end();
1967 
1968  for (iter=first; iter!=last; ++iter)
1969  {
1970  delete (*iter);
1971  }
1972 }
1973 
1974 // additional Declarations
1975 
1976 //-----------------------------------------------------------------------------
1977 // Function : Model::printOutInstances
1978 // Purpose : debugging tool.
1979 // Special Notes :
1980 // Scope : public
1981 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1982 // Creation Date : 03/21/2005
1983 //-----------------------------------------------------------------------------
1984 std::ostream &Model::printOutInstances(std::ostream &os) const
1985 {
1986  std::vector<Instance*>::const_iterator iter;
1987  std::vector<Instance*>::const_iterator first = instanceContainer.begin();
1988  std::vector<Instance*>::const_iterator last = instanceContainer.end();
1989 
1990  int i, isize;
1991  isize = instanceContainer.size();
1992 
1993  os << std::endl;
1994  os << "Number of MutIndNonLin instances: " << isize << std::endl;
1995  os << " name=\t\tmodelName\tParameters" << std::endl;
1996  for (i=0, iter=first; iter!=last; ++iter, ++i)
1997  {
1998  os << " " << i << ": " << (*iter)->getName() << "\t";
1999  os << getName();
2000  os << std::endl;
2001  }
2002 
2003  os << std::endl;
2004 
2005  return os;
2006 }
2007 
2008 //-----------------------------------------------------------------------------
2009 // Function : Model::forEachInstance
2010 // Purpose :
2011 // Special Notes :
2012 // Scope : public
2013 // Creator : David Baur
2014 // Creation Date : 2/4/2014
2015 //-----------------------------------------------------------------------------
2016 /// Apply a device instance "op" to all instances associated with this
2017 /// model
2018 ///
2019 /// @param[in] op Operator to apply to all instances.
2020 ///
2021 ///
2022 void Model::forEachInstance(DeviceInstanceOp &op) const /* override */
2023 {
2024  for (std::vector<Instance *>::const_iterator it = instanceContainer.begin(); it != instanceContainer.end(); ++it)
2025  op(*it);
2026 }
2027 
2028 
2029 //-----------------------------------------------------------------------------
2030 // Function : Master::updateState
2031 // Purpose :
2032 // Special Notes :
2033 // Scope : public
2034 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
2035 // Creation Date : 11/25/08
2036 //-----------------------------------------------------------------------------
2037 bool Master::updateState (double * solVec, double * staVec, double * stoVec)
2038 {
2039  bool bsuccess = true;
2040  bool tmpBool = true;
2041 
2042  for (InstanceVector::const_iterator it = getInstanceBegin(); it != getInstanceEnd(); ++it)
2043  {
2044  tmpBool = (*it)->updatePrimaryState ();
2045  bsuccess = bsuccess && tmpBool;
2046  }
2047 
2048  return bsuccess;
2049 }
2050 
2051 //-----------------------------------------------------------------------------
2052 // Function : Master::updateSecondaryState
2053 // Purpose :
2054 // Special Notes :
2055 // Scope : public
2056 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
2057 // Creation Date : 11/25/08
2058 //-----------------------------------------------------------------------------
2059 bool Master::updateSecondaryState (double * staDerivVec, double * stoVec)
2060 {
2061  bool bsuccess = true;
2062  bool tmpBool = true;
2063 
2064  for (InstanceVector::const_iterator it = getInstanceBegin(); it != getInstanceEnd(); ++it)
2065  {
2066  tmpBool = (*it)->updateSecondaryState ();
2067  bsuccess = bsuccess && tmpBool;
2068  }
2069 
2070  return bsuccess;
2071 }
2072 
2073 //-----------------------------------------------------------------------------
2074 // Function : Master::loadDAEVectors
2075 // Purpose :
2076 // Special Notes :
2077 // Scope : public
2078 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
2079 // Creation Date : 11/25/08
2080 //-----------------------------------------------------------------------------
2081 bool Master::loadDAEVectors (double * solVec, double * fVec, double *qVec, double * bVec, double * storeLeadF, double * storeLeadQ, double * leadF, double * leadQ, double * junctionV)
2082 {
2083  bool bsuccess = true;
2084  bool tmpBool = true;
2085 
2086  for (InstanceVector::const_iterator it = getInstanceBegin(); it != getInstanceEnd(); ++it)
2087  {
2088  tmpBool = (*it)->loadDAEFVector();
2089  bsuccess = bsuccess && tmpBool;
2090  tmpBool = (*it)->loadDAEQVector();
2091  bsuccess = bsuccess && tmpBool;
2092  }
2093 
2094  return bsuccess;
2095 }
2096 
2097 //-----------------------------------------------------------------------------
2098 // Function : Master::loadDAEMatrices
2099 // Purpose :
2100 // Special Notes :
2101 // Scope : public
2102 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
2103 // Creation Date : 11/25/08
2104 //-----------------------------------------------------------------------------
2105 bool Master::loadDAEMatrices (Linear::Matrix & dFdx, Linear::Matrix & dQdx)
2106 {
2107  bool bsuccess = true;
2108  bool tmpBool = true;
2109 
2110  for (InstanceVector::const_iterator it = getInstanceBegin(); it != getInstanceEnd(); ++it)
2111  {
2112  tmpBool = (*it)->loadDAEdFdx ();
2113  bsuccess = bsuccess && tmpBool;
2114  tmpBool = (*it)->loadDAEdQdx ();
2115  bsuccess = bsuccess && tmpBool;
2116  }
2117 
2118  return bsuccess;
2119 }
2120 
2121 Device *Traits::factory(const Configuration &configuration, const FactoryBlock &factory_block)
2122 {
2123 
2124  return new Master(configuration, factory_block, factory_block.solverState_, factory_block.deviceOptions_);
2125 }
2126 
2128 {
2130  .registerDevice("min", 1)
2131  .registerModelType("min", 1)
2132  .registerModelType("core", 1);
2133 }
2134 
2135 } // namespace MutIndNonLin
2136 } // namespace Device
2137 } // namespace Xyce
const InstanceName & getName() const
std::vector< std::string > inductorNames
void varTypes(std::vector< char > &varTypeVec)
const DeviceOptions & deviceOptions_
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
void registerStoreLIDs(const std::vector< int > &staLIDVecRef)
void addStateNode(Util::SymbolTable &symbol_table, int index, const InstanceName &instance_name, const std::string &lead_name)
const std::string & getSubcircuitName() const
Decodes the device name.
double pdt_
Previous delta time alpha/dt (Many devices)
Linear::Vector * currSolVectorPtr
Linear::Vector * nextSolVectorPtr
virtual bool loadDAEMatrices(Linear::Matrix &dFdx, Linear::Matrix &dQdx)
Populates the device's Jacobian object with these pointers.
bool given(const std::string &parameter_name) const
Pure virtual class to augment a linear system.
Devices and models are each named.
void addInternalNode(Util::SymbolTable &symbol_table, int index, const InstanceName &instance_name, const std::string &lead_name)
void registerStateLIDs(const std::vector< int > &staLIDVecRef)
std::vector< Instance * > instanceContainer
void setNumStoreVars(int num_store_vars)
InstanceVector::const_iterator getInstanceEnd() const
Returns an iterator to the ending of the vector of all instances created for this device...
#define AssertLIDs(cmp)
std::vector< std::vector< double > > LO
virtual void registerJacLIDs(const JacobianStamp &jacLIDVec)
double tnom
nominal temperature for device params.
InstanceVector::const_iterator getInstanceBegin() const
Returns an iterator to the beginning of the vector of all instances created for this device...
Teuchos::RCP< std::ofstream > outputFileStreamPtr
std::vector< Param > params
Parameters from the line.
void setParams(const std::vector< Param > &params)
const std::string & getName() const
std::vector< std::string > inductorsNode2
static Device * factory(const Configuration &configuration, const FactoryBlock &factory_block)
The FactoryBlock contains parameters needed by the device, instance and model creation functions...
const DeviceOptions & getDeviceOptions() const
bool outputPlotFiles(bool force_final_output)
void loadNodeSymbols(Util::SymbolTable &symbol_table) const
Populates and returns the store name map.
Linear::Vector * nextStaVectorPtr
static Config< T > & addConfiguration()
Adds the device to the Xyce device configuration.
static void loadInstanceParameters(ParametricData< Instance > &instance_parameters)
Linear::Matrix * dFdxMatrixPtr
virtual void forEachInstance(DeviceInstanceOp &op) const
Apply a device instance "op" to all instances associated with this model.
The Device class is an interface for device implementations.
Definition: N_DEV_Device.h:101
bool processInstanceParams()
processInstanceParams
void addStoreNode(Util::SymbolTable &symbol_table, int index, const InstanceName &instance_name, const std::string &lead_name)
virtual bool updateState(double *solVec, double *staVec, double *stoVec)
Updates the devices state information.
std::vector< std::string > inductorsNode1
const SolverState & solverState_
Class Configuration contains device configuration data.
#define CONSTTANH_THRESH
Definition: N_DEV_Const.h:117
virtual bool loadDAEVectors(double *solVec, double *fVec, double *qVec, double *bVec, double *storeLeadF, double *storeLeadQ, double *leadF, double *leadQ, double *junctionV)
Populates the device's ExternData object with these pointers.
void registerJacLIDs(const std::vector< std::vector< int > > &jacLIDVec)
const SolverState & getSolverState() const
#define M_PI
std::vector< int > inductorCurrentOffsets
Linear::Vector * nextStoVectorPtr
bool updateTemperature(const double &temp_tmp)
Instance(const Configuration &configuration, const InstanceBlock &IB, Model &Iiter, const FactoryBlock &factory_block)
virtual std::ostream & printOutInstances(std::ostream &os) const
std::string spiceInternalName(const InstanceName &instance_name, const std::string &lead)
virtual bool updateSecondaryState(double *staDeriv, double *stoVec)
Updates the devices secondary state information.
ModelBlock represents a .MODEL line from the netlist.
std::vector< std::string > couplingInductor
std::vector< std::vector< int > > jacStamp
Manages parameter binding for class C.
Definition: N_DEV_Pars.h:214
InstanceBlock represent a device instance line from the netlist.
double currTime_
DeviceEntity for expression time, breakpoints DeviceMgr for dependent parameters, breakpoints...
std::vector< Param > params
Linear::Matrix * dQdxMatrixPtr
const std::vector< std::vector< int > > & jacobianStamp() const
void registerLIDs(const std::vector< int > &intLIDVecRef, const std::vector< int > &extLIDVecRef)
std::vector< InductorInstanceData * > instanceData
Linear::Vector * nextStaDerivVectorPtr
void setModParams(const std::vector< Param > &params)
static void loadModelParameters(ParametricData< Model > &model_parameters)