Xyce  6.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
N_DEV_MutIndLin.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_MutIndLin.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.129.2.1 $
40 //
41 // Revision Date : $Date: 2014/08/26 22:11:58 $
42 //
43 // Current Owner : $Author: rlschie $
44 //-------------------------------------------------------------------------
45 #include <Xyce_config.h>
46 
47 // ---------- Standard Includes ----------
48 #include <N_UTL_Misc.h>
49 
50 #ifdef HAVE_CMATH
51 #include <cmath>
52 #else
53 #include <math.h>
54 #endif
55 
56 #include <algorithm>
57 #include <set>
58 
59 // ---------- Xyce Includes ----------
60 #include <N_DEV_DeviceOptions.h>
61 #include <N_DEV_DeviceMaster.h>
62 #include <N_DEV_ExternData.h>
63 #include <N_DEV_MatrixLoadData.h>
64 #include <N_DEV_MutIndLin.h>
65 #include <N_DEV_SolverState.h>
66 #include <N_DEV_Message.h>
67 #include <N_ERH_ErrorMgr.h>
68 
69 #include <N_LAS_Vector.h>
70 #include <N_LAS_Matrix.h>
71 
72 #include <N_UTL_Expression.h>
73 
74 namespace Xyce {
75 namespace Device {
76 
77 
78 //-----------------------------------------------------------------------------
79 // Function : InductorInstanceData::InductorInstanceData
80 // Purpose :
81 // Special Notes :
82 //
83 // Need to have a constructor for this data-only class
84 // or it is treated as in-line by some compilers (gcc 3.3 on the mac)
85 // Trying to inline this could cause problems
86 //
87 // Scope : public
88 // Creator : Rich Schiek
89 // Creation Date :
90 //-----------------------------------------------------------------------------
92  name(""), // name of inductor
93  L(0.0),
94  IC(0.0),
95  ICGiven(false),
96  baseL(0.0),
97  li_Pos(-1),
98  li_Neg(-1),
99  li_Branch(-1),
100  APosEquBraVarOffset(-1),
101  ANegEquBraVarOffset(-1),
102  ABraEquPosNodeOffset(-1),
103  ABraEquNegNodeOffset(-1),
104  ABraEquBraVarOffset(-1),
105  magOffset(-1),
106  vPosOffset(-1),
107  vNegOffset(-1),
109  // Pointers for dFdx entries:
110  f_PosEquBraVarPtr(0),
111  f_NegEquBraVarPtr(0),
112  f_BraEquPosNodePtr(0),
113  f_BraEquNegNodePtr(0),
114  f_BraEquBraVarPtr(0),
115  // offsets only needed in nonlinear application
116  f_magPtr(0),
117  f_vPosPtr(0),
118  f_vNegPtr(0),
119  q_PosEquBraVarPtr(0),
120  q_NegEquBraVarPtr(0),
121  q_BraEquPosNodePtr(0),
122  q_BraEquNegNodePtr(0),
123  q_BraEquBraVarPtr(0),
124  // offsets only needed in nonlinear application
125  q_magPtr(0),
126  q_vPosPtr(0),
127  q_vNegPtr(0)
128 #endif
129 {}
130 
131 namespace MutIndLin {
132 
133 
135 {
136 // Set up double precision variables:
137  p.addPar ("COUP_VAL", 1.0, false, ParameterType::NO_DEP,
140  U_NONE, CAT_NONE, "Coupling value");
141 
142  // Set up exceptions (ie variables that are not doubles):
143  p.addPar ("COUPLEDMutIndLin", std::vector<std::string>(), false, ParameterType::NO_DEP,
145  p.addPar ("COUPLEDINDUCTANCE", std::vector<double>(), false, ParameterType::NO_DEP,
147  p.addPar ("NODE1", std::vector<std::string>(), false, ParameterType::NO_DEP,
149  p.addPar ("NODE2", std::vector<std::string>(), false, ParameterType::NO_DEP,
151  p.addPar ("COUPLING", std::vector<double>(), false, ParameterType::SOLN_DEP,
152  &MutIndLin::Instance::couplingCoefficient, NULL, U_NONE, CAT_NONE, "Coupling coefficient");
153  p.addPar ("COUPLEDINDUCTOR", std::vector<std::string>(), false, ParameterType::NO_DEP,
155 }
156 
158 {
159  p.addPar ("TNOM", 27.0, false, ParameterType::NO_DEP,
161  NULL, U_DEGC, CAT_MATERIAL, "Reference temperature");
162 
163  p.addPar ("TC1", 0.0, false, ParameterType::NO_DEP,
165  NULL, U_NONE, CAT_MATERIAL, "First order temperature coeff.");
166 
167  p.addPar ("TC2", 0.0, false, ParameterType::NO_DEP,
169  NULL, U_NONE, CAT_MATERIAL, "Second order temperature coeff.");
170 }
171 
172 
173 //-----------------------------------------------------------------------------
174 // Function : Instance::Instance
175 // Purpose : constructor
176 // Special Notes :
177 // Scope : public
178 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
179 // Creation Date : 03/21/2005
180 //-----------------------------------------------------------------------------
182  const Configuration & configuration,
183  const InstanceBlock & IB,
184  Model & Iiter,
185  const FactoryBlock & factory_block)
186  : DeviceInstance(IB, configuration.getInstanceParameters(), factory_block),
187  model_(Iiter),
188  numInductors(0),
189  mutualCup(0.0),
190  mutualCupGiven(false),
191  scalingRHS(1.0),
192  temp(getDeviceOptions().temp.getImmutableValue<double>()),
193  tempGiven(false)
194 {
195  scalingRHS = 1.0;
196 
197  // set some default values. May be changed by processParams
198  numExtVars = 2;
199  numIntVars = 1;
200  numStateVars = 0;
201 
202  // Set params to constant default values:
203  setDefaultParams ();
204 
205  // Set params according to instance line and constant defaults from metadata:
206  setParams (IB.params);
207 
208  // now load the instance data vector
209  for( int i=0; i<inductorNames.size(); ++i )
210  {
211  InductorInstanceData * inductorData = new InductorInstanceData();
212  inductorData->name = inductorNames[i];
213  inductorData->L = inductorInductances[i];
214  inductorData->baseL = inductorInductances[i];
215  inductorData->ICGiven = false;
216  inductorData->inductorCurrentOffsets.resize( inductorNames.size() );
217  inductorData->f_inductorCurrentPtrs.resize( inductorNames.size() );
218  inductorData->q_inductorCurrentPtrs.resize( inductorNames.size() );
219  inductorData->depVarPairs.clear();
220  instanceData.push_back( inductorData );
221  }
222 
223  numInductors = instanceData.size();
224 
225  // set up the device connectivity map
226  // each simple inductor in this mutual inductor
227  // is maked as a connection (given a common, non-zero
228  // value in devConMap)
229  devConMap.resize(2*numInductors);
230  for(int i=0; i<numInductors; i++)
231  {
232  devConMap[i] = devConMap[i+1] = (i+1);
233  }
234 
235  // now assemble the mutual coupling coefficient matrix
236  // assume no coupling except to self, so just ones on dialgonal.
237  mutualCouplingCoef.resize( numInductors );
238  mutualCouplingCoefDerivs.resize( numInductors );
239  for( int i=0; i<numInductors; ++i )
240  {
241  mutualCouplingCoef[i].resize( numInductors );
242  fill( mutualCouplingCoef[i].begin(), mutualCouplingCoef[i].end(), 0.0 );
243  mutualCouplingCoef[i][i] = 1.0;
244 
245  mutualCouplingCoefDerivs[i].resize( numInductors );
246  fill( mutualCouplingCoefDerivs[i].begin(),
247  mutualCouplingCoefDerivs[i].end(), 0.0 );
248  }
249 
250  // take the vector couplingInductor in pairs to figure out
251  // the coupling value between inductors then insert it into
252  // the couplingCoefficient matrix.
253  indexPairs.resize(couplingInductor.size()/2);
254  for( int i=0; i<(couplingInductor.size()-1) ; i+=2 )
255  {
256  // from inductor names, derive indices.
257  int indexInductorOne = -1;
258  int indexInductorTwo = -1;
259  for( int j=0; j<inductorNames.size(); ++j)
260  {
261  if( couplingInductor[i] == inductorNames[j] )
262  {
263  indexInductorOne = j;
264  }
265  if( couplingInductor[i+1] == inductorNames[j] )
266  {
267  indexInductorTwo = j;
268  }
269  if( indexInductorOne != -1 && indexInductorTwo != -1 )
270  {
271  // stop searching if we found the answers.
272  break;
273  }
274  }
275  // Save the indices for later
276  indexPairs[i/2].first = indexInductorOne;
277  indexPairs[i/2].second = indexInductorTwo;
278 
279  // with the indices at hand, assign the value.
280  // note that couplingCoefficient can be of length 1 if all coefficients were the
281  // same. Catch that case here.
282  if( (i/2) < couplingCoefficient.size() )
283  {
284  mutualCouplingCoef[ indexInductorOne ][ indexInductorTwo ] =
285  mutualCouplingCoef[ indexInductorTwo ][ indexInductorOne ] = couplingCoefficient[ (i / 2) ];
286  }
287  else
288  {
289  mutualCouplingCoef[ indexInductorOne ][ indexInductorTwo ] =
290  mutualCouplingCoef[ indexInductorTwo ][ indexInductorOne ] = couplingCoefficient[ 0 ];
291  }
292  }
293 
294 #ifdef Xyce_DEBUG_DEVICE
295  if (getDeviceOptions().debugLevel > 0)
296  {
297  Xyce::dout() << "mutualCouplingCoef: " << std::endl;
298  for( int i=0; i<numInductors; ++i )
299  {
300  for( int j=0; j<numInductors; ++j )
301  {
302  Xyce::dout() << mutualCouplingCoef[i][j] << " ";
303  }
304  Xyce::dout() << std::endl;
305  }
306  }
307 #endif
308 
309  inductorCurrents.resize( numInductors );
310  dIdt.resize( numInductors );
311 
312  inductanceVals.resize( numInductors );
313  LOI.resize( numInductors );
314  LO.resize( numInductors );
315  for( int i=0; i<numInductors; ++i)
316  {
317  LO[i].resize( numInductors );
318  }
319 
321 
322  // Calculate any parameters specified as expressions:
324 
325  // calculate dependent (ie computed) params and check for errors:
326  processParams ();
327 
328  // update internal/external/state variable counts
331  numStateVars = 0;
332 
333  numCoupDepVars.resize(couplingCoefficient.size(),0);
334  expPtrs.resize(couplingCoefficient.size(),0); // null pointers
336 
337  // set up the jacobian stamp
338  // for an individual inductor the stamp would be:
339  //
340  // V1 V2 Ib
341  // kcl1 1
342  // kcl2 -1
343  // branch 1 -1 L/dt
344  //
345  // for a collection of these, the internal variable, branch equations,
346  // must be at the end of a given stamp row.
347  //
348  // So for N inductors the samp is:
349  //
350  // V1 V2 V3 V4 ... V2N I1 I2 ... IN
351  // kcl1 1
352  // kcl2 -1
353  // kcl3 1
354  // kcl4 -1
355  // branch1 1 -1 L/dt c ... c
356  // branch2 1 -1 c L/dt ... c
357  //
358  // where "c" is an induced current change.
359 
360  jacStamp.resize( 3 * numInductors );
361  for( int i=0; i< numInductors; ++i )
362  {
363  jacStamp[2*i].resize(1); // Vpos row
364  jacStamp[2*i+1].resize(1); // Vneg row
365  jacStamp[2*numInductors + i].resize(numInductors + 2); // Ibranch row
366 
367  jacStamp[2*i ][0] = 2*numInductors + i; // vpos-ibranch
368  jacStamp[2*i+1][0] = 2*numInductors + i; // vneg-ibranch
369  jacStamp[2*numInductors + i][0] = 2*i; // ibranch-vpos
370  jacStamp[2*numInductors + i][1] = 2*i + 1; // ibranch-vneg
371  for( int j=0; j<numInductors; ++j )
372  {
373  jacStamp[2*numInductors + i][j+2] = 2*numInductors + j;
374  // ibranch-ibranch[j]
375  }
376  }
377 
378  // Now we process our dependent variables.
379  // We need to bump up the jacstamp for all the dependencies on
380  // the coupling vars that are expressions. We also need to keep track of
381  // how many variables each inductor depends on.
382  std::vector<Depend>::const_iterator d;
383  std::vector<Depend>::const_iterator begin=getDependentParams().begin();
384  std::vector<Depend>::const_iterator end=getDependentParams().end();
385 
386  for (d=begin; d != end; ++d)
387  {
388 
389  if (d->name == "COUPLING" && d->vectorIndex != -1)
390  {
391  // sanity check, they should all have the first two conditions,
392  // and the next means we actually have work to do
393  if (d->n_vars >0)
394  {
395  // Keep track of the number of variables this coefficient depends on:
396  numCoupDepVars[d->vectorIndex] = d->n_vars;
397 
398  // Save the pointer to the expresison, we need to be able to evaluate it.
399  expPtrs[d->vectorIndex] = d->expr;
400 
401  // These are the two inductors that this coefficient couples
402  int indexInductorOne = indexPairs[d->vectorIndex].first;
403  int indexInductorTwo = indexPairs[d->vectorIndex].second;
404 
405  // The coupling coefficient introduces terms to both equations
406  int inductorOneBranchSize = jacStamp[2*numInductors+
407  indexInductorOne].size();
408  int inductorTwoBranchSize = jacStamp[2*numInductors+
409  indexInductorTwo].size();
410  jacStamp[2*numInductors+indexInductorOne].resize(inductorOneBranchSize+
411  d->n_vars);
412  jacStamp[2*numInductors+indexInductorTwo].resize(inductorTwoBranchSize+
413  d->n_vars);
414  for (int i=0; i<d->n_vars; i++)
415  {
416  jacStamp[2*numInductors+indexInductorOne][inductorOneBranchSize+i]
417  = d->lo_var+i+3*numInductors;
418  jacStamp[2*numInductors+indexInductorTwo][inductorTwoBranchSize+i]
419  = d->lo_var+i+3*numInductors;
420 
421  // We need to be able to determine which of our dependent var
422  // indices map onto which variable of what coupling coefficient
423  // when we go to do the jacobian loads. We do this so we don't
424  // need to keep all kinds of two-dimensional arrays like
425  // mutualCoupling that are just copies of the one-dimensonal
426  // arrays like couplingCoefficient.
427  instanceData[indexInductorOne]->depVarPairs.push_back(
428  std::pair<int,int>(d->vectorIndex,i));
429  instanceData[indexInductorTwo]->depVarPairs.push_back(
430  std::pair<int,int>(d->vectorIndex,i));
431  }
432  }
433  }
434  else
435  {
436  std::string msg="Error in mutual inductor constructor for ";
437  msg += getName().getEncodedName();
438  msg += ": Found a dependent parameter that fails sanity check.";
439  msg += "Name:" + d->name + " ";
440  msg += "Value: " + (d->expr)->get_expression() + " ";
441  N_ERH_ErrorMgr::report(N_ERH_ErrorMgr::DEV_FATAL, msg);
442  }
443  }
444 
445 
446  // We're now done processing all the dependent variables, so we know
447  // everything we need to know about expressions and derivative numbers.
448  for (int i=0; i<couplingCoefficient.size(); ++i)
449  {
451  }
452 
453 #ifdef Xyce_DEBUG_DEVICE
454  if (getDeviceOptions().debugLevel > 0)
455  {
456  Xyce::dout() << "Instance::Instance-----------------" << std::endl;
457  Xyce::dout() << "numExtVars = " << numExtVars << std::endl
458  << "numIntVars = " << numIntVars << std::endl
459  << "numStateVars = " << numStateVars << std::endl
460  << "numInductors = " << numInductors << std::endl
461  << "jacStamp = " << std::endl;
462  for( int i = 0; i<jacStamp.size(); ++i )
463  {
464  Xyce::dout() << "jacStamp[ " << i << " ] = { ";
465  for( int j=0; j<jacStamp[i].size(); ++j)
466  {
467  Xyce::dout() << jacStamp[i][j];
468  if( j != ( jacStamp[i].size() -1 ) )
469  {
470  Xyce::dout() << ", ";
471  }
472  }
473  Xyce::dout() << " }" << std::endl;
474  }
475  }
476 #endif
477 
478 }
479 
480 //-----------------------------------------------------------------------------
481 // Function : Instance::~Instance
482 // Purpose : destructor
483 // Special Notes :
484 // Scope : public
485 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
486 // Creation Date : 03/21/2005
487 //-----------------------------------------------------------------------------
489 {
490  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
491  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
492  for ( ; currentInductor != endInductor ; ++currentInductor)
493  {
494  delete *currentInductor;
495  }
496 }
497 
498 //-----------------------------------------------------------------------------
499 // Function : Instance::registerLIDs
500 // Purpose :
501 // Special Notes :
502 // Scope : public
503 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
504 // Creation Date : 03/21/2005
505 //-----------------------------------------------------------------------------
506 void Instance::registerLIDs(const std::vector<int> & intLIDVecRef,
507  const std::vector<int> & extLIDVecRef)
508 {
509  AssertLIDs(intLIDVecRef.size() == numIntVars);
510  AssertLIDs(extLIDVecRef.size() == numExtVars);
511 
512 // copy over the global ID lists.
513  intLIDVec = intLIDVecRef;
514  extLIDVec = extLIDVecRef;
515 
516  // Now use these lists to obtain the indices into the
517  // linear algebra entities. This assumes an order.
518  // For the matrix indices, first do the rows.
519  // get the current values of the inductances and currentOffsets
520  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
521  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
522  int i = 0;
523  int j = 0;
524  while( currentInductor != endInductor )
525  {
526  (*currentInductor)->li_Pos = extLIDVec[ i++ ];
527  (*currentInductor)->li_Neg = extLIDVec[ i++ ];
528  (*currentInductor)->li_Branch = intLIDVec[ j++ ];
529  ++currentInductor;
530  }
531 
532 #ifdef Xyce_DEBUG_DEVICE
533  if (getDeviceOptions().debugLevel > 0)
534  {
535  Xyce::dout() << "Instance::registerLIDs----------------------------" << std::endl;
536  currentInductor = instanceData.begin();
537  i=0;
538  while( currentInductor != endInductor )
539  {
540  Xyce::dout() << "Inductor [ " << i++ << " ] "
541  << " li_Pos = " << (*currentInductor)->li_Pos
542  << " li_Neg = " << (*currentInductor)->li_Neg
543  << " li_Branch = " << (*currentInductor)->li_Branch << std::endl;
544  ++currentInductor;
545  }
546  }
547 #endif
548 }
549 
550 //-----------------------------------------------------------------------------
551 // Function : DiodeInstance::getIntNameMap
552 // Purpose :
553 // Special Notes :
554 // Scope : public
555 // Creator : Dave Shirley, PSSI
556 // Creation Date : 04/26/06
557 //-----------------------------------------------------------------------------
558 std::map<int,std::string> & Instance::getIntNameMap ()
559 {
560  // set up the internal name map, if it hasn't been already.
561  if (intNameMap.empty ())
562  {
563  for( int i=0; i<instanceData.size() ; ++i)
564  {
565  intNameMap[ intLIDVec[i] ] = spiceInternalName(getName(), instanceData[i]->name + "_branch");
566  }
567  }
568 
569  return intNameMap;
570 }
571 
572 //-----------------------------------------------------------------------------
573 // Function : Instance::registerStateLIDs
574 // Purpose :
575 // Special Notes :
576 // Scope : public
577 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
578 // Creation Date : 03/21/2005
579 //-----------------------------------------------------------------------------
580 void Instance::registerStateLIDs( const std::vector<int> & staLIDVecRef )
581 {
582  AssertLIDs(staLIDVecRef.size() == numStateVars);
583 
584  // copy over the global ID lists.
585  staLIDVec = staLIDVecRef;
586 }
587 
588 //-----------------------------------------------------------------------------
589 // Function : BsrcInstance::getDepSolnVars
590 // Purpose :
591 // Special Notes :
592 // Scope : public
593 // Creator : Rob Hoekstra, SNL, Parallel Computational Sciences
594 // Creation Date : 06/06/01
595 //-----------------------------------------------------------------------------
596 const std::vector<std::string> & Instance::getDepSolnVars()
597 {
599 }
600 
601 //-----------------------------------------------------------------------------
602 // Function : Instance::jacobianStamp
603 // Purpose :
604 // Special Notes :
605 // Scope : public
606 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
607 // Creation Date : 03/21/2005
608 //-----------------------------------------------------------------------------
609 const std::vector< std::vector<int> > & Instance::jacobianStamp() const
610 {
611  return jacStamp;
612 }
613 
614 //-----------------------------------------------------------------------------
615 // Function : Instance::registerJacLIDs
616 // Purpose :
617 // Special Notes :
618 // Scope : public
619 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
620 // Creation Date : 03/21/2005
621 //-----------------------------------------------------------------------------
622 void Instance::registerJacLIDs( const std::vector< std::vector<int> > & jacLIDVec )
623 {
624  DeviceInstance::registerJacLIDs( jacLIDVec );
625 #ifdef Xyce_DEBUG_DEVICE
626  if (getDeviceOptions().debugLevel > 0)
627  {
628  Xyce::dout() << "Instance::registerJacLIDs ----------------------------" << std::endl;
629  }
630 #endif
631  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
632  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
633  int i = 0;
634  while( currentInductor != endInductor )
635  {
636  (*currentInductor)->APosEquBraVarOffset = jacLIDVec[ 2*i ][ 0 ];
637  (*currentInductor)->ANegEquBraVarOffset = jacLIDVec[ 2*i + 1 ][ 0 ];
638  (*currentInductor)->ABraEquPosNodeOffset = jacLIDVec[ 2*numInductors + i ][ 0 ];
639  (*currentInductor)->ABraEquNegNodeOffset = jacLIDVec[ 2*numInductors + i ][ 1 ];
640  for( int j=0; j<numInductors; ++j )
641  {
642  if( i == j )
643  {
644  (*currentInductor)->ABraEquBraVarOffset = jacLIDVec[ 2*numInductors + i ][ j + 2 ];
645  }
646  (*currentInductor)->inductorCurrentOffsets[ j ] = jacLIDVec[ 2*numInductors + i ][ j + 2 ];
647  }
648  // Now do the parts that are for the dependent variables
649  int numdepvars=(*currentInductor)->depVarPairs.size();
650  (*currentInductor)->ABraEquDepVarOffsets.resize(numdepvars);
651  for (int j = 0; j < numdepvars; ++j)
652  {
653  (*currentInductor)->ABraEquDepVarOffsets[j] =
654  jacLIDVec[2*numInductors+i][numInductors+2+j];
655  }
656  ++currentInductor;
657  ++i;
658  }
659 
660 #ifdef Xyce_DEBUG_DEVICE
661  if (getDeviceOptions().debugLevel > 0)
662  {
663  Xyce::dout() << "Instance::registerJacLIDs--------------------------" << std::endl;
664  currentInductor = instanceData.begin();
665  i=0;
666  while( currentInductor != endInductor )
667  {
668  Xyce::dout() << "Inductor [ " << i << " ] " << (*currentInductor)->name
669  << " APosEquBraVarOffset = " << (*currentInductor)->APosEquBraVarOffset
670  << " ANegEquBraVarOffset = " << (*currentInductor)->ANegEquBraVarOffset
671  << " ABraEquPosNodeOffset = " << (*currentInductor)->ABraEquPosNodeOffset
672  << " ABraEquNegNodeOffset = " << (*currentInductor)->ABraEquNegNodeOffset
673  << " ABraEquBraVarOffset = " << (*currentInductor)->ABraEquBraVarOffset << std::endl;
674  Xyce::dout() << "\tInductor branch offsets = { ";
675  for( int j=0; j<numInductors ; ++j )
676  {
677  Xyce::dout() << (*currentInductor)->inductorCurrentOffsets[ j ] << ", ";
678  }
679  Xyce::dout() << "} " << std::endl;
680  ++i;
681  ++currentInductor;
682  }
683  }
684 #endif
685 }
686 
687 //-----------------------------------------------------------------------------
688 // Function : Instance::setupPointers
689 // Purpose :
690 // Special Notes :
691 // Scope : public
692 // Creator : Eric Keiter, SNL
693 // Creation Date : 12/12/08
694 //-----------------------------------------------------------------------------
696 {
697 #ifndef Xyce_NONPOINTER_MATRIX_LOAD
698  N_LAS_Matrix & dFdx = *(extData.dFdxMatrixPtr);
699  N_LAS_Matrix & dQdx = *(extData.dQdxMatrixPtr);
700 
701  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
702  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
703  int i = 0;
704  while( currentInductor != endInductor )
705  {
706  // dFdx pointers:
707  (*currentInductor)->f_PosEquBraVarPtr =
708  &(dFdx[((*currentInductor)->li_Pos)] [((*currentInductor)->APosEquBraVarOffset)] );
709 
710  (*currentInductor)->f_NegEquBraVarPtr =
711  &(dFdx[((*currentInductor)->li_Neg)] [((*currentInductor)->ANegEquBraVarOffset)] );
712 
713  (*currentInductor)->f_BraEquPosNodePtr =
714  &(dFdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquPosNodeOffset)] );
715 
716  (*currentInductor)->f_BraEquNegNodePtr =
717  &(dFdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquNegNodeOffset)] );
718 
719  for( int j=0; j<numInductors; ++j )
720  {
721  if( i == j )
722  {
723  //(*currentInductor)->ABraEquBraVarOffset =
724  // Is this only used for JMat? (ie old DAE)?
725 
726  }
727  (*currentInductor)->f_inductorCurrentPtrs[ j ] =
728  &(dFdx[((*currentInductor)->li_Branch)][(*currentInductor)->inductorCurrentOffsets[j]] );
729  }
730  // Now do the parts that are for the dependent variables
731  int numdepvars=(*currentInductor)->depVarPairs.size();
732  (*currentInductor)->f_BraEquDepVarPtrs.resize(numdepvars);
733 
734  for (int j = 0; j < numdepvars; ++j)
735  {
736  (*currentInductor)->f_BraEquDepVarPtrs[j] =
737  &(dFdx[((*currentInductor)->li_Branch)][(*currentInductor)->ABraEquDepVarOffsets[j]] );
738 
739  }
740 
741  // dQdx pointers:
742  (*currentInductor)->q_PosEquBraVarPtr =
743  &(dQdx[((*currentInductor)->li_Pos)] [((*currentInductor)->APosEquBraVarOffset)] );
744 
745  (*currentInductor)->q_NegEquBraVarPtr =
746  &(dQdx[((*currentInductor)->li_Neg)] [((*currentInductor)->ANegEquBraVarOffset)] );
747 
748  (*currentInductor)->q_BraEquPosNodePtr =
749  &(dQdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquPosNodeOffset)] );
750 
751  (*currentInductor)->q_BraEquNegNodePtr =
752  &(dQdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquNegNodeOffset)] );
753 
754  for( int j=0; j<numInductors; ++j )
755  {
756  if( i == j )
757  {
758  //(*currentInductor)->ABraEquBraVarOffset =
759  // Is this only used for JMat? (ie old DAE)?
760 
761  }
762  (*currentInductor)->q_inductorCurrentPtrs[ j ] =
763  &(dQdx[((*currentInductor)->li_Branch)][(*currentInductor)->inductorCurrentOffsets[j]] );
764  }
765  // Now do the parts that are for the dependent variables
766  numdepvars=(*currentInductor)->depVarPairs.size();
767  (*currentInductor)->q_BraEquDepVarPtrs.resize(numdepvars);
768 
769  for (int j = 0; j < numdepvars; ++j)
770  {
771  (*currentInductor)->q_BraEquDepVarPtrs[j] =
772  &(dQdx[((*currentInductor)->li_Branch)][(*currentInductor)->ABraEquDepVarOffsets[j]] );
773 
774  }
775 
776  ++currentInductor;
777  ++i;
778  }
779 
780 
781 
782 #endif
783 }
784 
785 //-----------------------------------------------------------------------------
786 // Function : Instance::processParams
787 // Purpose :
788 // Special Notes :
789 // Scope : public
790 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
791 // Creation Date : 03/21/2005
792 //-----------------------------------------------------------------------------
794 {
795  int i, j;
796 
797  // take the vector couplingInductor in pairs to figure out
798  // the coupling value between inductors then insert it into
799  // the couplingCoefficient matrix.
800  j = indexPairs.size();
801  for( i=0; i<j ; ++i )
802  {
803  // note that couplingCoefficient can be of length 1 if all coefficients were the
804  // same. Catch that case here.
805  if( i < couplingCoefficient.size() )
806  {
807  mutualCouplingCoef[ indexPairs[i].first ][ indexPairs[i].second ] =
808  mutualCouplingCoef[ indexPairs[i].second ][ indexPairs[i].first ] = couplingCoefficient[i];
809  }
810  else
811  {
812  mutualCouplingCoef[ indexPairs[i].first ][ indexPairs[i].second ] =
813  mutualCouplingCoef[ indexPairs[i].second ][ indexPairs[i].first ] = couplingCoefficient[0];
814  }
815  }
816 #ifdef Xyce_DEBUG_DEVICE
818  {
819  Xyce::dout() << "processParams: mutualCouplingCoef: " << std::endl;
820  for( int i=0; i<numInductors; ++i )
821  {
822  for( int j=0; j<numInductors; ++j )
823  {
824  Xyce::dout() << mutualCouplingCoef[i][j] << " ";
825  }
826  Xyce::dout() << std::endl;
827  }
828  }
829 #endif
830 
831  // set the temperature related stuff.
833 
834  return true;
835 }
836 
837 //-----------------------------------------------------------------------------
838 // Function : Instance::updateTemperature
839 // Purpose :
840 // Special Notes :
841 // Scope : public
842 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
843 // Creation Date : 03/21/2005
844 //-----------------------------------------------------------------------------
845 bool Instance::updateTemperature ( const double & temp)
846 {
847  bool bsuccess = true;
848 
849  // current temp difference from reference temp.
850  double difference = temp - model_.tnom;
851 
852  std::vector< InductorInstanceData* >::iterator currentData = instanceData.begin();
853  while( currentData != instanceData.end() )
854  {
855  double factor = 1.0 + (model_.tempCoeff1)*difference +
856  (model_.tempCoeff2)*difference*difference;
857  (*currentData)->L = ((*currentData)->baseL)*factor;
858  ++currentData;
859  }
860 
861  // now that the inductances have changed we need to update the matrix.
863 
864  return bsuccess;
865 }
866 
867 //-----------------------------------------------------------------------------
868 // Function : Instance::updateIntermediateVars
869 // Purpose :
870 // Special Notes :
871 // Scope : public
872 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
873 // Creation Date : 03/21/2005
874 //-----------------------------------------------------------------------------
876 {
877  bool bsuccess = true;
878 
879  // This method is never called anymore
880 
881  return bsuccess;
882 }
883 
884 //-----------------------------------------------------------------------------
885 // Function : Instance::updateInductanceMatrix()
886 // Purpose : A matrix of inductances is used often enough that it
887 // calculated and stored as a member variable here
888 // If and inductance ever changes say from updating
889 // the temperature or a parameter udpate, then this
890 // routine must be called again.
891 // Special Notes :
892 // Scope : public
893 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
894 // Creation Date : 03/21/2005
895 //-----------------------------------------------------------------------------
897 {
898  std::vector< InductorInstanceData* >::iterator
899  currentInductor = instanceData.begin();
900  std::vector< InductorInstanceData* >::iterator
901  endInductor = instanceData.end();
902 
903  // collec the inductances
904  int i=0;
905  while( currentInductor != endInductor )
906  {
907  inductanceVals[ i ] = ((*currentInductor)->L);
908  ++i;
909  ++currentInductor;
910  }
911 
912  // compute the inductance matrix
913  for( i=0; i<numInductors; ++i)
914  {
915  for( int j=0; j<numInductors; ++j)
916  {
917  LO[i][j] = sqrt( inductanceVals[i]*inductanceVals[j] );
918  }
919  }
920 
921 }
922 //-----------------------------------------------------------------------------
923 // Function : Instance::updatePrimaryState
924 // Purpose :
925 // Special Notes :
926 // Scope : public
927 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
928 // Creation Date : 03/21/2005
929 //-----------------------------------------------------------------------------
931 {
932  // nothing happens in updateIntermediateVars() anymore
933  // bsuccess = updateIntermediateVars ();
934 
935 #ifdef Xyce_DEBUG_DEVICE
936  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
937  {
938  Xyce::dout() << "Instance::updatePrimaryState------------------" << std::endl
939  << "\tname = " << getName() << std::endl;
940  }
941 #endif
942 
943  double * solVec = extData.nextSolVectorRawPtr;
944 
945  // Evaluate the derivatives of all (dependent) coupling coefficients w.r.t
946  // their variables. We need these for both new and old DAE.
947  int ncoupcoef=couplingCoefficient.size();
948  for (int i=0; i<ncoupcoef; ++i)
949  {
950  if (expPtrs[i])
951  {
952  double junk;
953  expPtrs[i]->evaluate( junk, couplingCoefficientVarDerivs[i]);
954  }
955  }
956 
957  // get the currents in each inductor
958  std::vector< InductorInstanceData* >::iterator
959  currentInductor = instanceData.begin();
960  std::vector< InductorInstanceData* >::iterator
961  endInductor = instanceData.end();
962  int i = 0;
963  while( currentInductor != endInductor )
964  {
965  if( (getSolverState().dcopFlag) && ((*currentInductor)->ICGiven) )
966  {
967  inductorCurrents[ i ] = (*currentInductor)->IC;
968  }
969  else
970  {
971  inductorCurrents[ i ] = solVec[ (*currentInductor)->li_Branch ];
972  }
973  ++i;
974  ++currentInductor;
975  }
976 
977  return true;
978 }
979 
980 //-----------------------------------------------------------------------------
981 // Function : Instance::loadDeviceMask
982 //
983 // Purpose : Loads the zero elements of the device mask
984 //
985 // Special Notes : elements of the error vector associated with zero
986 // elements of the mask will not be included in weighted
987 // norms by the time integrator.
988 //
989 // Scope : public
990 // Creator : Richard Schiek, SNL, Electrical and Microsystems Modeling
991 // Creation Date : 02/06/07
992 //-----------------------------------------------------------------------------
994 {
995  bool returnVal=false;
996 #ifndef Xyce_NO_MUTINDLIN_MASK
997  N_LAS_Vector * maskVectorPtr = extData.deviceMaskVectorPtr;
998 
999  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1000  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1001  while( currentInductor != endInductor )
1002  {
1003  (*maskVectorPtr)[((*currentInductor)->li_Branch)] = 0.0;
1004  ++currentInductor;
1005  }
1006  returnVal = true;
1007 #endif
1008  return (returnVal);
1009 }
1010 
1011 //-----------------------------------------------------------------------------
1012 // Function : Instance::loadDAEQVector
1013 //
1014 // Purpose : Loads the Q-vector contributions for a single
1015 // Mutual Inductor instance.
1016 //
1017 // Special Notes : The "Q" vector is part of a standard DAE formalism in
1018 // which the system of equations is represented as:
1019 //
1020 // f(x) = dQ(x)/dt + F(x) - B(t) = 0
1021 //
1022 // Scope : public
1023 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1024 // Creation Date : 03/21/2005
1025 //-----------------------------------------------------------------------------
1027 {
1028 #ifdef Xyce_DEBUG_DEVICE
1029  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
1030  {
1031  Xyce::dout() << "Instance::loadDAEQVector---------------------------" << std::endl
1032  << "\tname = " << getName() << std::endl;
1033  }
1034 #endif
1035  double * qVec = extData.daeQVectorRawPtr;
1036 
1037  // calculate the following product
1038  // I = column vector of currents
1039  // L = row vector of inductances
1040  // LO = matrix = mutualCup * sqrt( L' * L )
1041  // LOI = column vector = mutualCup * sqrt( L' * L ) * I
1042  // LOI[1] = mutualCup * sqrt(L[1]*L[1])*I[1]) +
1043  // mutualCup * sqrt(L[1]*L[2])*I[2]) + ...
1044  // mutualCup * sqrt(L[1]*L[n])*I[n])
1045 
1046  for( int i = 0; i < numInductors; ++i )
1047  {
1048  LOI[ i ] = 0.0;
1049  for( int j = 0; j < numInductors; ++j )
1050  {
1051  LOI[i] += mutualCouplingCoef[i][j] * LO[i][j] * inductorCurrents[j];
1052  }
1053  }
1054 
1055  // loop over each inductor and load it's Q vector components
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  qVec[(*currentInductor)->li_Branch] += LOI[ i ];
1062  ++i;
1063  ++currentInductor;
1064  }
1065 
1066  return true;
1067 }
1068 
1069 //-----------------------------------------------------------------------------
1070 // Function : Instance::loadDAEFVector
1071 //
1072 // Purpose : Loads the F-vector contributions for a single
1073 // Mutual Inductor instance.
1074 //
1075 // Special Notes :
1076 //
1077 // Scope : public
1078 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1079 // Creation Date : 03/21/2005
1080 //-----------------------------------------------------------------------------
1082 {
1083 #ifdef Xyce_DEBUG_DEVICE
1084  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
1085  {
1086  Xyce::dout() << "Instance::loadDAEFVector---------------------------" << std::endl
1087  << "\tname = " << getName() << std::endl;
1088  }
1089 #endif
1090 
1091  double * fVec = extData.daeFVectorRawPtr;
1092  double * solVec = extData.nextSolVectorRawPtr;
1093 
1094  // loop over each inductor and load it's F vector components
1095  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1096  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1097  int i = 0;
1098 
1099  while( currentInductor != endInductor )
1100  {
1101  double current = solVec[(*currentInductor)->li_Branch];
1102  double vNodePos = solVec[(*currentInductor)->li_Pos];
1103  double vNodeNeg = solVec[(*currentInductor)->li_Neg];
1104  fVec[((*currentInductor)->li_Pos)] += scalingRHS * current;
1105  fVec[((*currentInductor)->li_Neg)] += -scalingRHS * current;
1106  fVec[((*currentInductor)->li_Branch)] += -(vNodePos - vNodeNeg);
1107 
1108 #ifdef Xyce_DEBUG_DEVICE
1109  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
1110  {
1111  Xyce::dout() << " Inductor = " << (*currentInductor)->name
1112  << "\tcurrent = " << current
1113  << "\tvNodePos = " << vNodePos
1114  << "\tvNodeNeg = " << vNodeNeg
1115  << std::endl;
1116  }
1117 #endif
1118 
1119  ++currentInductor;
1120  ++i;
1121  }
1122 
1123  return true;
1124 }
1125 
1126 //-----------------------------------------------------------------------------
1127 // Function : Instance::loadDAEdQdx
1128 //
1129 // Purpose : Loads the Q-vector contributions for a single
1130 // Mutual Inductor instance.
1131 // Scope : public
1132 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1133 // Creation Date : 03/21/2005
1134 //-----------------------------------------------------------------------------
1136 {
1137 #ifdef Xyce_DEBUG_DEVICE
1138  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
1139  {
1140  Xyce::dout() << "Instance::loadDAEdQdx--------------------------" << std::endl
1141  << "\tname = " << getName() << std::endl;
1142  }
1143 #endif
1144 
1145  // During the DCOP, the dQdx*pdt matrix is summed into the Jacobian,
1146  // even though the Q*pdt vector is not summed into the residual.
1147  //if (!getSolverState().dcopFlag)
1148  {
1149  N_LAS_Matrix & dQdx = *(extData.dQdxMatrixPtr);
1150 
1151  // loop over each inductor and load it's Q vector components
1152  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1153  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1154  int i=0;
1155  while( currentInductor != endInductor )
1156  {
1157  for( int j=0; j<numInductors; ++j )
1158  {
1159  dQdx[((*currentInductor)->li_Branch)]
1160  [(*currentInductor)->inductorCurrentOffsets[j]] += mutualCouplingCoef[i][j] * LO[i][j];
1161  }
1162  // finally do all the dependent variable terms
1163  int numdepterms=(*currentInductor)->depVarPairs.size();
1164  for (int j=0 ; j<numdepterms; ++j)
1165  {
1166  int coefficientNumber=(*currentInductor)->depVarPairs[j].first;
1167  int depVarNumber=(*currentInductor)->depVarPairs[j].second;
1168  int otherInductor;
1169 
1170  // indexPairs[coefficient] gives the two inductors coupled by that
1171  // coefficient. One of them is currentInductor because we saved that
1172  // coefficient number in our depVarPairs, so the other is the one
1173  // we need to know
1174 
1175  if (i==indexPairs[coefficientNumber].first)
1176  otherInductor=indexPairs[coefficientNumber].second;
1177  else
1178  otherInductor=indexPairs[coefficientNumber].first;
1179 
1180  dQdx[((*currentInductor)->li_Branch)][(*currentInductor)->ABraEquDepVarOffsets[j]] +=
1181  couplingCoefficientVarDerivs[coefficientNumber][depVarNumber]
1182  *LO[i][otherInductor]*inductorCurrents[otherInductor];
1183  }
1184 
1185  ++i;
1186  ++currentInductor;
1187  }
1188  }
1189 
1190  return true;
1191 }
1192 
1193 //-----------------------------------------------------------------------------
1194 // Function : Instance::loadDAEdFdx ()
1195 //
1196 // Purpose : Loads the F-vector contributions for a single
1197 // Mutual Inductor instance.
1198 //
1199 // Special Notes :
1200 //
1201 // Scope : public
1202 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1203 // Creation Date : 03/21/2005
1204 //-----------------------------------------------------------------------------
1206 {
1207  N_LAS_Matrix & dFdx = *(extData.dFdxMatrixPtr);
1208 
1209  // loop over each inductor and load it's dFdx components
1210  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1211  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1212  while( currentInductor != endInductor )
1213  {
1214  dFdx[((*currentInductor)->li_Pos)] [((*currentInductor)->APosEquBraVarOffset)] += scalingRHS;
1215  dFdx[((*currentInductor)->li_Neg)] [((*currentInductor)->ANegEquBraVarOffset)] += -scalingRHS;
1216  dFdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquPosNodeOffset)] += -1.0;
1217  dFdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquNegNodeOffset)] += 1.0;
1218 
1219  ++currentInductor;
1220  }
1221 
1222  return true;
1223 }
1224 
1225 //-----------------------------------------------------------------------------
1226 // Function : Instance::setIC
1227 // Purpose :
1228 // Special Notes :
1229 // Scope : public
1230 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1231 // Creation Date : 03/21/2005
1232 //-----------------------------------------------------------------------------
1234 {
1235  int i_bra_sol;
1236  int i_f_state;
1237 
1238  bool bsuccess = true;
1239  return bsuccess;
1240 }
1241 
1242 //-----------------------------------------------------------------------------
1243 // Function : Instance::varTypes
1244 // Purpose :
1245 // Special Notes :
1246 // Scope : public
1247 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1248 // Creation Date : 03/21/2005
1249 //-----------------------------------------------------------------------------
1250 void Instance::varTypes( std::vector<char> & varTypeVec )
1251 {
1252  varTypeVec.resize(numInductors);
1253  for(int i=0; i<numInductors; i++)
1254  {
1255  varTypeVec[i] = 'I';
1256  }
1257 }
1258 
1259 //-----------------------------------------------------------------------------
1260 // Function : Model::processParams
1261 // Purpose :
1262 // Special Notes :
1263 // Scope : public
1264 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1265 // Creation Date : 03/21/2005
1266 //-----------------------------------------------------------------------------
1268 {
1269  return true;
1270 }
1271 
1272 //----------------------------------------------------------------------------
1273 // Function : Model::processInstanceParams
1274 // Purpose :
1275 // Special Notes :
1276 // Scope : public
1277 // Creator : Dave Shirely, PSSI
1278 // Creation Date : 03/23/06
1279 //----------------------------------------------------------------------------
1281 {
1282  std::vector<Instance*>::iterator iter;
1283  std::vector<Instance*>::iterator first = instanceContainer.begin();
1284  std::vector<Instance*>::iterator last = instanceContainer.end();
1285 
1286  for (iter=first; iter!=last; ++iter)
1287  {
1288  (*iter)->processParams();
1289  }
1290 
1291  return true;
1292 }
1293 
1294 //-----------------------------------------------------------------------------
1295 // Function : Model::Model
1296 // Purpose : block constructor
1297 // Special Notes :
1298 // Scope : public
1299 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1300 // Creation Date : 03/21/2005
1301 //-----------------------------------------------------------------------------
1303  const Configuration & configuration,
1304  const ModelBlock & MB,
1305  const FactoryBlock & factory_block)
1306  : DeviceModel(MB, configuration.getModelParameters(), factory_block),
1307  tempCoeff1(0.0),
1308  tempCoeff2(0.0),
1309  tnom(getDeviceOptions().tnom)
1310 {
1311 
1312  // Set params to constant default values:
1313  setDefaultParams ();
1314 
1315  // Set params according to .model line and constant defaults from metadata:
1316  setModParams (MB.params);
1317 
1318  // Set any non-constant parameter defaults:
1319  if (!given("TNOM"))
1321 
1322  // Calculate any parameters specified as expressions:
1324 
1325  // calculate dependent (ie computed) params and check for errors:
1326  processParams ();
1327 }
1328 
1329 //-----------------------------------------------------------------------------
1330 // Function : Model::~Model
1331 // Purpose : destructor
1332 // Special Notes :
1333 // Scope : public
1334 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1335 // Creation Date : 03/21/2005
1336 //-----------------------------------------------------------------------------
1338 {
1339  std::vector<Instance*>::iterator iter;
1340  std::vector<Instance*>::iterator first = instanceContainer.begin();
1341  std::vector<Instance*>::iterator last = instanceContainer.end();
1342 
1343  for (iter=first; iter!=last; ++iter)
1344  {
1345  delete (*iter);
1346  }
1347 
1348 }
1349 
1350 // additional Declarations
1351 
1352 //-----------------------------------------------------------------------------
1353 // Function : Model::printOutInstances
1354 // Purpose : debugging tool.
1355 // Special Notes :
1356 // Scope : public
1357 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1358 // Creation Date : 03/21/2005
1359 //-----------------------------------------------------------------------------
1360 std::ostream &Model::printOutInstances(std::ostream &os) const
1361 {
1362  std::vector<Instance*>::const_iterator iter;
1363  std::vector<Instance*>::const_iterator first = instanceContainer.begin();
1364  std::vector<Instance*>::const_iterator last = instanceContainer.end();
1365 
1366  int i, isize;
1367  isize = instanceContainer.size();
1368 
1369  os << std::endl;
1370  os << "Number of MutIndLin instances: " << isize << std::endl;
1371  os << " name=\t\tmodelName\tParameters" << std::endl;
1372  for (i=0, iter=first; iter!=last; ++iter, ++i)
1373  {
1374  os << " " << i << ": " << (*iter)->getName() << "\t";
1375  os << getName();
1376  os << std::endl;
1377  }
1378 
1379  os << std::endl;
1380 
1381  return os;
1382 }
1383 
1384 //-----------------------------------------------------------------------------
1385 // Function : Model::forEachInstance
1386 // Purpose :
1387 // Special Notes :
1388 // Scope : public
1389 // Creator : David Baur
1390 // Creation Date : 2/4/2014
1391 //-----------------------------------------------------------------------------
1392 /// Apply a device instance "op" to all instances associated with this
1393 /// model
1394 ///
1395 /// @param[in] op Operator to apply to all instances.
1396 ///
1397 ///
1398 void Model::forEachInstance(DeviceInstanceOp &op) const /* override */
1399 {
1400  for (std::vector<Instance *>::const_iterator it = instanceContainer.begin(); it != instanceContainer.end(); ++it)
1401  op(*it);
1402 }
1403 
1404 
1405 // MutIndLin Master functions:
1406 
1407 //-----------------------------------------------------------------------------
1408 // Function : Master::updateState
1409 // Purpose :
1410 // Special Notes :
1411 // Scope : public
1412 // Creator : Eric Keiter, SNL
1413 // Creation Date : 12/12/08
1414 //-----------------------------------------------------------------------------
1415 bool Master::updateState (double * solVec, double * staVec, double * stoVec)
1416 {
1417  for (InstanceVector::const_iterator it = getInstanceBegin(); it != getInstanceEnd(); ++it)
1418  {
1419  Instance & inst = *(*it);
1420 
1421  // Evaluate the derivatives of all (dependent) coupling coefficients w.r.t
1422  // their variables. We need these for both new and old DAE.
1423  int ncoupcoef=inst.couplingCoefficient.size();
1424  for (int i=0; i<ncoupcoef; ++i)
1425  {
1426  if (inst.expPtrs[i])
1427  {
1428  double junk;
1429  inst.expPtrs[i]->evaluate( junk, inst.couplingCoefficientVarDerivs[i]);
1430  }
1431  }
1432 
1433  // get the currents in each inductor
1434  std::vector< InductorInstanceData* >::iterator
1435  currentInductor = inst.instanceData.begin();
1436  std::vector< InductorInstanceData* >::iterator
1437  endInductor = inst.instanceData.end();
1438  {
1439  int i = 0;
1440  while( currentInductor != endInductor )
1441  {
1442  if( (getSolverState().dcopFlag) && ((*currentInductor)->ICGiven) )
1443  {
1444  inst.inductorCurrents[ i ] = (*currentInductor)->IC;
1445  }
1446  else
1447  {
1448  inst.inductorCurrents[ i ] = solVec[ (*currentInductor)->li_Branch ];
1449  }
1450  ++i;
1451  ++currentInductor;
1452  }
1453  }
1454  }
1455 
1456  return true;
1457 }
1458 
1459 //-----------------------------------------------------------------------------
1460 // Function : Master::loadDAEVectors
1461 // Purpose :
1462 // Special Notes :
1463 // Scope : public
1464 // Creator : Eric Keiter, SNL
1465 // Creation Date : 12/12/08
1466 //-----------------------------------------------------------------------------
1467 bool Master::loadDAEVectors (double * solVec, double * fVec, double *qVec, double * bVec, double * storeLeadF, double * storeLeadQ)
1468 {
1469  for (InstanceVector::const_iterator it = getInstanceBegin(); it != getInstanceEnd(); ++it)
1470  {
1471  Instance & inst = *(*it);
1472 
1473  // F-vector:
1474  // loop over each inductor and load it's F vector components
1475  std::vector< InductorInstanceData* >::iterator currentInductor = inst.instanceData.begin();
1476  std::vector< InductorInstanceData* >::iterator endInductor = inst.instanceData.end();
1477 
1478  while( currentInductor != endInductor )
1479  {
1480  double current = solVec[(*currentInductor)->li_Branch];
1481  double vNodePos = solVec[(*currentInductor)->li_Pos];
1482  double vNodeNeg = solVec[(*currentInductor)->li_Neg];
1483 
1484  fVec[((*currentInductor)->li_Pos)] += inst.scalingRHS * current;
1485  fVec[((*currentInductor)->li_Neg)] += -inst.scalingRHS * current;
1486  fVec[((*currentInductor)->li_Branch)] += -(vNodePos - vNodeNeg);
1487 
1488  ++currentInductor;
1489  }
1490 
1491  // Q-vector:
1492  // calculate the following product
1493  // I = column vector of currents
1494  // L = row vector of inductances
1495  // LO = matrix = mutualCup * sqrt( L' * L )
1496  // LOI = column vector = mutualCup * sqrt( L' * L ) * I
1497  // LOI[1] = mutualCup * sqrt(L[1]*L[1])*I[1]) +
1498  // mutualCup * sqrt(L[1]*L[2])*I[2]) + ...
1499  // mutualCup * sqrt(L[1]*L[n])*I[n])
1500 
1501  for( int i = 0; i < inst.numInductors; ++i )
1502  {
1503  inst.LOI[ i ] = 0.0;
1504  for( int j = 0; j < inst.numInductors; ++j )
1505  {
1506  inst.LOI[i] += inst.mutualCouplingCoef[i][j] * inst.LO[i][j] * inst.inductorCurrents[j];
1507  }
1508  }
1509 
1510  // loop over each inductor and load it's Q vector components
1511  currentInductor = inst.instanceData.begin();
1512  endInductor = inst.instanceData.end();
1513  int li=0;
1514  while( currentInductor != endInductor )
1515  {
1516  qVec[(*currentInductor)->li_Branch] += inst.LOI[ li ];
1517  ++li;
1518  ++currentInductor;
1519  }
1520  }
1521  return true;
1522 }
1523 
1524 #ifndef Xyce_NONPOINTER_MATRIX_LOAD
1525 //-----------------------------------------------------------------------------
1526 // Function : Master::loadDAEMatrices
1527 // Purpose :
1528 // Special Notes :
1529 // Scope : public
1530 // Creator : Eric Keiter, SNL
1531 // Creation Date : 12/12/08
1532 //-----------------------------------------------------------------------------
1533 bool Master::loadDAEMatrices (N_LAS_Matrix & dFdx, N_LAS_Matrix & dQdx)
1534 {
1535  for (InstanceVector::const_iterator it = getInstanceBegin(); it != getInstanceEnd(); ++it)
1536  {
1537  Instance & inst = *(*it);
1538  // loop over each inductor and load it's dFdx components
1539  std::vector< InductorInstanceData* >::iterator currentInductor = inst.instanceData.begin();
1540  std::vector< InductorInstanceData* >::iterator endInductor = inst.instanceData.end();
1541  while( currentInductor != endInductor )
1542  {
1543  *((*currentInductor)->f_PosEquBraVarPtr) += inst.scalingRHS;
1544  *((*currentInductor)->f_NegEquBraVarPtr) += -inst.scalingRHS;
1545  *((*currentInductor)->f_BraEquPosNodePtr) += -1.0;
1546  *((*currentInductor)->f_BraEquNegNodePtr) += 1.0;
1547 
1548  ++currentInductor;
1549  }
1550 
1551  // During the DCOP, the dQdx*pdt matrix is summed into the Jacobian,
1552  // even though the Q*pdt vector is not summed into the residual.
1553  //if (!getSolverState().dcopFlag)
1554  {
1555  // loop over each inductor and load it's Q vector components
1556  currentInductor = inst.instanceData.begin();
1557  endInductor = inst.instanceData.end();
1558  int li=0;
1559  while( currentInductor != endInductor )
1560  {
1561  for( int j=0; j<inst.numInductors; ++j )
1562  {
1563  *((*currentInductor)->q_inductorCurrentPtrs[j]) += inst.mutualCouplingCoef[li][j] * inst.LO[li][j];
1564  }
1565  // finally do all the dependent variable terms
1566  int numdepterms=(*currentInductor)->depVarPairs.size();
1567  for (int j=0 ; j<numdepterms; ++j)
1568  {
1569  int coefficientNumber=(*currentInductor)->depVarPairs[j].first;
1570  int depVarNumber=(*currentInductor)->depVarPairs[j].second;
1571  int otherInductor;
1572 
1573  // indexPairs[coefficient] gives the two inductors coupled by that
1574  // coefficient. One of them is currentInductor because we saved that
1575  // coefficient number in our depVarPairs, so the other is the one
1576  // we need to know
1577 
1578  if (li==inst.indexPairs[coefficientNumber].first)
1579  otherInductor=inst.indexPairs[coefficientNumber].second;
1580  else
1581  otherInductor=inst.indexPairs[coefficientNumber].second;
1582 
1583  *((*currentInductor)->q_BraEquDepVarPtrs[j]) +=
1584  inst.couplingCoefficientVarDerivs[coefficientNumber][depVarNumber]
1585  *inst.LO[li][otherInductor]*inst.inductorCurrents[otherInductor];
1586  }
1587 
1588  ++li;
1589  ++currentInductor;
1590  }
1591  }
1592  }
1593 
1594  return true;
1595 }
1596 
1597 #else
1598 //-----------------------------------------------------------------------------
1599 // Function : Master::loadDAEMatrices
1600 // Purpose :
1601 // Special Notes :
1602 // Scope : public
1603 // Creator : Eric Keiter, SNL
1604 // Creation Date : 12/12/08
1605 //-----------------------------------------------------------------------------
1606 bool Master::loadDAEMatrices (N_LAS_Matrix & dFdx, N_LAS_Matrix & dQdx)
1607 {
1608  int sizeInstances = instanceContainer_.size();
1609  for (int i=0; i<sizeInstances; ++i)
1610  {
1611  Instance & inst = *(instanceContainer_.at(i));
1612  // loop over each inductor and load it's dFdx components
1613  std::vector< InductorInstanceData* >::iterator currentInductor = inst.instanceData.begin();
1614  std::vector< InductorInstanceData* >::iterator endInductor = inst.instanceData.end();
1615  while( currentInductor != endInductor )
1616  {
1617  dFdx[((*currentInductor)->li_Pos)] [((*currentInductor)->APosEquBraVarOffset)] += inst.scalingRHS;
1618  dFdx[((*currentInductor)->li_Neg)] [((*currentInductor)->ANegEquBraVarOffset)] += -inst.scalingRHS;
1619  dFdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquPosNodeOffset)] += -1.0;
1620  dFdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquNegNodeOffset)] += 1.0;
1621 
1622  ++currentInductor;
1623  }
1624 
1625  // During the DCOP, the dQdx*pdt matrix is summed into the Jacobian,
1626  // even though the Q*pdt vector is not summed into the residual.
1627  //if (!getSolverState().dcopFlag)
1628  {
1629  // loop over each inductor and load it's Q vector components
1630  currentInductor = inst.instanceData.begin();
1631  endInductor = inst.instanceData.end();
1632  int i=0;
1633  while( currentInductor != endInductor )
1634  {
1635  for( int j=0; j<inst.numInductors; ++j )
1636  {
1637  dQdx[((*currentInductor)->li_Branch)]
1638  [(*currentInductor)->inductorCurrentOffsets[j]] += inst.mutualCouplingCoef[i][j] * inst.LO[i][j];
1639  }
1640  // finally do all the dependent variable terms
1641  int numdepterms=(*currentInductor)->depVarPairs.size();
1642  for (int j=0 ; j<numdepterms; ++j)
1643  {
1644  int coefficientNumber=(*currentInductor)->depVarPairs[j].first;
1645  int depVarNumber=(*currentInductor)->depVarPairs[j].second;
1646  int otherInductor;
1647 
1648  // indexPairs[coefficient] gives the two inductors coupled by that
1649  // coefficient. One of them is currentInductor because we saved that
1650  // coefficient number in our depVarPairs, so the other is the one
1651  // we need to know
1652 
1653  if (i==inst.indexPairs[coefficientNumber].first)
1654  otherInductor=inst.indexPairs[coefficientNumber].second;
1655  else
1656  otherInductor=inst.indexPairs[coefficientNumber].second;
1657 
1658  dQdx[((*currentInductor)->li_Branch)][(*currentInductor)->ABraEquDepVarOffsets[j]] +=
1659  inst.couplingCoefficientVarDerivs[coefficientNumber][depVarNumber]
1660  *inst.LO[i][otherInductor]*inst.inductorCurrents[otherInductor];
1661  }
1662 
1663  ++i;
1664  ++currentInductor;
1665  }
1666  }
1667  }
1668  return true;
1669 }
1670 #endif
1671 
1672 Device *Traits::factory(const Configuration &configuration, const FactoryBlock &factory_block)
1673 {
1674 
1675  return new Master(configuration, factory_block, factory_block.solverState_, factory_block.deviceOptions_);
1676 }
1677 
1679 {
1681  .registerDevice("mil", 1)
1682  .registerModelType("mil", 1);
1683 }
1684 
1685 } // namespace MutIndLin
1686 } // namespace Device
1687 } // namespace Xyce