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.53 $
40 //
41 // Revision Date : $Date: 2014/05/19 15:49:14 $
42 //
43 // Current Owner : $Author: dgbaur $
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().getEncodedName() );
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  UserError(*this) << "Could not open file " << filename << " for output of state variables";
365  }
366  else {
367  (*outputFileStreamPtr).setf(std::ios::scientific, std::ios::floatfield );
368  (*outputFileStreamPtr).width(20);
369  (*outputFileStreamPtr).precision(12);
370  }
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.get() && outputFileStreamPtr->is_open() )
510  {
511  outputFileStreamPtr->close();
512  }
513 
514  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
515  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
516  for ( ; currentInductor != endInductor ; ++currentInductor)
517  {
518  if (*currentInductor != NULL)
519  {
520  delete *currentInductor;
521  *currentInductor = NULL;
522  }
523  }
524  instanceData.clear();
525 }
526 
527 //-----------------------------------------------------------------------------
528 // Function : Instance::registerLIDs
529 // Purpose :
530 // Special Notes :
531 // Scope : public
532 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
533 // Creation Date : 03/21/2005
534 //-----------------------------------------------------------------------------
535 void Instance::registerLIDs(const std::vector<int> & intLIDVecRef,
536  const std::vector<int> & extLIDVecRef)
537 {
538  AssertLIDs(intLIDVecRef.size() == numIntVars);
539  AssertLIDs(extLIDVecRef.size() == numExtVars);
540 
541  // copy over the global ID lists.
542  intLIDVec = intLIDVecRef;
543  extLIDVec = extLIDVecRef;
544 
545  // Now use these lists to obtain the indices into the
546  // linear algebra entities. This assumes an order.
547  // For the matrix indices, first do the rows.
548  // get the current values of the inductances and currentOffsets
549  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
550  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
551  int i = 0;
552  int j = 0;
553  while( currentInductor != endInductor )
554  {
555  (*currentInductor)->li_Pos = extLIDVec[ i++ ];
556  (*currentInductor)->li_Neg = extLIDVec[ i++ ];
557  (*currentInductor)->li_Branch = intLIDVec[ j++ ];
558  currentInductor++;
559  }
560 
561  if(includeDeltaM)
562  {
563  // now get deltaHapp and deltaM
564  //li_deltaHappVar = intLIDVec[ j++ ];
565  li_deltaMagVar = intLIDVec[ j++ ];
566  }
567 
568 #ifdef Xyce_DEBUG_DEVICE
569  if (getDeviceOptions().debugLevel > 0)
570  {
571  Xyce::dout() << "Instance::registerLIDs------------------------" << std::endl;
572  currentInductor = instanceData.begin();
573  i=0;
574  while( currentInductor != endInductor )
575  {
576  Xyce::dout() << "Inductor [ " << i++ << " ] "
577  << " li_Pos = " << (*currentInductor)->li_Pos
578  << " li_Neg = " << (*currentInductor)->li_Neg
579  << " li_Branch = " << (*currentInductor)->li_Branch << std::endl;
580  currentInductor++;
581  }
582  if(includeDeltaM)
583  {
584  Xyce::dout() << " li_deltaMagVar = " << li_deltaMagVar << std::endl;
585  }
586  }
587 #endif
588 }
589 
590 //-----------------------------------------------------------------------------
591 // Function : Instance::getIntNameMap
592 // Purpose :
593 // Special Notes :
594 // Scope : public
595 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
596 // Creation Date : 03/21/2005
597 //-----------------------------------------------------------------------------
598 std::map<int,std::string> & Instance::getIntNameMap()
599 {
600  // set up the internal name map, if it hasn't been already.
601  if (intNameMap.empty ())
602  {
603  for (std::vector< InductorInstanceData* >::const_iterator it = instanceData.begin(), end = instanceData.end(); it != end; ++it)
604  {
605  intNameMap[ (*it)->li_Branch ] = spiceInternalName(getName(), (*it)->name + "_branch");
606  }
607  }
608 
609  return intNameMap;
610 }
611 
612 
613 //-----------------------------------------------------------------------------
614 // Function : Instance::registerStateLIDs
615 // Purpose :
616 // Special Notes :
617 // Scope : public
618 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
619 // Creation Date : 03/21/2005
620 //-----------------------------------------------------------------------------
621 void Instance::registerStateLIDs( const std::vector<int> & staLIDVecRef )
622 {
623  AssertLIDs(staLIDVecRef.size() == numStateVars);
624 
625  // copy over the global ID lists.
626  staLIDVec = staLIDVecRef;
627 
628 #ifdef Xyce_DEBUG_DEVICE
629  if (getDeviceOptions().debugLevel > 0)
630  {
631  Xyce::dout() << "Instance::registerStateLIDs-------------------" << std::endl;
632  }
633 #endif
634 }
635 
636 //-----------------------------------------------------------------------------
637 // Function : Instance::jacobianStamp
638 // Purpose :
639 // Special Notes :
640 // Scope : public
641 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
642 // Creation Date : 03/21/2005
643 //-----------------------------------------------------------------------------
644 const std::vector< std::vector<int> > & Instance::jacobianStamp() const
645 {
646  return jacStamp;
647 }
648 
649 //-----------------------------------------------------------------------------
650 // Function : Instance::registerJacLIDs
651 // Purpose :
652 // Special Notes :
653 // Scope : public
654 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
655 // Creation Date : 03/21/2005
656 //-----------------------------------------------------------------------------
657 void Instance::registerJacLIDs( const std::vector< std::vector<int> > & jacLIDVec )
658 {
659  DeviceInstance::registerJacLIDs( jacLIDVec );
660 #ifdef Xyce_DEBUG_DEVICE
661  if (getDeviceOptions().debugLevel > 0)
662  {
663  Xyce::dout() << "Instance::registerJacLIDs ----------------------------" << std::endl;
664 
665  Xyce::dout() << "jacLIDVec = " << std::endl;
666  for( int i = 0; i<jacStamp.size(); ++i )
667  {
668  Xyce::dout() << "jacLIDVec[ " << i << " ] = { ";
669  for( int j=0; j<jacLIDVec[i].size(); ++j)
670  {
671  Xyce::dout() << jacLIDVec[i][j];
672  if( j != ( jacLIDVec[i].size() -1 ) )
673  {
674  Xyce::dout() << ", ";
675  }
676  }
677  Xyce::dout() << " }" << std::endl;
678  }
679  }
680 #endif
681 
682  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
683  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
684  int i = 0;
685  while( currentInductor != endInductor )
686  {
687  (*currentInductor)->APosEquBraVarOffset = jacLIDVec[ 2*i ][ 0 ];
688  (*currentInductor)->ANegEquBraVarOffset = jacLIDVec[ 2*i + 1 ][ 0 ];
689  (*currentInductor)->vPosOffset = jacLIDVec[ 2*numInductors + i ][ 0 ];
690  (*currentInductor)->vNegOffset = jacLIDVec[ 2*numInductors + i ][ 1 ];
691 
692  (*currentInductor)->ABraEquPosNodeOffset = jacLIDVec[ 2*numInductors + i ][ 0 ];
693  (*currentInductor)->ABraEquNegNodeOffset = jacLIDVec[ 2*numInductors + i ][ 1 ];
694  for( int j=0; j<numInductors; ++j )
695  {
696  if( i == j )
697  {
698  (*currentInductor)->ABraEquBraVarOffset = jacLIDVec[ 2*numInductors + i ][ j + 2 ];
699  }
700  (*currentInductor)->inductorCurrentOffsets[ j ] = jacLIDVec[ 2*numInductors + i ][ j + 2 ];
701  }
702 
703  currentInductor++;
704  i++;
705  }
706 
707  if(includeDeltaM)
708  {
709  // now get the deltaM
710  for( int i=0; i<numInductors; i++)
711  {
712  deltaMEquInductorOffsets[i] = jacLIDVec[ 3*numInductors ][ i ];
713  }
714  deltaMEquDeltaMOffset = jacLIDVec[ 3*numInductors ][ numInductors ];
715  }
716 
717 #ifdef Xyce_DEBUG_DEVICE
718  if (getDeviceOptions().debugLevel > 0)
719  {
720  currentInductor = instanceData.begin();
721  i=0;
722  while( currentInductor != endInductor )
723  {
724  Xyce::dout() << "Inductor [ " << i << " ] " << (*currentInductor)->name << std::endl
725  << " APosEquBraVarOffset = " << (*currentInductor)->APosEquBraVarOffset << std::endl
726  << " ANegEquBraVarOffset = " << (*currentInductor)->ANegEquBraVarOffset << std::endl
727  << " vPosOffset = " << (*currentInductor)->vPosOffset << std::endl
728  << " vNegOffset = " << (*currentInductor)->vNegOffset << std::endl
729  << " ABraEquPosNodeOffset = " << (*currentInductor)->ABraEquPosNodeOffset << std::endl
730  << " ABraEquNegNodeOffset = " << (*currentInductor)->ABraEquNegNodeOffset << std::endl
731  << " ABraEquBraVarOffset = " << (*currentInductor)->ABraEquBraVarOffset << std::endl
732  << " magOffset = " << (*currentInductor)->magOffset << std::endl;
733  Xyce::dout() << "\tInductor branch offsets = { ";
734  for( int j=0; j<numInductors ; ++j )
735  {
736  Xyce::dout() << (*currentInductor)->inductorCurrentOffsets[ j ] << ", ";
737  }
738  Xyce::dout() << "} " << std::endl;
739  i++;
740  currentInductor++;
741  }
742  Xyce::dout() << std::endl;
743 
744  if(includeDeltaM)
745  {
746  Xyce::dout() << "deltaMEquInductorOffsets = ";
747  for( int i=0; i<numInductors; i++ )
748  {
749  Xyce::dout() << deltaMEquInductorOffsets[i] << " ";
750  }
751  Xyce::dout() //<< "deltaMEquDeltaHappOffset = " << deltaMEquDeltaHappOffset
752  << " deltaMEquDeltaMOffset = " << deltaMEquDeltaMOffset << std::endl;
753  }
754  }
755 #endif
756 }
757 
758 
759 //-----------------------------------------------------------------------------
760 // Function : Instance::processParams
761 // Purpose :
762 // Special Notes :
763 // Scope : public
764 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
765 // Creation Date : 03/21/2005
766 //-----------------------------------------------------------------------------
768 {
769  // now set the temperature related stuff.
770  if (tempGiven)
771  {
773  }
774 
775  return true;
776 }
777 
778 //-----------------------------------------------------------------------------
779 // Function : Instance::updateTemperature
780 // Purpose :
781 // Special Notes :
782 // Scope : public
783 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
784 // Creation Date : 03/21/2005
785 //-----------------------------------------------------------------------------
786 bool Instance::updateTemperature ( const double & temp)
787 {
788  bool bsuccess = true;
789 
790  // current temp difference from reference temp.
791  double difference = temp - model_.tnom;
792 
793  std::vector< InductorInstanceData* >::iterator currentData = instanceData.begin();
794  while( currentData != instanceData.end() )
795  {
796  double factor = 1.0 + (model_.tempCoeff1)*difference +
797  (model_.tempCoeff2)*difference*difference;
798  (*currentData)->L = ((*currentData)->baseL)*factor;
799  currentData++;
800  }
801 
802  // now that the inductances have changed we need to update the matrix.
804 
805  return bsuccess;
806 }
807 
808 //-----------------------------------------------------------------------------
809 // Function : Instance::updateIntermediateVars
810 // Purpose : updates a set of common variables used by RHS and jacobian
811 // Special Notes :
812 // Scope : public
813 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
814 // Creation Date : 03/21/2005
815 //-----------------------------------------------------------------------------
817 {
818  bool bsuccess = true;
819 #ifdef Xyce_DEBUG_DEVICE
820  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
821  {
822  Xyce::dout() << "Instance::updateIntermediateVars" << std::endl;
823  }
824 #endif
825  N_LAS_Vector & solVector = *(extData.nextSolVectorPtr);
826  N_LAS_Vector & lastSolVector = *(extData.lastSolVectorPtr);
827 
828  // some parameters in the model class that we will use often
829  const double A = model_.A;
830  const double Alpha = model_.Alpha;
831  const double Area = model_.Area;
832  const double BetaH = model_.BetaH;
833  const double BetaM = model_.BetaM;
834  const double C = model_.C;
835  const double DeltaV = model_.DeltaV;
836  const double Gap = model_.Gap;
837  const double Ms = model_.Ms;
838  const double Kirr = model_.Kirr;
839  const double Path = model_.Path;
840  const double Vinf = model_.Vinf;
841 
842  double latestMag; //NoMag = solVector[ li_MagVar ];
843 
844  //sum of currents through the inductors
845  branchCurrentSum = 0.0;
846  for(int i=0; i<numInductors; i++)
847  {
848  branchCurrentSum += solVector[instanceData[i]->li_Branch] * inductanceVals[ i ];
849  }
850 
851  latestMag = MagVar + MagVarUpdate;
852 
853  // used in voltage drop over first inductor
854  double V1Pos = solVector[(instanceData[0])->li_Pos];
855  double V1Neg = solVector[(instanceData[0])->li_Neg];
856 
857  double qV = (DeltaV / Vinf) * (V1Pos - V1Neg);
858 
859  double tanh_qV = 0.0;
860  if (fabs(qV) < CONSTTANH_THRESH)
861  {
862  tanh_qV = tanh(qV);
863  }
864  else if (qV < 0.0)
865  {
866  tanh_qV = -1.0;
867  }
868  else
869  {
870  tanh_qV = 1.0;
871  }
872 
873  double Happ = branchCurrentSum / Path;
874 
875 #ifdef MS_FACTORING2
876  double H = Happ - (Gap / Path) * latestMag * Ms;
877  double He = H + Alpha * latestMag * Ms;
878 #else
879  double H = Happ - (Gap / Path) * latestMag;
880  double He = H + Alpha * latestMag;
881 #endif
882 
883  double Heo = BetaH*A;
884 
885  // terms that come up frequently
886  double gap_path = Gap / Path;
887  double He2 = He*He;
888  double Heo2 = Heo*Heo;
889  double sq_Heo2He2 = sqrt(Heo2 + He2);
890 
891  double delM0 = model_.BetaM * Ms;
892  double Man = Ms * He / ( A + sq_Heo2He2 );
893 #ifdef MS_FACTORING2
894  double delM = Man - latestMag*Ms;
895 #else
896  double delM = Man - latestMag;
897 #endif
898 
899  double delM2 = delM*delM;
900  double delM02 = delM0*delM0;
901  double sq_delM02delM2 = sqrt( delM02 + delM2 );
902 
903  double Pold = P;
904 
905  #ifdef MS_FACTORING2
906  double Mirrp = (delM * tanh_qV + sq_delM02delM2 ) / (2*( Kirr- Alpha * sq_delM02delM2));
907  double Manp = Ms * (A + Heo2/sq_Heo2He2) / pow(A + sq_Heo2He2, 2.0);
908  P = ( C * Manp + (1 - C) * Mirrp) / ((1 + (gap_path - Alpha) * C * Manp + gap_path * (1-C) * Mirrp)*Ms);
909  #else
910  double Mirrp = (delM * tanh_qV + sq_delM02delM2 ) / (2*( Kirr- Alpha * sq_delM02delM2));
911  double Manp = Ms*(A + Heo2/sq_Heo2He2) / pow(A + sq_Heo2He2, 2.0);
912  P = ( C * Manp + (1 - C) * Mirrp) / (1 + (gap_path - Alpha) * C * Manp + gap_path * (1-C) * Mirrp);
913 #endif
914 
915 
916  // at this point we have P so now we can update mag.
917  //
918  // The problem is that if deltaM is too big, then we need to shrink
919  // the time step. One way to control this is to set a max time
920  // step. But what we really need to do is calculate deltaM and
921  // then if it's over some fraction of Ms then turn on the limiting
922  // flag (or bail on the step but I think turning on limiting is
923  // safer and if we hit maxItter with it on then we'll get that step
924  // rejected.
925 
926  if( useRKIntegration )
927  {
928  // use 4th order runga-kutta to estimate MagVarUpdate
930  MagVarUpdate = stepLen * ( PFunctionHistory[0] +
931  2*PFunctionHistory[1] +
932  2*PFunctionHistory[2] +
933  P) / 6;
934  }
935  else
936  {
937  // forward euler method
938  MagVarUpdate = P * (branchCurrentSum - oldBranchCurrentSum) / model_.Path;
939 
940  // trap
941  double MagVarUpdateWithTrap = 0.5 * (P + PPreviousStep) * (branchCurrentSum - oldBranchCurrentSum) / model_.Path;
942  origFlag = true;
943  if( fabs( MagVarUpdate ) > 0.25 * Ms )
944  {
945  // step was too big, so
946  // turn on limiting
947  origFlag = false;
948  }
949  }
950 
951  latestMag = MagVar + MagVarUpdate;
952 
953  if(includeDeltaM)
954  {
955  // in this case we're scaling MagVarUpdate by Ms because it's being solved
956  // with the full system and big changes in
957  MagVarUpdate /=model_.Ms;
958  }
959 
960  double dP_Denom = 1.0 + (gap_path - Alpha)*C*Manp + gap_path * (1.0-C) * Mirrp;
961 
962  // now get dP_dI for each inductor
963  for( int i=0; i<numInductors; i++)
964  {
965 
966  dHe_dI[ i ] = inductanceVals[ i ] / Path;
967 
968  dManp_dI[i] = ( -Ms * He / (pow(A + sq_Heo2He2, 2.0)*sq_Heo2He2)) *
969  ( (Heo2 / (Heo2 + He2)) + (2.0*(A + Heo2 / sq_Heo2He2)/(A+sq_Heo2He2)) ) * dHe_dI[i];
970 
971  ddelM_dI[i] = (Ms / (A + sq_Heo2He2)) * (1.0 - He2/((A + sq_Heo2He2)*sq_Heo2He2)) * dHe_dI[i];
972 
973  dMirrp_dI[i] = (1.0/(2.0*(Kirr - Alpha*sq_delM02delM2))) *
974  (tanh_qV + delM/sq_delM02delM2 +
975  (2.0*Alpha*delM*(delM*tanh_qV +
976  sq_delM02delM2)/(2.0*(Kirr-Alpha*sq_delM02delM2)*sq_delM02delM2))) * ddelM_dI[i];
977 
978  dP_dI[i] = (1.0/dP_Denom) * (C * dManp_dI[i] + (1.0-C) * dMirrp_dI[i]) -
979  ( (C*Manp + (1.0-C)*Mirrp)/pow(dP_Denom,2.0) ) *
980  ( (gap_path - Alpha)*C*dManp_dI[i] + gap_path*(1.0-C)*dMirrp_dI[i] );
981 
982  }
983 
984  return true;
985 }
986 
987 //-----------------------------------------------------------------------------
988 // Function : Instance::updateInductanceMatrix()
989 // Purpose : A matrix of inductances is used often enough that it
990 // calculated and stored as a member variable here
991 // If and inductance ever changes say from updating
992 // the temperature or a parameter udpate, then this
993 // routine must be called again.
994 // Special Notes :
995 // Scope : public
996 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
997 // Creation Date : 03/21/2005
998 //-----------------------------------------------------------------------------
1000 {
1001  std::vector< InductorInstanceData* >::iterator
1002  currentInductor = instanceData.begin();
1003  std::vector< InductorInstanceData* >::iterator
1004  endInductor = instanceData.end();
1005 
1006  // collec the inductances
1007  int i=0;
1008  while( currentInductor != endInductor )
1009  {
1010  inductanceVals[ i ] = (*currentInductor)->L;
1011  i++;
1012  currentInductor++;
1013  }
1014 
1015  double Area = model_.Area;
1016  double Path = model_.Path;
1017 
1018  // compute the inductance matrix
1019  for( i=0; i<numInductors; ++i)
1020  {
1021  for( int j=0; j<numInductors; ++j)
1022  {
1023  // 4.0e-7 * M_PI is a magnetic constant, the permeability of free space [Henries/m]
1024  LO[i][j] = mutualCup * 4.0e-7 * M_PI * (Area / Path) * inductanceVals[i] * inductanceVals[j];
1025  }
1026  }
1027 
1028 }
1029 
1030 //-----------------------------------------------------------------------------
1031 // Function : Instance::updatePrimaryState
1032 // Purpose :
1033 // Special Notes :
1034 // Scope : public
1035 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1036 // Creation Date : 03/21/2005
1037 //-----------------------------------------------------------------------------
1039 {
1040  bool bsuccess = true;
1041 #ifdef Xyce_DEBUG_DEVICE
1042  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
1043  {
1044  Xyce::dout() << "Instance::updatePrimaryState---------------" << std::endl
1045  << "\tname = " << getName() << std::endl;
1046  }
1047 #endif
1048  // udate dependent parameters
1050 
1051 #if 0
1052  // don't need to do this as we're not using the state vector
1053  N_LAS_Vector & solVector = *(extData.nextSolVectorPtr);
1054  N_LAS_Vector & staVector = *(extData.nextStaVectorPtr);
1055 
1056  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1057  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1058  int i = 0;
1059  while( currentInductor != endInductor )
1060  {
1061  double current = solVector[ ( (*currentInductor)->li_Branch) ];
1062  if( (getSolverState().dcopFlag) && ((*currentInductor)->ICGiven) )
1063  {
1064  current = (*currentInductor)->IC;
1065  }
1066  // place this value for the charge in the state vector.
1067  staVector[((*currentInductor)->li_currentState)] = current;
1068  currentInductor++;
1069  i++;
1070  }
1071 #endif
1072 
1073  return bsuccess;
1074 }
1075 
1076 //-----------------------------------------------------------------------------
1077 // Function : Instance::updateSecondaryState
1078 // Purpose :
1079 // Special Notes :
1080 // Scope : public
1081 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1082 // Creation Date : 03/21/2005
1083 //-----------------------------------------------------------------------------
1085 {
1086  bool bsuccess = true;
1087  return bsuccess;
1088 }
1089 
1090 //-----------------------------------------------------------------------------
1091 // Function : Instance::acceptStep
1092 // Purpose : This function updates the value of MagVar
1093 //
1094 // Scope : public
1095 // Creator : Rich Schiek, SNL, Electrical Systems Modeling
1096 // Creation Date : 01/25/2011
1097 //-----------------------------------------------------------------------------
1099 {
1100  if (!getSolverState().dcopFlag)
1101  {
1102  if(includeDeltaM)
1103  {
1105  }
1106  else
1107  {
1108  MagVar += MagVarUpdate;
1109  }
1111  PPreviousStep = P;
1112  if( fabs(MagVar) > 2*model_.Ms )
1113  {
1114  MagVar = 0.0;
1115  }
1116 
1117  if( useRKIntegration )
1118  {
1119  // fill in history for RK integration of dM/dH
1120  for(int i=0; i<2; i++)
1121  {
1124  }
1127  }
1129  }
1130 }
1131 
1132 //-----------------------------------------------------------------------------
1133 // Function : Instance::loadDeviceMask
1134 //
1135 // Purpose : Loads the zero elements of the device mask
1136 //
1137 // Special Notes : elements of the error vector associated with zero
1138 // elements of the mask will not be included in weighted
1139 // norms by the time integrator.
1140 //
1141 // Scope : public
1142 // Creator : Tom Russo, SNL, Electrical and Microsystems Modeling
1143 // Creation Date : 01/19/07
1144 //-----------------------------------------------------------------------------
1146 {
1147  bool returnVal=false;
1148 #if 0
1149  // ifndef Xyce_NO_MUTIND_MASK
1150  N_LAS_Vector * maskVectorPtr = extData.deviceMaskVectorPtr;
1151 
1152  (*maskVectorPtr)[li_MagVar] = 0.0;
1153  returnVal = true;
1154 #endif
1155  return (returnVal);
1156 }
1157 
1158 //-----------------------------------------------------------------------------
1159 // Function : Instance::loadDAEQVector
1160 //
1161 // Purpose : Loads the Q-vector contributions for a single
1162 // instance.
1163 //
1164 // Special Notes : The "Q" vector is part of a standard DAE formalism in
1165 // which the system of equations is represented as:
1166 //
1167 // f(x) = dQ(x)/dt + F(x) - B(t) = 0
1168 //
1169 // Scope : public
1170 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1171 // Creation Date : 03/21/2005
1172 //-----------------------------------------------------------------------------
1174 {
1175  bool bsuccess = true;
1176 
1177 #ifdef Xyce_DEBUG_DEVICE
1178  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
1179  {
1180  Xyce::dout() << "Instance::loadDAEQVector------------------------" << std::endl
1181  << "\tname = " << getName() << std::endl;
1182  }
1183 #endif
1184 
1185  N_LAS_Vector * daeQVecPtr = extData.daeQVectorPtr;
1186  N_LAS_Vector & solVector = *(extData.nextSolVectorPtr);
1187 
1188  // update LOI -- the following product
1189  // I = column vector of currents
1190  // L = row vector of inductances
1191  // LO = matrix = mutualCup * sqrt( L' * L )
1192  // LOI = column vector = mutualCup * sqrt( L' * L ) * I
1193  // LOI[1] = mutualCup * sqrt(L[1]*L[1])*I[1]) +
1194  // mutualCup * sqrt(L[1]*L[2])*I[2]) + ...
1195  // mutualCup * sqrt(L[1]*L[n])*I[n])
1196 
1197  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1198  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1199  int i = 0;
1200  while( currentInductor != endInductor )
1201  {
1202  if( (getSolverState().dcopFlag) && (*currentInductor)->ICGiven == true )
1203  {
1204  inductorCurrents[ i ] = (*currentInductor)->IC;
1205  }
1206  else
1207  {
1208  inductorCurrents[ i ] = solVector[ (*currentInductor)->li_Branch ];
1209  }
1210  i++;
1211  currentInductor++;
1212  }
1213 
1214  for( i = 0; i < numInductors; ++i )
1215  {
1216  LOI[ i ] = 0;
1217  for( int j = 0; j < numInductors; ++j )
1218  {
1219  LOI[i] += LO[i][j] * inductorCurrents[j];
1220  }
1221  }
1222 
1223  // loop over each inductor and load it's Q vector components
1224  // and each inductor's contribution to the R equ.
1225  currentInductor = instanceData.begin();
1226  endInductor = instanceData.end();
1227  i = 0;
1228  while( currentInductor != endInductor )
1229  {
1230 
1231  (*daeQVecPtr)[((*currentInductor)->li_Branch)] += LOI[ i ];
1232  i++;
1233  currentInductor++;
1234  }
1235 
1236  return bsuccess;
1237 }
1238 
1239 
1240 //-----------------------------------------------------------------------------
1241 // Function : Instance::loadDAEFVector
1242 //
1243 // Purpose : Loads the F-vector contributions for a single
1244 // instance.
1245 //
1246 // Special Notes : See the special notes for loadDAEFVector.
1247 //
1248 // Same as loadRHS, but without the capacitor
1249 // currents.
1250 //
1251 // Scope : public
1252 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1253 // Creation Date : 03/21/2005
1254 //-----------------------------------------------------------------------------
1256 {
1257  bool bsuccess=true;
1258 
1259 #ifdef Xyce_DEBUG_DEVICE
1260  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
1261  {
1262  Xyce::dout() << "Instance::loadDAEFVector------------------------" << std::endl
1263  << "\tname = " << getName() << std::endl;
1264  }
1265 #endif
1266 
1267  N_LAS_Vector * daeFVecPtr = extData.daeFVectorPtr;
1268  N_LAS_Vector & solVector = *(extData.nextSolVectorPtr);
1269 
1270  double Gap = model_.Gap;
1271  double Path = model_.Path;
1272 
1273 
1274  // used in scaling the branch equation;
1275  double mid = 1.0 + (1.0 - (Gap / Path))*P;
1276 
1277  // loop over each inductor and load it's F vector components
1278  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1279  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1280  int i=0;
1281  while( currentInductor != endInductor )
1282  {
1283  double current = solVector[(*currentInductor)->li_Branch];
1284  double vNodePos = solVector[(*currentInductor)->li_Pos];
1285  double vNodeNeg = solVector[(*currentInductor)->li_Neg];
1286 
1287 
1288  (*daeFVecPtr)[((*currentInductor)->li_Pos)] += current;
1289 
1290  (*daeFVecPtr)[((*currentInductor)->li_Neg)] += -current;
1291 
1292  (*daeFVecPtr)[((*currentInductor)->li_Branch)] += -((vNodePos - vNodeNeg)/mid);
1293 
1294  currentInductor++;
1295  i++;
1296  }
1297 
1298  if(includeDeltaM)
1299  {
1300  // the deltaHapp equation
1301  //(*daeFVecPtr)[li_deltaHappVar] += solVector[li_deltaHappVar];
1302  //(*daeFVecPtr)[li_deltaHappVar] -= HappVarUpdate;
1303 
1304  // the deltaM equation
1305  (*daeFVecPtr)[li_deltaMagVar] += solVector[li_deltaMagVar];
1306  (*daeFVecPtr)[li_deltaMagVar] -= MagVarUpdate;
1307  }
1308 
1309  return bsuccess;
1310 }
1311 
1312 //-----------------------------------------------------------------------------
1313 // Function : Instance::loadDAEdQdx
1314 //
1315 // Purpose : Loads the Q-vector contributions for a single
1316 // instance.
1317 // Scope : public
1318 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1319 // Creation Date : 03/21/2005
1320 //-----------------------------------------------------------------------------
1322 {
1323  bool bsuccess = true;
1324  int i;
1325 
1326 #ifdef Xyce_DEBUG_DEVICE
1327  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
1328  {
1329  Xyce::dout() << "Instance::loadDAEdQdx-----------------------" << std::endl
1330  << "\tname = " << getName() << std::endl;
1331  }
1332 #endif
1333  N_LAS_Matrix * dQdxMatPtr = extData.dQdxMatrixPtr;
1334 
1335  // loop over each inductor and load it's Q vector components
1336  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1337  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1338  i = 0;
1339  while( currentInductor != endInductor )
1340  {
1341  for( int j=0; j<numInductors; ++j )
1342  {
1343  (*dQdxMatPtr)[((*currentInductor)->li_Branch)]
1344  [(*currentInductor)->inductorCurrentOffsets[j]] += LO[i][j];
1345  }
1346  i++;
1347  currentInductor++;
1348  }
1349 
1350  return bsuccess;
1351 }
1352 
1353 
1354 
1355 //-----------------------------------------------------------------------------
1356 // Function : Instance::loadDAEdFdx ()
1357 //
1358 // Purpose : Loads the F-vector contributions for a single
1359 // instance.
1360 //
1361 // Special Notes :
1362 //
1363 // Scope : public
1364 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1365 // Creation Date : 03/21/2005
1366 //-----------------------------------------------------------------------------
1368 {
1369  bool bsuccess = true;
1370 
1371 #ifdef Xyce_DEBUG_DEVICE
1372  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
1373  {
1374  Xyce::dout() << "Instance::loadDAEdFdx----------------------" << std::endl
1375  << "\tname = " << getName() << std::endl;
1376  }
1377 #endif
1378 
1379  N_LAS_Vector & solVector = *(extData.nextSolVectorPtr);
1380  N_LAS_Vector & lastSolVector = *(extData.lastSolVectorPtr);
1381  N_LAS_Matrix * dFdxMatPtr = extData.dFdxMatrixPtr;
1382 
1383  // pull these parameters up from the model class to make it easier
1384  // to view the equations.
1385  const double Gap = model_.Gap;
1386  const double Path = model_.Path;
1387 
1388 
1389  // loop over each inductor and load it's dFdx components
1390 #ifdef MS_FACTORING2
1391  double mid = 1.0 + (1.0 - (Gap/Path))*P*model_.Ms;
1392 #else
1393  double mid = 1.0 + (1.0 - (Gap/Path))*P;
1394 #endif
1395 
1396  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1397  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1398  while( currentInductor != endInductor )
1399  {
1400  // do the normal work for an inductor
1401  (*dFdxMatPtr)[((*currentInductor)->li_Pos)] [((*currentInductor)->APosEquBraVarOffset)] += 1.0;
1402  (*dFdxMatPtr)[((*currentInductor)->li_Neg)] [((*currentInductor)->ANegEquBraVarOffset)] += -1.0;
1403  (*dFdxMatPtr)[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquPosNodeOffset)] += -1.0/mid;
1404  (*dFdxMatPtr)[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquNegNodeOffset)] += 1.0/mid;
1405 
1406  double delV = solVector[(*currentInductor)->li_Pos] - solVector[(*currentInductor)->li_Neg];
1407 
1408  for( int j = 0; j<numInductors; ++j )
1409  {
1410 
1411  (*dFdxMatPtr)[((*currentInductor)->li_Branch)][(*currentInductor)->inductorCurrentOffsets[j]] +=
1412  delV * (1.0 - (Gap/Path)) * dP_dI[j]/(mid*mid);
1413 
1414  }
1415 
1416  currentInductor++;
1417  }
1418 
1419  if(includeDeltaM)
1420  {
1421  (*dFdxMatPtr)[li_deltaMagVar][deltaMEquDeltaMOffset] = 1.0;
1422 
1423  for( int i=0; i<numInductors; i++ )
1424  {
1425  (*dFdxMatPtr)[li_deltaMagVar][deltaMEquInductorOffsets[i]] =
1426  -((inductanceVals[i]*(solVector[ instanceData[i]->li_Branch ] - lastSolVector[ instanceData[i]->li_Branch ] ) * dP_dI[i])
1427  + P*inductanceVals[i])/(model_.Path*model_.Ms );
1428  }
1429  }
1430 
1431  return bsuccess;
1432 }
1433 
1434 //-----------------------------------------------------------------------------
1435 // Function : Instance::outputPlotFiles
1436 // Purpose : If requested by the use in the model statement,
1437 // this routine outputs values of the internal
1438 // state variables M, H and R to a file
1439 // named "Inductor_name.dat". File is opened
1440 // and closed in the contructor and destructor.
1441 // Special Notes :
1442 // Scope : public
1443 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1444 // Creation Date : 03/21/2005
1445 //-----------------------------------------------------------------------------
1447 {
1448  bool bsuccess = true;
1449  if( outputStateVarsFlag && outputFileStreamPtr.get() && (*outputFileStreamPtr) )
1450  {
1451  N_LAS_Vector & solVector = *(extData.nextSolVectorPtr);
1452 
1453 #ifdef MS_FACTORING2
1454  double latestMag = MagVar*model_.Ms;
1455 #else
1456  double latestMag = MagVar;
1457 #endif
1458 
1459  if( includeDeltaM )
1460  {
1461  (*outputFileStreamPtr)
1462  << getSolverState().currTime << " "
1463  << latestMag
1464  << std::endl;
1465  }
1466  else
1467  {
1468  (*outputFileStreamPtr)
1469  << getSolverState().currTime << " "
1470  << latestMag
1471  << std::endl;
1472  }
1473 
1474  }
1475  return bsuccess;
1476 }
1477 
1478 
1479 //-----------------------------------------------------------------------------
1480 // Function : Instance::setIC
1481 // Purpose :
1482 // Special Notes :
1483 // Scope : public
1484 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1485 // Creation Date : 03/21/2005
1486 //-----------------------------------------------------------------------------
1488 {
1489  int i_bra_sol;
1490  int i_f_state;
1491 
1492  bool bsuccess = true;
1493  return bsuccess;
1494 }
1495 
1496 //-----------------------------------------------------------------------------
1497 // Function : Instance::varTypes
1498 // Purpose :
1499 // Special Notes :
1500 // Scope : public
1501 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1502 // Creation Date : 03/21/2005
1503 //-----------------------------------------------------------------------------
1504 void Instance::varTypes( std::vector<char> & varTypeVec )
1505 {
1506  varTypeVec.resize(numInductors);
1507  for(int i=0; i<numInductors; i++)
1508  {
1509  varTypeVec[i] = 'I';
1510  }
1511 }
1512 
1513 
1514 
1515 //-----------------------------------------------------------------------------
1516 // Function : Instance::Pcalc
1517 // Purpose :
1518  // this is a templated function for a complicated term P(M,I_1... I_n) that relates
1519  // the magnetic saturation of the mutual indcutor to the individual currents
1520  // through the inductors. We'll need dP_dM and this tempated function automates
1521  // that calculation via Sacado
1522 // Special Notes :
1523 // Scope : public
1524 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1525 // Creation Date : 10/13/2011
1526 //-----------------------------------------------------------------------------
1527 template <typename ScalarT>
1529  const ScalarT & Mag, const ScalarT & CurrentSum,
1530  const ScalarT & Vpos, const ScalarT & Vneg)
1531  // independent variable M
1532  // independent variable Sum(n_i * I_i) windings * current through each inductor
1533  // independent variable Vpos, Vneg -- voltage drop over first inductor
1534 {
1535  // some parameters in the model class that we will use often
1536  const double A = model_.A;
1537  const double Alpha = model_.Alpha;
1538  const double Area = model_.Area;
1539  const double BetaH = model_.BetaH;
1540  const double BetaM = model_.BetaM;
1541  const double C = model_.C;
1542  const double DeltaV = model_.DeltaV;
1543  const double Gap = model_.Gap;
1544  const double Ms = model_.Ms;
1545  const double Kirr = model_.Kirr;
1546  const double Path = model_.Path;
1547  const double Vinf = model_.Vinf;
1548 
1549  ScalarT qV = (DeltaV / Vinf) * (Vpos - Vneg);
1550 
1551  ScalarT tanh_qV = 0.0;
1552  if (fabs(qV) < CONSTTANH_THRESH)
1553  {
1554  ScalarT tanh_qV = tanh(qV);
1555  }
1556  else if (qV < 0.0)
1557  {
1558  tanh_qV = -1.0;
1559  }
1560  else
1561  {
1562  tanh_qV = 1.0;
1563  }
1564 
1565  ScalarT Happ = CurrentSum / Path;
1566 
1567 #ifdef MS_FACTORING2
1568  ScalarT H = Happ - (Gap / Path) * Mag * Ms;
1569  ScalarT He = H + Alpha * Mag * Ms;
1570 #else
1571  ScalarT H = Happ - (Gap / Path) * Mag;
1572  ScalarT He = H + Alpha * Mag;
1573 #endif
1574 
1575  ScalarT Heo = BetaH*A;
1576 
1577  // terms that come up frequently
1578  ScalarT gap_path = Gap / Path;
1579  ScalarT He2 = He*He;
1580  ScalarT Heo2 = Heo*Heo;
1581  ScalarT sq_Heo2He2 = sqrt(Heo2 + He2);
1582 
1583  ScalarT delM0 = model_.BetaM * Ms;
1584  ScalarT Man = Ms * He / ( A + sq_Heo2He2 );
1585 #ifdef MS_FACTORING2
1586  ScalarT delM = Man - Mag*Ms;
1587 #else
1588  ScalarT delM = Man - Mag;
1589 #endif
1590 
1591  ScalarT delM2 = delM*delM;
1592  ScalarT delM02 = delM0*delM0;
1593  ScalarT sq_delM02delM2 = sqrt( delM02 + delM2 );
1594 
1595 #ifdef MS_FACTORING2
1596  ScalarT Mirrp = (delM * tanh_qV + sq_delM02delM2 ) / (2*( Kirr- Alpha * sq_delM02delM2));
1597  ScalarT Manp = Ms * (A + Heo2/sq_Heo2He2) / pow(A + sq_Heo2He2, 2.0);
1598  ScalarT Pval = ( C * Manp + (1 - C) * Mirrp) / ((1 + (gap_path - Alpha) * C * Manp + gap_path * (1-C) * Mirrp)*Ms);
1599 #else
1600  ScalarT Mirrp = (delM * tanh_qV + sq_delM02delM2 ) / (2*( Kirr- Alpha * sq_delM02delM2));
1601  ScalarT Manp = Ms*(A + Heo2/sq_Heo2He2) / pow(A + sq_Heo2He2, 2.0);
1602  ScalarT Pval = ( C * Manp + (1 - C) * Mirrp) / (1 + (gap_path - Alpha) * C * Manp + gap_path * (1-C) * Mirrp);
1603 #endif
1604 
1605  return Pval;
1606 } // end of Pcalc() function
1607 
1608 
1609 //-----------------------------------------------------------------------------
1610 // Function : Model::processParams
1611 // Purpose :
1612 // Special Notes :
1613 // Scope : public
1614 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1615 // Creation Date : 03/21/2005
1616 //-----------------------------------------------------------------------------
1618 {
1619  return true;
1620 }
1621 
1622 //----------------------------------------------------------------------------
1623 // Function : Model::processInstanceParams
1624 // Purpose :
1625 // Special Notes :
1626 // Scope : public
1627 // Creator : Dave Shirely, PSSI
1628 // Creation Date : 03/23/06
1629 //----------------------------------------------------------------------------
1631 {
1632  std::vector<Instance*>::iterator iter;
1633  std::vector<Instance*>::iterator first = instanceContainer.begin();
1634  std::vector<Instance*>::iterator last = instanceContainer.end();
1635 
1636  for (iter=first; iter!=last; ++iter)
1637  {
1638  (*iter)->processParams();
1639  }
1640 
1641  return true;
1642 }
1643 
1644 //-----------------------------------------------------------------------------
1645 // Function : Model::Model
1646 // Purpose : block constructor
1647 // Special Notes :
1648 // Scope : public
1649 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1650 // Creation Date : 03/21/2005
1651 //-----------------------------------------------------------------------------
1653  const Configuration & configuration,
1654  const ModelBlock & MB,
1655  const FactoryBlock & factory_block)
1656  : DeviceModel(MB, configuration.getModelParameters(), factory_block),
1657  A(0.0),
1658  Alpha(0.0),
1659  Area(0.0),
1660  BetaH(0.0),
1661  BetaM(0.0),
1662  C(0.0),
1663  DeltaV(0.0),
1664  Gap(0.0),
1665  Kirr(0.0),
1666  Ms(0.0),
1667  Path(0.0),
1668  Vinf(0.0),
1669  tempCoeff1(0.0),
1670  tempCoeff2(0.0),
1671  outputStateVars(0),
1672  useRKIntegration(0),
1673  includeDeltaM(0),
1674  tnom(getDeviceOptions().tnom)
1675 {
1676  setLevel(2);
1677 
1678 
1679  // Set params to constant default values:
1680  setDefaultParams ();
1681 
1682  // Set params according to .model line and constant defaults from metadata:
1683  setModParams (MB.params);
1684 
1685  // scale gap, path and area from cm and cm^2 to m and m^2
1686  Gap *= 1.0e-2;
1687  Path *= 1.0e-2;
1688  Area *= 1.0e-4;
1689 
1690  // Set any non-constant parameter defaults:
1691 
1692  if (!given("TNOM"))
1693  {
1695  }
1696 
1697  // Calculate any parameters specified as expressions:
1699 
1700  // calculate dependent (ie computed) params and check for errors:
1701  processParams ();
1702 }
1703 
1704 //-----------------------------------------------------------------------------
1705 // Function : Model::~Model
1706 // Purpose : destructor
1707 // Special Notes :
1708 // Scope : public
1709 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1710 // Creation Date : 03/21/2005
1711 //-----------------------------------------------------------------------------
1713 {
1714  std::vector<Instance*>::iterator iter;
1715  std::vector<Instance*>::iterator first = instanceContainer.begin();
1716  std::vector<Instance*>::iterator last = instanceContainer.end();
1717 
1718  for (iter=first; iter!=last; ++iter)
1719  {
1720  delete (*iter);
1721  }
1722 }
1723 
1724 // additional Declarations
1725 
1726 //-----------------------------------------------------------------------------
1727 // Function : Model::printOutInstances
1728 // Purpose : debugging tool.
1729 // Special Notes :
1730 // Scope : public
1731 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1732 // Creation Date : 03/21/2005
1733 //-----------------------------------------------------------------------------
1734 std::ostream &Model::printOutInstances(std::ostream &os) const
1735 {
1736  std::vector<Instance*>::const_iterator iter;
1737  std::vector<Instance*>::const_iterator first = instanceContainer.begin();
1738  std::vector<Instance*>::const_iterator last = instanceContainer.end();
1739 
1740  int i, isize;
1741  isize = instanceContainer.size();
1742 
1743  os << std::endl;
1744  os << "Number of MutIndNonLin instances: " << isize << std::endl;
1745  os << " name=\t\tmodelName\tParameters" << std::endl;
1746  for (i=0, iter=first; iter!=last; ++iter, ++i)
1747  {
1748  os << " " << i << ": " << (*iter)->getName() << "\t";
1749  os << getName();
1750  os << std::endl;
1751  }
1752 
1753  os << std::endl;
1754 
1755  return os;
1756 }
1757 
1758 //-----------------------------------------------------------------------------
1759 // Function : Model::forEachInstance
1760 // Purpose :
1761 // Special Notes :
1762 // Scope : public
1763 // Creator : David Baur
1764 // Creation Date : 2/4/2014
1765 //-----------------------------------------------------------------------------
1766 /// Apply a device instance "op" to all instances associated with this
1767 /// model
1768 ///
1769 /// @param[in] op Operator to apply to all instances.
1770 ///
1771 ///
1772 void Model::forEachInstance(DeviceInstanceOp &op) const /* override */
1773 {
1774  for (std::vector<Instance *>::const_iterator it = instanceContainer.begin(); it != instanceContainer.end(); ++it)
1775  op(*it);
1776 }
1777 
1778 
1779 Device *Traits::factory(const Configuration &configuration, const FactoryBlock &factory_block)
1780 {
1781 
1782  return new DeviceMaster<Traits>(configuration, factory_block, factory_block.solverState_, factory_block.deviceOptions_);
1783 }
1784 
1786 {
1788  .registerDevice("min", 2)
1789  .registerModelType("min", 2);
1790 }
1791 
1792 } // namespace MutIndNonLin2
1793 } // namespace Device
1794 } // namespace Xyce