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