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