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