My Project  0.0.16
QUCS Mapping
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
transient.cpp
Go to the documentation of this file.
1 /*
2  * transient.cpp - transient helper class implementation
3  *
4  * Copyright (C) 2004, 2006 Stefan Jahn <stefan@lkcc.org>
5  *
6  * This is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2, or (at your option)
9  * any later version.
10  *
11  * This software is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with this package; see the file COPYING. If not, write to
18  * the Free Software Foundation, Inc., 51 Franklin Street - Fifth Floor,
19  * Boston, MA 02110-1301, USA.
20  *
21  * $Id: transient.cpp 1825 2011-03-11 20:42:14Z ela $
22  *
23  */
24 
25 #if HAVE_CONFIG_H
26 # include <config.h>
27 #endif
28 
29 #include <stdio.h>
30 #include <stdlib.h>
31 #include <string.h>
32 
33 #include "object.h"
34 #include "complex.h"
35 #include "circuit.h"
36 #include "net.h"
37 #include "tvector.h"
38 #include "tmatrix.h"
39 #include "eqnsys.h"
40 #include "transient.h"
41 
42 #define COEFFDEBUG 0
43 #define FIXEDCOEFF 0
44 
45 // Defines where the equivalent admittance coefficient is going to be stored.
46 #define COEFF_G 0
47 
48 using namespace transient;
49 
50 /* The function calculates the integration coefficient for numerical
51  integration methods. Supported methods are: Gear (order 1-6),
52  Trapezoidal, backward Euler and Adams-Moulton (order 1-6). */
53 void transient::calcCorrectorCoeff (int Method, int order,
54  nr_double_t * coefficients,
55  nr_double_t * delta) {
56 
57  tmatrix<nr_double_t> A (order + 1);
58  tvector<nr_double_t> x (order + 1);
59  tvector<nr_double_t> b (order + 1);
62 
63  switch (Method) {
64  case INTEGRATOR_GEAR: // GEAR order 1 to 6
65  {
66 #if FIXEDCOEFF
67  int i, r, c;
68  // right hand side vector
69  for (i = 0; i < order + 1; i++) b.set (i, 1);
70  for (i = 1; i < order + 1; i++) {
71  A.set (i, 0, i); // first column
72  A.set (0, i, 1); // first row
73  }
74  for (c = 1; c <= order - 1; c++) {
75  nr_double_t entry = -c;
76  for (r = 1; r <= order; r++) {
77  A.set (r, c + 1, entry);
78  entry *= -c;
79  }
80  }
81  e.passEquationSys (&A, &x, &b);
82  e.solve ();
83 
84  // vector x consists of b_{-1}, a_{0}, a_{1} ... a_{k-1} right here
85 #if COEFFDEBUG
86  logprint (LOG_STATUS, "DEBUG: Gear order %d:", order);
87  for (i = 0; i < x.getRows (); i++) {
88  logprint (LOG_STATUS, " %g", x.get (i));
89  }
90  logprint (LOG_STATUS, "\n");
91 #endif
92  nr_double_t k = x.get (0);
93  coefficients[COEFF_G] = 1 / delta[0] / k;
94  for (i = 1; i <= order; i++) {
95  coefficients[i] = - 1 / delta[0] / k * x.get (i);
96  }
97 #else /* !FIXEDCOEFF */
98  int c, r;
99  // right hand side vector
100  b.set (1, -1 / delta[0]);
101  // first row
102  for (c = 0; c < order + 1; c++) A.set (0, c, 1);
103  nr_double_t f, a;
104  for (f = 0, c = 0; c < order; c++) {
105  f += delta[c];
106  for (a = 1, r = 0; r < order; r++) {
107  a *= f / delta[0];
108  A.set (r + 1, c + 1, a);
109  }
110  }
111  e.passEquationSys (&A, &x, &b);
112  e.solve ();
113  for (r = 0; r <= order; r++) coefficients[r] = x.get (r);
114 #endif /* !FIXEDCOEFF */
115  }
116  break;
117  case INTEGRATOR_EULER: // BACKWARD EULER
118  coefficients[COEFF_G] = 1 / delta[0];
119  coefficients[1] = - 1 / delta[0];
120  break;
121  case INTEGRATOR_TRAPEZOIDAL: // TRAPEZOIDAL (bilinear)
122  coefficients[COEFF_G] = 2 / delta[0];
123  coefficients[1] = - 2 / delta[0];
124  break;
125  case INTEGRATOR_ADAMSMOULTON: // ADAMS-MOULTON order 1 to 6
126  {
127  int i, r, c;
128  // right hand side vector
129  for (i = 0; i < order + 1; i++) b.set (i, 1);
130  for (i = 1; i < order + 1; i++) {
131  A.set (i, 1, i); // second column
132  A.set (1, i, 1); // second row
133  }
134  A.set (0, 0, 1);
135  for (c = 1; c <= order - 2; c++) {
136  nr_double_t entry = -c;
137  for (r = 2; r <= order; r++) {
138  A.set (r, c + 2, r * entry);
139  entry *= -c;
140  }
141  }
142  e.passEquationSys (&A, &x, &b);
143  e.solve ();
144 
145  // vector x consists of a_{0}, b_{-1}, b_{0} ... b_{k-2} right here
146 #if COEFFDEBUG
147  logprint (LOG_STATUS, "DEBUG: Moulton order %d:", order);
148  for (i = 0; i < x.getRows (); i++) {
149  logprint (LOG_STATUS, " %g", x.get (i));
150  }
151  logprint (LOG_STATUS, "\n");
152 #endif
153  nr_double_t k = x.get (1);
154  coefficients[COEFF_G] = 1 / delta[0] / k;
155  coefficients[1] = -x.get (0) / delta[0] / k;
156  for (i = 2; i <= order; i++) {
157  coefficients[i] = -x.get (i) / k;
158  }
159  }
160  break;
161  }
162 }
163 
164 /* The function calculates the integration coefficient for numerical
165  integration methods. Supported methods are: Adams-Bashford (order
166  1-6), forward Euler and explicit Gear (order 1-6). */
167 void transient::calcPredictorCoeff (int Method, int order,
168  nr_double_t * coefficients,
169  nr_double_t * delta) {
170 
171  tmatrix<nr_double_t> A (order + 1);
172  tvector<nr_double_t> x (order + 1);
173  tvector<nr_double_t> b (order + 1);
176 
177  switch (Method) {
178  case INTEGRATOR_GEAR: // explicit GEAR order 1 to 6
179  {
180  int c, r;
181  // right hand side vector
182  b.set (0, 1);
183  // first row
184  for (c = 0; c < order + 1; c++) A.set (0, c, 1);
185  nr_double_t f, a;
186  for (f = 0, c = 0; c < order + 1; c++) {
187  f += delta[c];
188  for (a = 1, r = 0; r < order; r++) {
189  a *= f / delta[0];
190  A.set (r + 1, c, a);
191  }
192  }
193  e.passEquationSys (&A, &x, &b);
194  e.solve ();
195  for (r = 0; r <= order; r++) coefficients[r] = x.get (r);
196  }
197  break;
198  case INTEGRATOR_ADAMSBASHFORD: // ADAMS-BASHFORD order 1 to 6
199  {
200  int i, r, c;
201  // right hand side vector
202  for (i = 0; i < order + 1; i++) b.set (i, 1);
203  for (i = 1; i < order + 1; i++) A.set (1, i, 1); // second row
204  A.set (0, 0, 1);
205  for (c = 1; c <= order - 1; c++) {
206  nr_double_t entry = -c;
207  for (r = 2; r <= order; r++) {
208  A.set (r, c + 1, r * entry);
209  entry *= -c;
210  }
211  }
212  e.passEquationSys (&A, &x, &b);
213  e.solve ();
214 
215  // vector x consists of a_{0}, b_{0}, b_{1} ... b_{k-1} right here
216 #if COEFFDEBUG
217  logprint (LOG_STATUS, "DEBUG: Bashford order %d:", order);
218  for (i = 0; i < x.getRows (); i++) {
219  logprint (LOG_STATUS, " %g", x.get (i));
220  }
221  logprint (LOG_STATUS, "\n");
222 #endif
223  coefficients[COEFF_G] = x.get (0);
224  for (i = 1; i <= order; i++) {
225  coefficients[i] = x.get (i) * delta[0];
226  }
227 #if !FIXEDCOEFF
228  if (order == 2) {
229  nr_double_t f = - delta[0] / (2 * delta[1]);
230  coefficients[0] = 1;
231  coefficients[1] = (1 - f) * delta[0];
232  coefficients[2] = f * delta[0];
233  }
234 #endif
235  }
236  break;
237  case INTEGRATOR_EULER: // FORWARD EULER
238  coefficients[COEFF_G] = 1;
239  coefficients[1] = delta[0];
240  break;
241  }
242 }
243 
244 // Loads the equivalent conductance.
245 void transient::getConductance (integrator * c, nr_double_t cap,
246  nr_double_t& geq) {
247  nr_double_t * coeff = c->getCoefficients ();
248  geq = cap * coeff[COEFF_G];
249 }
250 
251 // This is the implicit Euler integrator.
252 void transient::integrateEuler (integrator * c, int qstate, nr_double_t cap,
253  nr_double_t& geq, nr_double_t& ceq) {
254  nr_double_t * coeff = c->getCoefficients ();
255  int cstate = qstate + 1;
256  nr_double_t cur;
257  geq = cap * coeff[COEFF_G];
258  ceq = c->getState (qstate, 1) * coeff[1];
259  cur = c->getState (qstate) * coeff[COEFF_G] + ceq;
260  c->setState (cstate, cur);
261 }
262 
263 // Trapezoidal integrator.
264 void transient::integrateBilinear (integrator * c, int qstate, nr_double_t cap,
265  nr_double_t& geq, nr_double_t& ceq) {
266  nr_double_t * coeff = c->getCoefficients ();
267  int cstate = qstate + 1;
268  nr_double_t cur;
269  geq = cap * coeff[COEFF_G];
270  ceq = c->getState (qstate, 1) * coeff[1] - c->getState (cstate, 1);
271  cur = c->getState (qstate) * coeff[COEFF_G] + ceq;
272  c->setState (cstate, cur);
273 }
274 
275 // Integrator using the Gear coefficients.
276 void transient::integrateGear (integrator * c, int qstate, nr_double_t cap,
277  nr_double_t& geq, nr_double_t& ceq) {
278  nr_double_t * coeff = c->getCoefficients ();
279  int i, cstate = qstate + 1;
280  nr_double_t cur;
281  geq = cap * coeff[COEFF_G];
282  for (ceq = 0, i = 1; i <= c->getOrder (); i++) {
283  ceq += c->getState (qstate, i) * coeff[i];
284  }
285  cur = c->getState (qstate) * coeff[COEFF_G] + ceq;
286  c->setState (cstate, cur);
287 }
288 
289 // Integrator using the Adams-Moulton coefficients.
290 void transient::integrateMoulton (integrator * c, int qstate, nr_double_t cap,
291  nr_double_t& geq, nr_double_t& ceq) {
292  nr_double_t * coeff = c->getCoefficients ();
293  int i, cstate = qstate + 1;
294  nr_double_t cur;
295  geq = cap * coeff[COEFF_G];
296  ceq = c->getState (qstate, 1) * coeff[1];
297  for (i = 2; i <= c->getOrder (); i++) {
298  ceq += c->getState (cstate, i - 1) * coeff[i];
299  }
300  cur = c->getState (qstate) * coeff[COEFF_G] + ceq;
301  c->setState (cstate, cur);
302 }
303 
304 /* The function applies the appropriate integration function to the
305  given circuit object. */
307  switch (Method) {
308  case INTEGRATOR_GEAR:
310  break;
313  break;
314  case INTEGRATOR_EULER:
316  break;
319  break;
320  default:
321  c->setIntegration (NULL);
322  break;
323  }
325 }
326 
327 /* Returns an appropriate integrator type identifier and the maximum
328  order depending on the given string argument. */
329 int transient::correctorType (char * Method, int& MaxOrder) {
330  if (!strcmp (Method, "Gear")) {
331  if (MaxOrder > 6) MaxOrder = 6;
332  if (MaxOrder < 1) MaxOrder = 1;
333  return INTEGRATOR_GEAR;
334  }
335  else if (!strcmp (Method, "Trapezoidal")) {
336  MaxOrder = 2;
337  return INTEGRATOR_TRAPEZOIDAL;
338  }
339  else if (!strcmp (Method, "Euler")) {
340  MaxOrder = 1;
341  return INTEGRATOR_EULER;
342  }
343  else if (!strcmp (Method, "AdamsMoulton")) {
344  if (MaxOrder > 6) MaxOrder = 6;
345  if (MaxOrder < 1) MaxOrder = 1;
347  }
348  else if (!strcmp (Method, "AdamsBashford")) {
349  if (MaxOrder > 6) MaxOrder = 6;
350  if (MaxOrder < 1) MaxOrder = 1;
352  }
353  return INTEGRATOR_UNKNOWN;
354 }
355 
356 /* The function returns the appropriate predictor integration method
357  for the given corrector method and adjusts the order of the
358  predictor as well based on the given corrector method. */
359 int transient::predictorType (int corrMethod, int corrOrder, int& predOrder) {
360  int predMethod = INTEGRATOR_UNKNOWN;
361  switch (corrMethod) {
362  case INTEGRATOR_GEAR:
363  predMethod = INTEGRATOR_GEAR;
364  break;
366  predMethod = INTEGRATOR_ADAMSBASHFORD;
367  break;
369  predMethod = INTEGRATOR_ADAMSBASHFORD;
370  break;
371  case INTEGRATOR_EULER:
372  predMethod = INTEGRATOR_EULER;
373  break;
374  }
375  predOrder = corrOrder;
376  return predMethod;
377 }
378 
379 // Structure defining integration algorithm for each possible order.
381  int Method;
382  int integratorType[6];
383  nr_double_t corrErrorConstant[6];
384  nr_double_t predErrorConstant[6];
385 };
386 
387 static struct integration_types_t integration_types[] = {
389  { INTEGRATOR_EULER },
390  { -1.0/2 },
391  { +1.0/2 }
392  },
395  { -1.0/2, -1.0/12 },
396  { +1.0/2, +5.0/12 }
397  },
398  { INTEGRATOR_GEAR,
400  INTEGRATOR_GEAR, INTEGRATOR_GEAR, INTEGRATOR_GEAR },
401  { -1.0/2, -2.0/9, -3.0/22, -12.0/125, -10.0/137, -20.0/343 },
402  { +1.0, +1.0, +1.0, +1.0, +1.0, +1.0 }
403  },
407  INTEGRATOR_ADAMSMOULTON, INTEGRATOR_ADAMSMOULTON },
408  { -1.0/2, -1.0/12, -1.0/24, -19.0/720, -3.0/160, -863.0/60480 },
409  { +1.0/2, +1.0/12, +1.0/24, +19.0/720, +3.0/160, +863.0/60480 }
410  },
414  INTEGRATOR_ADAMSBASHFORD, INTEGRATOR_ADAMSBASHFORD },
415  { -1.0/2, -5.0/12, -3.0/8, -251.0/720, -95.0/288, -19087.0/60480 },
416  { +1.0/2, +5.0/12, +3.0/8, +251.0/720, +95.0/288, +19087.0/60480 }
417  }
418 };
419 
420 /* The function returns the appropriate integration type for the given
421  corrector integration type and order. */
422 int transient::correctorType (int Method, int order) {
423  return integration_types[Method].integratorType[order - 1];
424 }
425 
426 // Returns the error constant for the given corrector.
427 nr_double_t transient::getCorrectorError (int Method, int order) {
428  return integration_types[Method].corrErrorConstant[order - 1];
429 }
430 
431 // Returns the error constant for the given predictor.
432 nr_double_t transient::getPredictorError (int Method, int order) {
433  return integration_types[Method].predErrorConstant[order - 1];
434 }