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.124.2.2 $
40 //
41 // Revision Date : $Date: 2014/03/06 23:33:43 $
42 //
43 // Current Owner : $Author: tvrusso $
44 //-------------------------------------------------------------------------
45 #include <Xyce_config.h>
46 
47 // ---------- Standard Includes ----------
48 #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() ; 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<sDepend>::iterator d;
383  std::vector<sDepend>::iterator begin=dependentParams.begin();
384  std::vector<sDepend>::iterator end=dependentParams.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();
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  std::string tmpstr;
564  for( int i=0; i<instanceData.size() ; ++i)
565  {
566  tmpstr = getName() + "_" + instanceData[i]->name + "_branch";
567  spiceInternalName (tmpstr);
568  intNameMap[ intLIDVec[i] ] = tmpstr;
569  }
570  }
571 
572  return intNameMap;
573 }
574 
575 //-----------------------------------------------------------------------------
576 // Function : Instance::registerStateLIDs
577 // Purpose :
578 // Special Notes :
579 // Scope : public
580 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
581 // Creation Date : 03/21/2005
582 //-----------------------------------------------------------------------------
583 void Instance::registerStateLIDs( const std::vector<int> & staLIDVecRef )
584 {
585  AssertLIDs(staLIDVecRef.size() == numStateVars);
586 
587  // copy over the global ID lists.
588  staLIDVec = staLIDVecRef;
589 }
590 
591 //-----------------------------------------------------------------------------
592 // Function : BsrcInstance::getDepSolnVars
593 // Purpose :
594 // Special Notes :
595 // Scope : public
596 // Creator : Rob Hoekstra, SNL, Parallel Computational Sciences
597 // Creation Date : 06/06/01
598 //-----------------------------------------------------------------------------
599 const std::vector<std::string> & Instance::getDepSolnVars()
600 {
602 }
603 
604 //-----------------------------------------------------------------------------
605 // Function : Instance::jacobianStamp
606 // Purpose :
607 // Special Notes :
608 // Scope : public
609 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
610 // Creation Date : 03/21/2005
611 //-----------------------------------------------------------------------------
612 const std::vector< std::vector<int> > & Instance::jacobianStamp() const
613 {
614  return jacStamp;
615 }
616 
617 //-----------------------------------------------------------------------------
618 // Function : Instance::registerJacLIDs
619 // Purpose :
620 // Special Notes :
621 // Scope : public
622 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
623 // Creation Date : 03/21/2005
624 //-----------------------------------------------------------------------------
625 void Instance::registerJacLIDs( const std::vector< std::vector<int> > & jacLIDVec )
626 {
627  DeviceInstance::registerJacLIDs( jacLIDVec );
628 #ifdef Xyce_DEBUG_DEVICE
629  if (getDeviceOptions().debugLevel > 0)
630  {
631  Xyce::dout() << "Instance::registerJacLIDs ----------------------------" << std::endl;
632  }
633 #endif
634  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
635  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
636  int i = 0;
637  while( currentInductor != endInductor )
638  {
639  (*currentInductor)->APosEquBraVarOffset = jacLIDVec[ 2*i ][ 0 ];
640  (*currentInductor)->ANegEquBraVarOffset = jacLIDVec[ 2*i + 1 ][ 0 ];
641  (*currentInductor)->ABraEquPosNodeOffset = jacLIDVec[ 2*numInductors + i ][ 0 ];
642  (*currentInductor)->ABraEquNegNodeOffset = jacLIDVec[ 2*numInductors + i ][ 1 ];
643  for( int j=0; j<numInductors; ++j )
644  {
645  if( i == j )
646  {
647  (*currentInductor)->ABraEquBraVarOffset = jacLIDVec[ 2*numInductors + i ][ j + 2 ];
648  }
649  (*currentInductor)->inductorCurrentOffsets[ j ] = jacLIDVec[ 2*numInductors + i ][ j + 2 ];
650  }
651  // Now do the parts that are for the dependent variables
652  int numdepvars=(*currentInductor)->depVarPairs.size();
653  (*currentInductor)->ABraEquDepVarOffsets.resize(numdepvars);
654  for (int j = 0; j < numdepvars; ++j)
655  {
656  (*currentInductor)->ABraEquDepVarOffsets[j] =
657  jacLIDVec[2*numInductors+i][numInductors+2+j];
658  }
659  ++currentInductor;
660  ++i;
661  }
662 
663 #ifdef Xyce_DEBUG_DEVICE
664  if (getDeviceOptions().debugLevel > 0)
665  {
666  Xyce::dout() << "Instance::registerJacLIDs--------------------------" << std::endl;
667  currentInductor = instanceData.begin();
668  i=0;
669  while( currentInductor != endInductor )
670  {
671  Xyce::dout() << "Inductor [ " << i << " ] " << (*currentInductor)->name
672  << " APosEquBraVarOffset = " << (*currentInductor)->APosEquBraVarOffset
673  << " ANegEquBraVarOffset = " << (*currentInductor)->ANegEquBraVarOffset
674  << " ABraEquPosNodeOffset = " << (*currentInductor)->ABraEquPosNodeOffset
675  << " ABraEquNegNodeOffset = " << (*currentInductor)->ABraEquNegNodeOffset
676  << " ABraEquBraVarOffset = " << (*currentInductor)->ABraEquBraVarOffset << std::endl;
677  Xyce::dout() << "\tInductor branch offsets = { ";
678  for( int j=0; j<numInductors ; ++j )
679  {
680  Xyce::dout() << (*currentInductor)->inductorCurrentOffsets[ j ] << ", ";
681  }
682  Xyce::dout() << "} " << std::endl;
683  ++i;
684  ++currentInductor;
685  }
686  }
687 #endif
688 }
689 
690 //-----------------------------------------------------------------------------
691 // Function : Instance::setupPointers
692 // Purpose :
693 // Special Notes :
694 // Scope : public
695 // Creator : Eric Keiter, SNL
696 // Creation Date : 12/12/08
697 //-----------------------------------------------------------------------------
699 {
700 #ifndef Xyce_NONPOINTER_MATRIX_LOAD
701  N_LAS_Matrix & dFdx = *(extData.dFdxMatrixPtr);
702  N_LAS_Matrix & dQdx = *(extData.dQdxMatrixPtr);
703 
704  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
705  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
706  int i = 0;
707  while( currentInductor != endInductor )
708  {
709  // dFdx pointers:
710  (*currentInductor)->f_PosEquBraVarPtr =
711  &(dFdx[((*currentInductor)->li_Pos)] [((*currentInductor)->APosEquBraVarOffset)] );
712 
713  (*currentInductor)->f_NegEquBraVarPtr =
714  &(dFdx[((*currentInductor)->li_Neg)] [((*currentInductor)->ANegEquBraVarOffset)] );
715 
716  (*currentInductor)->f_BraEquPosNodePtr =
717  &(dFdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquPosNodeOffset)] );
718 
719  (*currentInductor)->f_BraEquNegNodePtr =
720  &(dFdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquNegNodeOffset)] );
721 
722  for( int j=0; j<numInductors; ++j )
723  {
724  if( i == j )
725  {
726  //(*currentInductor)->ABraEquBraVarOffset =
727  // Is this only used for JMat? (ie old DAE)?
728 
729  }
730  (*currentInductor)->f_inductorCurrentPtrs[ j ] =
731  &(dFdx[((*currentInductor)->li_Branch)][(*currentInductor)->inductorCurrentOffsets[j]] );
732  }
733  // Now do the parts that are for the dependent variables
734  int numdepvars=(*currentInductor)->depVarPairs.size();
735  (*currentInductor)->f_BraEquDepVarPtrs.resize(numdepvars);
736 
737  for (int j = 0; j < numdepvars; ++j)
738  {
739  (*currentInductor)->f_BraEquDepVarPtrs[j] =
740  &(dFdx[((*currentInductor)->li_Branch)][(*currentInductor)->ABraEquDepVarOffsets[j]] );
741 
742  }
743 
744  // dQdx pointers:
745  (*currentInductor)->q_PosEquBraVarPtr =
746  &(dQdx[((*currentInductor)->li_Pos)] [((*currentInductor)->APosEquBraVarOffset)] );
747 
748  (*currentInductor)->q_NegEquBraVarPtr =
749  &(dQdx[((*currentInductor)->li_Neg)] [((*currentInductor)->ANegEquBraVarOffset)] );
750 
751  (*currentInductor)->q_BraEquPosNodePtr =
752  &(dQdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquPosNodeOffset)] );
753 
754  (*currentInductor)->q_BraEquNegNodePtr =
755  &(dQdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquNegNodeOffset)] );
756 
757  for( int j=0; j<numInductors; ++j )
758  {
759  if( i == j )
760  {
761  //(*currentInductor)->ABraEquBraVarOffset =
762  // Is this only used for JMat? (ie old DAE)?
763 
764  }
765  (*currentInductor)->q_inductorCurrentPtrs[ j ] =
766  &(dQdx[((*currentInductor)->li_Branch)][(*currentInductor)->inductorCurrentOffsets[j]] );
767  }
768  // Now do the parts that are for the dependent variables
769  numdepvars=(*currentInductor)->depVarPairs.size();
770  (*currentInductor)->q_BraEquDepVarPtrs.resize(numdepvars);
771 
772  for (int j = 0; j < numdepvars; ++j)
773  {
774  (*currentInductor)->q_BraEquDepVarPtrs[j] =
775  &(dQdx[((*currentInductor)->li_Branch)][(*currentInductor)->ABraEquDepVarOffsets[j]] );
776 
777  }
778 
779  ++currentInductor;
780  ++i;
781  }
782 
783 
784 
785 #endif
786 }
787 
788 //-----------------------------------------------------------------------------
789 // Function : Instance::processParams
790 // Purpose :
791 // Special Notes :
792 // Scope : public
793 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
794 // Creation Date : 03/21/2005
795 //-----------------------------------------------------------------------------
797 {
798  int i, j;
799 
800  // take the vector couplingInductor in pairs to figure out
801  // the coupling value between inductors then insert it into
802  // the couplingCoefficient matrix.
803  j = indexPairs.size();
804  for( i=0; i<j ; ++i )
805  {
806  // note that couplingCoefficient can be of length 1 if all coefficients were the
807  // same. Catch that case here.
808  if( i < couplingCoefficient.size() )
809  {
810  mutualCouplingCoef[ indexPairs[i].first ][ indexPairs[i].second ] =
811  mutualCouplingCoef[ indexPairs[i].second ][ indexPairs[i].first ] = couplingCoefficient[i];
812  }
813  else
814  {
815  mutualCouplingCoef[ indexPairs[i].first ][ indexPairs[i].second ] =
816  mutualCouplingCoef[ indexPairs[i].second ][ indexPairs[i].first ] = couplingCoefficient[0];
817  }
818  }
819 #ifdef Xyce_DEBUG_DEVICE
821  {
822  Xyce::dout() << "processParams: mutualCouplingCoef: " << std::endl;
823  for( int i=0; i<numInductors; ++i )
824  {
825  for( int j=0; j<numInductors; ++j )
826  {
827  Xyce::dout() << mutualCouplingCoef[i][j] << " ";
828  }
829  Xyce::dout() << std::endl;
830  }
831  }
832 #endif
833 
834  // set the temperature related stuff.
836 
837  return true;
838 }
839 
840 //-----------------------------------------------------------------------------
841 // Function : Instance::updateTemperature
842 // Purpose :
843 // Special Notes :
844 // Scope : public
845 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
846 // Creation Date : 03/21/2005
847 //-----------------------------------------------------------------------------
848 bool Instance::updateTemperature ( const double & temp)
849 {
850  bool bsuccess = true;
851 
852  // current temp difference from reference temp.
853  double difference = temp - model_.tnom;
854 
855  std::vector< InductorInstanceData* >::iterator currentData = instanceData.begin();
856  while( currentData != instanceData.end() )
857  {
858  double factor = 1.0 + (model_.tempCoeff1)*difference +
859  (model_.tempCoeff2)*difference*difference;
860  (*currentData)->L = ((*currentData)->baseL)*factor;
861  ++currentData;
862  }
863 
864  // now that the inductances have changed we need to update the matrix.
866 
867  return bsuccess;
868 }
869 
870 //-----------------------------------------------------------------------------
871 // Function : Instance::updateIntermediateVars
872 // Purpose :
873 // Special Notes :
874 // Scope : public
875 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
876 // Creation Date : 03/21/2005
877 //-----------------------------------------------------------------------------
879 {
880  bool bsuccess = true;
881 
882  // This method is never called anymore
883 
884  return bsuccess;
885 }
886 
887 //-----------------------------------------------------------------------------
888 // Function : Instance::updateInductanceMatrix()
889 // Purpose : A matrix of inductances is used often enough that it
890 // calculated and stored as a member variable here
891 // If and inductance ever changes say from updating
892 // the temperature or a parameter udpate, then this
893 // routine must be called again.
894 // Special Notes :
895 // Scope : public
896 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
897 // Creation Date : 03/21/2005
898 //-----------------------------------------------------------------------------
900 {
901  std::vector< InductorInstanceData* >::iterator
902  currentInductor = instanceData.begin();
903  std::vector< InductorInstanceData* >::iterator
904  endInductor = instanceData.end();
905 
906  // collec the inductances
907  int i=0;
908  while( currentInductor != endInductor )
909  {
910  inductanceVals[ i ] = ((*currentInductor)->L);
911  ++i;
912  ++currentInductor;
913  }
914 
915  // compute the inductance matrix
916  for( i=0; i<numInductors; ++i)
917  {
918  for( int j=0; j<numInductors; ++j)
919  {
920  LO[i][j] = sqrt( inductanceVals[i]*inductanceVals[j] );
921  }
922  }
923 
924 }
925 //-----------------------------------------------------------------------------
926 // Function : Instance::updatePrimaryState
927 // Purpose :
928 // Special Notes :
929 // Scope : public
930 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
931 // Creation Date : 03/21/2005
932 //-----------------------------------------------------------------------------
934 {
935  // nothing happens in updateIntermediateVars() anymore
936  // bsuccess = updateIntermediateVars ();
937 
938 #ifdef Xyce_DEBUG_DEVICE
939  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
940  {
941  Xyce::dout() << "Instance::updatePrimaryState------------------" << std::endl
942  << "\tname = " << getName() << std::endl;
943  }
944 #endif
945 
946  double * solVec = extData.nextSolVectorRawPtr;
947 
948  // Evaluate the derivatives of all (dependent) coupling coefficients w.r.t
949  // their variables. We need these for both new and old DAE.
950  int ncoupcoef=couplingCoefficient.size();
951  for (int i=0; i<ncoupcoef; ++i)
952  {
953  if (expPtrs[i])
954  {
955  double junk;
956  expPtrs[i]->evaluate( junk, couplingCoefficientVarDerivs[i]);
957  }
958  }
959 
960  // get the currents in each inductor
961  std::vector< InductorInstanceData* >::iterator
962  currentInductor = instanceData.begin();
963  std::vector< InductorInstanceData* >::iterator
964  endInductor = instanceData.end();
965  int i = 0;
966  while( currentInductor != endInductor )
967  {
968  if( (getSolverState().dcopFlag) && ((*currentInductor)->ICGiven) )
969  {
970  inductorCurrents[ i ] = (*currentInductor)->IC;
971  }
972  else
973  {
974  inductorCurrents[ i ] = solVec[ (*currentInductor)->li_Branch ];
975  }
976  ++i;
977  ++currentInductor;
978  }
979 
980  return true;
981 }
982 
983 //-----------------------------------------------------------------------------
984 // Function : Instance::loadDeviceMask
985 //
986 // Purpose : Loads the zero elements of the device mask
987 //
988 // Special Notes : elements of the error vector associated with zero
989 // elements of the mask will not be included in weighted
990 // norms by the time integrator.
991 //
992 // Scope : public
993 // Creator : Richard Schiek, SNL, Electrical and Microsystems Modeling
994 // Creation Date : 02/06/07
995 //-----------------------------------------------------------------------------
997 {
998  bool returnVal=false;
999 #ifndef Xyce_NO_MUTINDLIN_MASK
1000  N_LAS_Vector * maskVectorPtr = extData.deviceMaskVectorPtr;
1001 
1002  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1003  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1004  while( currentInductor != endInductor )
1005  {
1006  (*maskVectorPtr)[((*currentInductor)->li_Branch)] = 0.0;
1007  ++currentInductor;
1008  }
1009  returnVal = true;
1010 #endif
1011  return (returnVal);
1012 }
1013 
1014 //-----------------------------------------------------------------------------
1015 // Function : Instance::loadDAEQVector
1016 //
1017 // Purpose : Loads the Q-vector contributions for a single
1018 // Mutual Inductor instance.
1019 //
1020 // Special Notes : The "Q" vector is part of a standard DAE formalism in
1021 // which the system of equations is represented as:
1022 //
1023 // f(x) = dQ(x)/dt + F(x) - B(t) = 0
1024 //
1025 // Scope : public
1026 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1027 // Creation Date : 03/21/2005
1028 //-----------------------------------------------------------------------------
1030 {
1031 #ifdef Xyce_DEBUG_DEVICE
1032  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
1033  {
1034  Xyce::dout() << "Instance::loadDAEQVector---------------------------" << std::endl
1035  << "\tname = " << getName() << std::endl;
1036  }
1037 #endif
1038  double * qVec = extData.daeQVectorRawPtr;
1039 
1040  // calculate the following product
1041  // I = column vector of currents
1042  // L = row vector of inductances
1043  // LO = matrix = mutualCup * sqrt( L' * L )
1044  // LOI = column vector = mutualCup * sqrt( L' * L ) * I
1045  // LOI[1] = mutualCup * sqrt(L[1]*L[1])*I[1]) +
1046  // mutualCup * sqrt(L[1]*L[2])*I[2]) + ...
1047  // mutualCup * sqrt(L[1]*L[n])*I[n])
1048 
1049  for( int i = 0; i < numInductors; ++i )
1050  {
1051  LOI[ i ] = 0.0;
1052  for( int j = 0; j < numInductors; ++j )
1053  {
1054  LOI[i] += mutualCouplingCoef[i][j] * LO[i][j] * inductorCurrents[j];
1055  }
1056  }
1057 
1058  // loop over each inductor and load it's Q vector components
1059  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1060  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1061  int i=0;
1062  while( currentInductor != endInductor )
1063  {
1064  qVec[(*currentInductor)->li_Branch] += LOI[ i ];
1065  ++i;
1066  ++currentInductor;
1067  }
1068 
1069  return true;
1070 }
1071 
1072 //-----------------------------------------------------------------------------
1073 // Function : Instance::loadDAEFVector
1074 //
1075 // Purpose : Loads the F-vector contributions for a single
1076 // Mutual Inductor instance.
1077 //
1078 // Special Notes :
1079 //
1080 // Scope : public
1081 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1082 // Creation Date : 03/21/2005
1083 //-----------------------------------------------------------------------------
1085 {
1086 #ifdef Xyce_DEBUG_DEVICE
1087  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
1088  {
1089  Xyce::dout() << "Instance::loadDAEFVector---------------------------" << std::endl
1090  << "\tname = " << getName() << std::endl;
1091  }
1092 #endif
1093 
1094  double * fVec = extData.daeFVectorRawPtr;
1095  double * solVec = extData.nextSolVectorRawPtr;
1096 
1097  // loop over each inductor and load it's F vector components
1098  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1099  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1100  int i = 0;
1101 
1102  while( currentInductor != endInductor )
1103  {
1104  double current = solVec[(*currentInductor)->li_Branch];
1105  double vNodePos = solVec[(*currentInductor)->li_Pos];
1106  double vNodeNeg = solVec[(*currentInductor)->li_Neg];
1107  fVec[((*currentInductor)->li_Pos)] += scalingRHS * current;
1108  fVec[((*currentInductor)->li_Neg)] += -scalingRHS * current;
1109  fVec[((*currentInductor)->li_Branch)] += -(vNodePos - vNodeNeg);
1110 
1111 #ifdef Xyce_DEBUG_DEVICE
1112  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
1113  {
1114  Xyce::dout() << " Inductor = " << (*currentInductor)->name
1115  << "\tcurrent = " << current
1116  << "\tvNodePos = " << vNodePos
1117  << "\tvNodeNeg = " << vNodeNeg
1118  << std::endl;
1119  }
1120 #endif
1121 
1122  ++currentInductor;
1123  ++i;
1124  }
1125 
1126  return true;
1127 }
1128 
1129 //-----------------------------------------------------------------------------
1130 // Function : Instance::loadDAEdQdx
1131 //
1132 // Purpose : Loads the Q-vector contributions for a single
1133 // Mutual Inductor instance.
1134 // Scope : public
1135 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1136 // Creation Date : 03/21/2005
1137 //-----------------------------------------------------------------------------
1139 {
1140 #ifdef Xyce_DEBUG_DEVICE
1141  if (getDeviceOptions().debugLevel > 0 && getSolverState().debugTimeFlag)
1142  {
1143  Xyce::dout() << "Instance::loadDAEdQdx--------------------------" << std::endl
1144  << "\tname = " << getName() << std::endl;
1145  }
1146 #endif
1147 
1148  // During the DCOP, the dQdx*pdt matrix is summed into the Jacobian,
1149  // even though the Q*pdt vector is not summed into the residual.
1150  //if (!getSolverState().dcopFlag)
1151  {
1152  N_LAS_Matrix & dQdx = *(extData.dQdxMatrixPtr);
1153 
1154  // loop over each inductor and load it's Q vector components
1155  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1156  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1157  int i=0;
1158  while( currentInductor != endInductor )
1159  {
1160  for( int j=0; j<numInductors; ++j )
1161  {
1162  dQdx[((*currentInductor)->li_Branch)]
1163  [(*currentInductor)->inductorCurrentOffsets[j]] += mutualCouplingCoef[i][j] * LO[i][j];
1164  }
1165  // finally do all the dependent variable terms
1166  int numdepterms=(*currentInductor)->depVarPairs.size();
1167  for (int j=0 ; j<numdepterms; ++j)
1168  {
1169  int coefficientNumber=(*currentInductor)->depVarPairs[j].first;
1170  int depVarNumber=(*currentInductor)->depVarPairs[j].second;
1171  int otherInductor;
1172 
1173  // indexPairs[coefficient] gives the two inductors coupled by that
1174  // coefficient. One of them is currentInductor because we saved that
1175  // coefficient number in our depVarPairs, so the other is the one
1176  // we need to know
1177 
1178  if (i==indexPairs[coefficientNumber].first)
1179  otherInductor=indexPairs[coefficientNumber].second;
1180  else
1181  otherInductor=indexPairs[coefficientNumber].first;
1182 
1183  dQdx[((*currentInductor)->li_Branch)][(*currentInductor)->ABraEquDepVarOffsets[j]] +=
1184  couplingCoefficientVarDerivs[coefficientNumber][depVarNumber]
1185  *LO[i][otherInductor]*inductorCurrents[otherInductor];
1186  }
1187 
1188  ++i;
1189  ++currentInductor;
1190  }
1191  }
1192 
1193  return true;
1194 }
1195 
1196 //-----------------------------------------------------------------------------
1197 // Function : Instance::loadDAEdFdx ()
1198 //
1199 // Purpose : Loads the F-vector contributions for a single
1200 // Mutual Inductor instance.
1201 //
1202 // Special Notes :
1203 //
1204 // Scope : public
1205 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1206 // Creation Date : 03/21/2005
1207 //-----------------------------------------------------------------------------
1209 {
1210  N_LAS_Matrix & dFdx = *(extData.dFdxMatrixPtr);
1211 
1212  // loop over each inductor and load it's dFdx components
1213  std::vector< InductorInstanceData* >::iterator currentInductor = instanceData.begin();
1214  std::vector< InductorInstanceData* >::iterator endInductor = instanceData.end();
1215  while( currentInductor != endInductor )
1216  {
1217  dFdx[((*currentInductor)->li_Pos)] [((*currentInductor)->APosEquBraVarOffset)] += scalingRHS;
1218  dFdx[((*currentInductor)->li_Neg)] [((*currentInductor)->ANegEquBraVarOffset)] += -scalingRHS;
1219  dFdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquPosNodeOffset)] += -1.0;
1220  dFdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquNegNodeOffset)] += 1.0;
1221 
1222  ++currentInductor;
1223  }
1224 
1225  return true;
1226 }
1227 
1228 //-----------------------------------------------------------------------------
1229 // Function : Instance::setIC
1230 // Purpose :
1231 // Special Notes :
1232 // Scope : public
1233 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1234 // Creation Date : 03/21/2005
1235 //-----------------------------------------------------------------------------
1237 {
1238  int i_bra_sol;
1239  int i_f_state;
1240 
1241  bool bsuccess = true;
1242  return bsuccess;
1243 }
1244 
1245 //-----------------------------------------------------------------------------
1246 // Function : Instance::varTypes
1247 // Purpose :
1248 // Special Notes :
1249 // Scope : public
1250 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1251 // Creation Date : 03/21/2005
1252 //-----------------------------------------------------------------------------
1253 void Instance::varTypes( std::vector<char> & varTypeVec )
1254 {
1255  varTypeVec.resize(numInductors);
1256  for(int i=0; i<numInductors; i++)
1257  {
1258  varTypeVec[i] = 'I';
1259  }
1260 }
1261 
1262 //-----------------------------------------------------------------------------
1263 // Function : Model::processParams
1264 // Purpose :
1265 // Special Notes :
1266 // Scope : public
1267 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1268 // Creation Date : 03/21/2005
1269 //-----------------------------------------------------------------------------
1271 {
1272  return true;
1273 }
1274 
1275 //----------------------------------------------------------------------------
1276 // Function : Model::processInstanceParams
1277 // Purpose :
1278 // Special Notes :
1279 // Scope : public
1280 // Creator : Dave Shirely, PSSI
1281 // Creation Date : 03/23/06
1282 //----------------------------------------------------------------------------
1284 {
1285  std::vector<Instance*>::iterator iter;
1286  std::vector<Instance*>::iterator first = instanceContainer.begin();
1287  std::vector<Instance*>::iterator last = instanceContainer.end();
1288 
1289  for (iter=first; iter!=last; ++iter)
1290  {
1291  (*iter)->processParams();
1292  }
1293 
1294  return true;
1295 }
1296 
1297 //-----------------------------------------------------------------------------
1298 // Function : Model::Model
1299 // Purpose : block constructor
1300 // Special Notes :
1301 // Scope : public
1302 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1303 // Creation Date : 03/21/2005
1304 //-----------------------------------------------------------------------------
1306  const Configuration & configuration,
1307  const ModelBlock & MB,
1308  const FactoryBlock & factory_block)
1309  : DeviceModel(MB, configuration.getModelParameters(), factory_block),
1310  tempCoeff1(0.0),
1311  tempCoeff2(0.0),
1312  tnom(getDeviceOptions().tnom)
1313 {
1314 
1315  // Set params to constant default values:
1316  setDefaultParams ();
1317 
1318  // Set params according to .model line and constant defaults from metadata:
1319  setModParams (MB.params);
1320 
1321  // Set any non-constant parameter defaults:
1322  if (!given("TNOM"))
1324 
1325  // Calculate any parameters specified as expressions:
1327 
1328  // calculate dependent (ie computed) params and check for errors:
1329  processParams ();
1330 }
1331 
1332 //-----------------------------------------------------------------------------
1333 // Function : Model::~Model
1334 // Purpose : destructor
1335 // Special Notes :
1336 // Scope : public
1337 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1338 // Creation Date : 03/21/2005
1339 //-----------------------------------------------------------------------------
1341 {
1342  std::vector<Instance*>::iterator iter;
1343  std::vector<Instance*>::iterator first = instanceContainer.begin();
1344  std::vector<Instance*>::iterator last = instanceContainer.end();
1345 
1346  for (iter=first; iter!=last; ++iter)
1347  {
1348  delete (*iter);
1349  }
1350 
1351 }
1352 
1353 // additional Declarations
1354 
1355 //-----------------------------------------------------------------------------
1356 // Function : Model::printOutInstances
1357 // Purpose : debugging tool.
1358 // Special Notes :
1359 // Scope : public
1360 // Creator : Rich Schiek, SNL, Parallel Computational Sciences
1361 // Creation Date : 03/21/2005
1362 //-----------------------------------------------------------------------------
1363 std::ostream &Model::printOutInstances(std::ostream &os) const
1364 {
1365  std::vector<Instance*>::const_iterator iter;
1366  std::vector<Instance*>::const_iterator first = instanceContainer.begin();
1367  std::vector<Instance*>::const_iterator last = instanceContainer.end();
1368 
1369  int i, isize;
1370  isize = instanceContainer.size();
1371 
1372  os << std::endl;
1373  os << "Number of MutIndLin instances: " << isize << std::endl;
1374  os << " name=\t\tmodelName\tParameters" << std::endl;
1375  for (i=0, iter=first; iter!=last; ++iter, ++i)
1376  {
1377  os << " " << i << ": " << (*iter)->getName() << "\t";
1378  os << getName();
1379  os << std::endl;
1380  }
1381 
1382  os << std::endl;
1383 
1384  return os;
1385 }
1386 
1387 //-----------------------------------------------------------------------------
1388 // Function : Model::forEachInstance
1389 // Purpose :
1390 // Special Notes :
1391 // Scope : public
1392 // Creator : David Baur
1393 // Creation Date : 2/4/2014
1394 //-----------------------------------------------------------------------------
1395 /// Apply a device instance "op" to all instances associated with this
1396 /// model
1397 ///
1398 /// @param[in] op Operator to apply to all instances.
1399 ///
1400 ///
1401 void Model::forEachInstance(DeviceInstanceOp &op) const /* override */
1402 {
1403  for (std::vector<Instance *>::const_iterator it = instanceContainer.begin(); it != instanceContainer.end(); ++it)
1404  op(*it);
1405 }
1406 
1407 
1408 // MutIndLin Master functions:
1409 
1410 //-----------------------------------------------------------------------------
1411 // Function : Master::updateState
1412 // Purpose :
1413 // Special Notes :
1414 // Scope : public
1415 // Creator : Eric Keiter, SNL
1416 // Creation Date : 12/12/08
1417 //-----------------------------------------------------------------------------
1418 bool Master::updateState (double * solVec, double * staVec, double * stoVec)
1419 {
1420  for (InstanceVector::const_iterator it = getInstanceBegin(); it != getInstanceEnd(); ++it)
1421  {
1422  Instance & inst = *(*it);
1423 
1424  // Evaluate the derivatives of all (dependent) coupling coefficients w.r.t
1425  // their variables. We need these for both new and old DAE.
1426  int ncoupcoef=inst.couplingCoefficient.size();
1427  for (int i=0; i<ncoupcoef; ++i)
1428  {
1429  if (inst.expPtrs[i])
1430  {
1431  double junk;
1432  inst.expPtrs[i]->evaluate( junk, inst.couplingCoefficientVarDerivs[i]);
1433  }
1434  }
1435 
1436  // get the currents in each inductor
1437  std::vector< InductorInstanceData* >::iterator
1438  currentInductor = inst.instanceData.begin();
1439  std::vector< InductorInstanceData* >::iterator
1440  endInductor = inst.instanceData.end();
1441  {
1442  int i = 0;
1443  while( currentInductor != endInductor )
1444  {
1445  if( (getSolverState().dcopFlag) && ((*currentInductor)->ICGiven) )
1446  {
1447  inst.inductorCurrents[ i ] = (*currentInductor)->IC;
1448  }
1449  else
1450  {
1451  inst.inductorCurrents[ i ] = solVec[ (*currentInductor)->li_Branch ];
1452  }
1453  ++i;
1454  ++currentInductor;
1455  }
1456  }
1457  }
1458 
1459  return true;
1460 }
1461 
1462 //-----------------------------------------------------------------------------
1463 // Function : Master::loadDAEVectors
1464 // Purpose :
1465 // Special Notes :
1466 // Scope : public
1467 // Creator : Eric Keiter, SNL
1468 // Creation Date : 12/12/08
1469 //-----------------------------------------------------------------------------
1470 bool Master::loadDAEVectors (double * solVec, double * fVec, double *qVec, double * storeLeadF, double * storeLeadQ)
1471 {
1472  for (InstanceVector::const_iterator it = getInstanceBegin(); it != getInstanceEnd(); ++it)
1473  {
1474  Instance & inst = *(*it);
1475 
1476  // F-vector:
1477  // loop over each inductor and load it's F vector components
1478  std::vector< InductorInstanceData* >::iterator currentInductor = inst.instanceData.begin();
1479  std::vector< InductorInstanceData* >::iterator endInductor = inst.instanceData.end();
1480 
1481  while( currentInductor != endInductor )
1482  {
1483  double current = solVec[(*currentInductor)->li_Branch];
1484  double vNodePos = solVec[(*currentInductor)->li_Pos];
1485  double vNodeNeg = solVec[(*currentInductor)->li_Neg];
1486 
1487  fVec[((*currentInductor)->li_Pos)] += inst.scalingRHS * current;
1488  fVec[((*currentInductor)->li_Neg)] += -inst.scalingRHS * current;
1489  fVec[((*currentInductor)->li_Branch)] += -(vNodePos - vNodeNeg);
1490 
1491  ++currentInductor;
1492  }
1493 
1494  // Q-vector:
1495  // calculate the following product
1496  // I = column vector of currents
1497  // L = row vector of inductances
1498  // LO = matrix = mutualCup * sqrt( L' * L )
1499  // LOI = column vector = mutualCup * sqrt( L' * L ) * I
1500  // LOI[1] = mutualCup * sqrt(L[1]*L[1])*I[1]) +
1501  // mutualCup * sqrt(L[1]*L[2])*I[2]) + ...
1502  // mutualCup * sqrt(L[1]*L[n])*I[n])
1503 
1504  for( int i = 0; i < inst.numInductors; ++i )
1505  {
1506  inst.LOI[ i ] = 0.0;
1507  for( int j = 0; j < inst.numInductors; ++j )
1508  {
1509  inst.LOI[i] += inst.mutualCouplingCoef[i][j] * inst.LO[i][j] * inst.inductorCurrents[j];
1510  }
1511  }
1512 
1513  // loop over each inductor and load it's Q vector components
1514  currentInductor = inst.instanceData.begin();
1515  endInductor = inst.instanceData.end();
1516  int li=0;
1517  while( currentInductor != endInductor )
1518  {
1519  qVec[(*currentInductor)->li_Branch] += inst.LOI[ li ];
1520  ++li;
1521  ++currentInductor;
1522  }
1523  }
1524  return true;
1525 }
1526 
1527 #ifndef Xyce_NONPOINTER_MATRIX_LOAD
1528 //-----------------------------------------------------------------------------
1529 // Function : Master::loadDAEMatrices
1530 // Purpose :
1531 // Special Notes :
1532 // Scope : public
1533 // Creator : Eric Keiter, SNL
1534 // Creation Date : 12/12/08
1535 //-----------------------------------------------------------------------------
1536 bool Master::loadDAEMatrices (N_LAS_Matrix & dFdx, N_LAS_Matrix & dQdx)
1537 {
1538  for (InstanceVector::const_iterator it = getInstanceBegin(); it != getInstanceEnd(); ++it)
1539  {
1540  Instance & inst = *(*it);
1541  // loop over each inductor and load it's dFdx components
1542  std::vector< InductorInstanceData* >::iterator currentInductor = inst.instanceData.begin();
1543  std::vector< InductorInstanceData* >::iterator endInductor = inst.instanceData.end();
1544  while( currentInductor != endInductor )
1545  {
1546  *((*currentInductor)->f_PosEquBraVarPtr) += inst.scalingRHS;
1547  *((*currentInductor)->f_NegEquBraVarPtr) += -inst.scalingRHS;
1548  *((*currentInductor)->f_BraEquPosNodePtr) += -1.0;
1549  *((*currentInductor)->f_BraEquNegNodePtr) += 1.0;
1550 
1551  ++currentInductor;
1552  }
1553 
1554  // During the DCOP, the dQdx*pdt matrix is summed into the Jacobian,
1555  // even though the Q*pdt vector is not summed into the residual.
1556  //if (!getSolverState().dcopFlag)
1557  {
1558  // loop over each inductor and load it's Q vector components
1559  currentInductor = inst.instanceData.begin();
1560  endInductor = inst.instanceData.end();
1561  int li=0;
1562  while( currentInductor != endInductor )
1563  {
1564  for( int j=0; j<inst.numInductors; ++j )
1565  {
1566  *((*currentInductor)->q_inductorCurrentPtrs[j]) += inst.mutualCouplingCoef[li][j] * inst.LO[li][j];
1567  }
1568  // finally do all the dependent variable terms
1569  int numdepterms=(*currentInductor)->depVarPairs.size();
1570  for (int j=0 ; j<numdepterms; ++j)
1571  {
1572  int coefficientNumber=(*currentInductor)->depVarPairs[j].first;
1573  int depVarNumber=(*currentInductor)->depVarPairs[j].second;
1574  int otherInductor;
1575 
1576  // indexPairs[coefficient] gives the two inductors coupled by that
1577  // coefficient. One of them is currentInductor because we saved that
1578  // coefficient number in our depVarPairs, so the other is the one
1579  // we need to know
1580 
1581  if (li==inst.indexPairs[coefficientNumber].first)
1582  otherInductor=inst.indexPairs[coefficientNumber].second;
1583  else
1584  otherInductor=inst.indexPairs[coefficientNumber].second;
1585 
1586  *((*currentInductor)->q_BraEquDepVarPtrs[j]) +=
1587  inst.couplingCoefficientVarDerivs[coefficientNumber][depVarNumber]
1588  *inst.LO[li][otherInductor]*inst.inductorCurrents[otherInductor];
1589  }
1590 
1591  ++li;
1592  ++currentInductor;
1593  }
1594  }
1595  }
1596 
1597  return true;
1598 }
1599 
1600 #else
1601 //-----------------------------------------------------------------------------
1602 // Function : Master::loadDAEMatrices
1603 // Purpose :
1604 // Special Notes :
1605 // Scope : public
1606 // Creator : Eric Keiter, SNL
1607 // Creation Date : 12/12/08
1608 //-----------------------------------------------------------------------------
1609 bool Master::loadDAEMatrices (N_LAS_Matrix & dFdx, N_LAS_Matrix & dQdx)
1610 {
1611  int sizeInstances = instanceContainer_.size();
1612  for (int i=0; i<sizeInstances; ++i)
1613  {
1614  Instance & inst = *(instanceContainer_.at(i));
1615  // loop over each inductor and load it's dFdx components
1616  std::vector< InductorInstanceData* >::iterator currentInductor = inst.instanceData.begin();
1617  std::vector< InductorInstanceData* >::iterator endInductor = inst.instanceData.end();
1618  while( currentInductor != endInductor )
1619  {
1620  dFdx[((*currentInductor)->li_Pos)] [((*currentInductor)->APosEquBraVarOffset)] += inst.scalingRHS;
1621  dFdx[((*currentInductor)->li_Neg)] [((*currentInductor)->ANegEquBraVarOffset)] += -inst.scalingRHS;
1622  dFdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquPosNodeOffset)] += -1.0;
1623  dFdx[((*currentInductor)->li_Branch)][((*currentInductor)->ABraEquNegNodeOffset)] += 1.0;
1624 
1625  ++currentInductor;
1626  }
1627 
1628  // During the DCOP, the dQdx*pdt matrix is summed into the Jacobian,
1629  // even though the Q*pdt vector is not summed into the residual.
1630  //if (!getSolverState().dcopFlag)
1631  {
1632  // loop over each inductor and load it's Q vector components
1633  currentInductor = inst.instanceData.begin();
1634  endInductor = inst.instanceData.end();
1635  int i=0;
1636  while( currentInductor != endInductor )
1637  {
1638  for( int j=0; j<inst.numInductors; ++j )
1639  {
1640  dQdx[((*currentInductor)->li_Branch)]
1641  [(*currentInductor)->inductorCurrentOffsets[j]] += inst.mutualCouplingCoef[i][j] * inst.LO[i][j];
1642  }
1643  // finally do all the dependent variable terms
1644  int numdepterms=(*currentInductor)->depVarPairs.size();
1645  for (int j=0 ; j<numdepterms; ++j)
1646  {
1647  int coefficientNumber=(*currentInductor)->depVarPairs[j].first;
1648  int depVarNumber=(*currentInductor)->depVarPairs[j].second;
1649  int otherInductor;
1650 
1651  // indexPairs[coefficient] gives the two inductors coupled by that
1652  // coefficient. One of them is currentInductor because we saved that
1653  // coefficient number in our depVarPairs, so the other is the one
1654  // we need to know
1655 
1656  if (i==inst.indexPairs[coefficientNumber].first)
1657  otherInductor=inst.indexPairs[coefficientNumber].second;
1658  else
1659  otherInductor=inst.indexPairs[coefficientNumber].second;
1660 
1661  dQdx[((*currentInductor)->li_Branch)][(*currentInductor)->ABraEquDepVarOffsets[j]] +=
1662  inst.couplingCoefficientVarDerivs[coefficientNumber][depVarNumber]
1663  *inst.LO[i][otherInductor]*inst.inductorCurrents[otherInductor];
1664  }
1665 
1666  ++i;
1667  ++currentInductor;
1668  }
1669  }
1670  }
1671  return true;
1672 }
1673 #endif
1674 
1675 Device *Traits::factory(const Configuration &configuration, const FactoryBlock &factory_block)
1676 {
1677 
1678  return new Master(configuration, factory_block, factory_block.solverState_, factory_block.deviceOptions_);
1679 }
1680 
1682 {
1684  .registerDevice("mil", 1)
1685  .registerModelType("mil", 1);
1686 }
1687 
1688 } // namespace MutIndLin
1689 } // namespace Device
1690 } // namespace Xyce