Qucs-core  0.0.18
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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$
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 namespace qucs {
49 
50 using namespace transient;
51 
52 /* The function calculates the integration coefficient for numerical
53  integration methods. Supported methods are: Gear (order 1-6),
54  Trapezoidal, backward Euler and Adams-Moulton (order 1-6). */
55 void transient::calcCorrectorCoeff (int Method, int order,
56  nr_double_t * coefficients,
57  nr_double_t * delta) {
58 
59  tmatrix<nr_double_t> A (order + 1);
60  tvector<nr_double_t> x (order + 1);
61  tvector<nr_double_t> b (order + 1);
62  eqnsys<nr_double_t> e;
63  e.setAlgo (ALGO_LU_DECOMPOSITION);
64 
65  switch (Method) {
66  case INTEGRATOR_GEAR: // GEAR order 1 to 6
67  {
68 #if FIXEDCOEFF
69  int i, r, c;
70  // right hand side vector
71  for (i = 0; i < order + 1; i++) b.set (i, 1);
72  for (i = 1; i < order + 1; i++) {
73  A.set (i, 0, i); // first column
74  A.set (0, i, 1); // first row
75  }
76  for (c = 1; c <= order - 1; c++) {
77  nr_double_t entry = -c;
78  for (r = 1; r <= order; r++) {
79  A.set (r, c + 1, entry);
80  entry *= -c;
81  }
82  }
83  e.passEquationSys (&A, &x, &b);
84  e.solve ();
85 
86  // vector x consists of b_{-1}, a_{0}, a_{1} ... a_{k-1} right here
87 #if COEFFDEBUG
88  logprint (LOG_STATUS, "DEBUG: Gear order %d:", order);
89  for (i = 0; i < x.getRows (); i++) {
90  logprint (LOG_STATUS, " %g", x.get (i));
91  }
92  logprint (LOG_STATUS, "\n");
93 #endif
94  nr_double_t k = x.get (0);
95  coefficients[COEFF_G] = 1 / delta[0] / k;
96  for (i = 1; i <= order; i++) {
97  coefficients[i] = - 1 / delta[0] / k * x.get (i);
98  }
99 #else /* !FIXEDCOEFF */
100  int c, r;
101  // right hand side vector
102  b.set (1, -1 / delta[0]);
103  // first row
104  for (c = 0; c < order + 1; c++) A.set (0, c, 1);
105  nr_double_t f, a;
106  for (f = 0, c = 0; c < order; c++) {
107  f += delta[c];
108  for (a = 1, r = 0; r < order; r++) {
109  a *= f / delta[0];
110  A.set (r + 1, c + 1, a);
111  }
112  }
113  e.passEquationSys (&A, &x, &b);
114  e.solve ();
115  for (r = 0; r <= order; r++) coefficients[r] = x.get (r);
116 #endif /* !FIXEDCOEFF */
117  }
118  break;
119  case INTEGRATOR_EULER: // BACKWARD EULER
120  coefficients[COEFF_G] = 1 / delta[0];
121  coefficients[1] = - 1 / delta[0];
122  break;
123  case INTEGRATOR_TRAPEZOIDAL: // TRAPEZOIDAL (bilinear)
124  coefficients[COEFF_G] = 2 / delta[0];
125  coefficients[1] = - 2 / delta[0];
126  break;
127  case INTEGRATOR_ADAMSMOULTON: // ADAMS-MOULTON order 1 to 6
128  {
129  int i, r, c;
130  // right hand side vector
131  for (i = 0; i < order + 1; i++) b.set (i, 1);
132  for (i = 1; i < order + 1; i++) {
133  A.set (i, 1, i); // second column
134  A.set (1, i, 1); // second row
135  }
136  A.set (0, 0, 1);
137  for (c = 1; c <= order - 2; c++) {
138  nr_double_t entry = -c;
139  for (r = 2; r <= order; r++) {
140  A.set (r, c + 2, r * entry);
141  entry *= -c;
142  }
143  }
144  e.passEquationSys (&A, &x, &b);
145  e.solve ();
146 
147  // vector x consists of a_{0}, b_{-1}, b_{0} ... b_{k-2} right here
148 #if COEFFDEBUG
149  logprint (LOG_STATUS, "DEBUG: Moulton order %d:", order);
150  for (i = 0; i < x.getRows (); i++) {
151  logprint (LOG_STATUS, " %g", x.get (i));
152  }
153  logprint (LOG_STATUS, "\n");
154 #endif
155  nr_double_t k = x.get (1);
156  coefficients[COEFF_G] = 1 / delta[0] / k;
157  coefficients[1] = -x.get (0) / delta[0] / k;
158  for (i = 2; i <= order; i++) {
159  coefficients[i] = -x.get (i) / k;
160  }
161  }
162  break;
163  }
164 }
165 
166 /* The function calculates the integration coefficient for numerical
167  integration methods. Supported methods are: Adams-Bashford (order
168  1-6), forward Euler and explicit Gear (order 1-6). */
169 void transient::calcPredictorCoeff (int Method, int order,
170  nr_double_t * coefficients,
171  nr_double_t * delta) {
172 
173  tmatrix<nr_double_t> A (order + 1);
174  tvector<nr_double_t> x (order + 1);
175  tvector<nr_double_t> b (order + 1);
176  eqnsys<nr_double_t> e;
177  e.setAlgo (ALGO_LU_DECOMPOSITION);
178 
179  switch (Method) {
180  case INTEGRATOR_GEAR: // explicit GEAR order 1 to 6
181  {
182  int c, r;
183  // right hand side vector
184  b.set (0, 1);
185  // first row
186  for (c = 0; c < order + 1; c++) A.set (0, c, 1);
187  nr_double_t f, a;
188  for (f = 0, c = 0; c < order + 1; c++) {
189  f += delta[c];
190  for (a = 1, r = 0; r < order; r++) {
191  a *= f / delta[0];
192  A.set (r + 1, c, a);
193  }
194  }
195  e.passEquationSys (&A, &x, &b);
196  e.solve ();
197  for (r = 0; r <= order; r++) coefficients[r] = x.get (r);
198  }
199  break;
200  case INTEGRATOR_ADAMSBASHFORD: // ADAMS-BASHFORD order 1 to 6
201  {
202  int i, r, c;
203  // right hand side vector
204  for (i = 0; i < order + 1; i++) b.set (i, 1);
205  for (i = 1; i < order + 1; i++) A.set (1, i, 1); // second row
206  A.set (0, 0, 1);
207  for (c = 1; c <= order - 1; c++) {
208  nr_double_t entry = -c;
209  for (r = 2; r <= order; r++) {
210  A.set (r, c + 1, r * entry);
211  entry *= -c;
212  }
213  }
214  e.passEquationSys (&A, &x, &b);
215  e.solve ();
216 
217  // vector x consists of a_{0}, b_{0}, b_{1} ... b_{k-1} right here
218 #if COEFFDEBUG
219  logprint (LOG_STATUS, "DEBUG: Bashford order %d:", order);
220  for (i = 0; i < x.getRows (); i++) {
221  logprint (LOG_STATUS, " %g", x.get (i));
222  }
223  logprint (LOG_STATUS, "\n");
224 #endif
225  coefficients[COEFF_G] = x.get (0);
226  for (i = 1; i <= order; i++) {
227  coefficients[i] = x.get (i) * delta[0];
228  }
229 #if !FIXEDCOEFF
230  if (order == 2) {
231  nr_double_t f = - delta[0] / (2 * delta[1]);
232  coefficients[0] = 1;
233  coefficients[1] = (1 - f) * delta[0];
234  coefficients[2] = f * delta[0];
235  }
236 #endif
237  }
238  break;
239  case INTEGRATOR_EULER: // FORWARD EULER
240  coefficients[COEFF_G] = 1;
241  coefficients[1] = delta[0];
242  break;
243  }
244 }
245 
246 // Loads the equivalent conductance.
247 void transient::getConductance (integrator * c, nr_double_t cap,
248  nr_double_t& geq) {
249  nr_double_t * coeff = c->getCoefficients ();
250  geq = cap * coeff[COEFF_G];
251 }
252 
253 // This is the implicit Euler integrator.
254 void transient::integrateEuler (integrator * c, int qstate, nr_double_t cap,
255  nr_double_t& geq, nr_double_t& ceq) {
256  nr_double_t * coeff = c->getCoefficients ();
257  int cstate = qstate + 1;
258  nr_double_t cur;
259  geq = cap * coeff[COEFF_G];
260  ceq = c->getState (qstate, 1) * coeff[1];
261  cur = c->getState (qstate) * coeff[COEFF_G] + ceq;
262  c->setState (cstate, cur);
263 }
264 
265 // Trapezoidal integrator.
266 void transient::integrateBilinear (integrator * c, int qstate, nr_double_t cap,
267  nr_double_t& geq, nr_double_t& ceq) {
268  nr_double_t * coeff = c->getCoefficients ();
269  int cstate = qstate + 1;
270  nr_double_t cur;
271  geq = cap * coeff[COEFF_G];
272  ceq = c->getState (qstate, 1) * coeff[1] - c->getState (cstate, 1);
273  cur = c->getState (qstate) * coeff[COEFF_G] + ceq;
274  c->setState (cstate, cur);
275 }
276 
277 // Integrator using the Gear coefficients.
278 void transient::integrateGear (integrator * c, int qstate, nr_double_t cap,
279  nr_double_t& geq, nr_double_t& ceq) {
280  nr_double_t * coeff = c->getCoefficients ();
281  int i, cstate = qstate + 1;
282  nr_double_t cur;
283  geq = cap * coeff[COEFF_G];
284  for (ceq = 0, i = 1; i <= c->getOrder (); i++) {
285  ceq += c->getState (qstate, i) * coeff[i];
286  }
287  cur = c->getState (qstate) * coeff[COEFF_G] + ceq;
288  c->setState (cstate, cur);
289 }
290 
291 // Integrator using the Adams-Moulton coefficients.
292 void transient::integrateMoulton (integrator * c, int qstate, nr_double_t cap,
293  nr_double_t& geq, nr_double_t& ceq) {
294  nr_double_t * coeff = c->getCoefficients ();
295  int i, cstate = qstate + 1;
296  nr_double_t cur;
297  geq = cap * coeff[COEFF_G];
298  ceq = c->getState (qstate, 1) * coeff[1];
299  for (i = 2; i <= c->getOrder (); i++) {
300  ceq += c->getState (cstate, i - 1) * coeff[i];
301  }
302  cur = c->getState (qstate) * coeff[COEFF_G] + ceq;
303  c->setState (cstate, cur);
304 }
305 
306 /* The function applies the appropriate integration function to the
307  given circuit object. */
308 void transient::setIntegrationMethod (circuit * c, int Method) {
309  switch (Method) {
310  case INTEGRATOR_GEAR:
311  c->setIntegration (integrateGear);
312  break;
314  c->setIntegration (integrateBilinear);
315  break;
316  case INTEGRATOR_EULER:
317  c->setIntegration (integrateEuler);
318  break;
320  c->setIntegration (integrateMoulton);
321  break;
322  default:
323  c->setIntegration (NULL);
324  break;
325  }
326  c->setConductance (getConductance);
327 }
328 
329 /* Returns an appropriate integrator type identifier and the maximum
330  order depending on the given string argument. */
331 int transient::correctorType (char * Method, int& MaxOrder) {
332  if (!strcmp (Method, "Gear")) {
333  if (MaxOrder > 6) MaxOrder = 6;
334  if (MaxOrder < 1) MaxOrder = 1;
335  return INTEGRATOR_GEAR;
336  }
337  else if (!strcmp (Method, "Trapezoidal")) {
338  MaxOrder = 2;
339  return INTEGRATOR_TRAPEZOIDAL;
340  }
341  else if (!strcmp (Method, "Euler")) {
342  MaxOrder = 1;
343  return INTEGRATOR_EULER;
344  }
345  else if (!strcmp (Method, "AdamsMoulton")) {
346  if (MaxOrder > 6) MaxOrder = 6;
347  if (MaxOrder < 1) MaxOrder = 1;
349  }
350  else if (!strcmp (Method, "AdamsBashford")) {
351  if (MaxOrder > 6) MaxOrder = 6;
352  if (MaxOrder < 1) MaxOrder = 1;
354  }
355  return INTEGRATOR_UNKNOWN;
356 }
357 
358 /* The function returns the appropriate predictor integration method
359  for the given corrector method and adjusts the order of the
360  predictor as well based on the given corrector method. */
361 int transient::predictorType (int corrMethod, int corrOrder, int& predOrder) {
362  int predMethod = INTEGRATOR_UNKNOWN;
363  switch (corrMethod) {
364  case INTEGRATOR_GEAR:
365  predMethod = INTEGRATOR_GEAR;
366  break;
368  predMethod = INTEGRATOR_ADAMSBASHFORD;
369  break;
371  predMethod = INTEGRATOR_ADAMSBASHFORD;
372  break;
373  case INTEGRATOR_EULER:
374  predMethod = INTEGRATOR_EULER;
375  break;
376  }
377  predOrder = corrOrder;
378  return predMethod;
379 }
380 
381 // Structure defining integration algorithm for each possible order.
383  int Method;
384  int integratorType[6];
385  nr_double_t corrErrorConstant[6];
386  nr_double_t predErrorConstant[6];
387 };
388 
391  { INTEGRATOR_EULER },
392  { -1.0/2 },
393  { +1.0/2 }
394  },
397  { -1.0/2, -1.0/12 },
398  { +1.0/2, +5.0/12 }
399  },
400  { INTEGRATOR_GEAR,
402  INTEGRATOR_GEAR, INTEGRATOR_GEAR, INTEGRATOR_GEAR },
403  { -1.0/2, -2.0/9, -3.0/22, -12.0/125, -10.0/137, -20.0/343 },
404  { +1.0, +1.0, +1.0, +1.0, +1.0, +1.0 }
405  },
409  INTEGRATOR_ADAMSMOULTON, INTEGRATOR_ADAMSMOULTON },
410  { -1.0/2, -1.0/12, -1.0/24, -19.0/720, -3.0/160, -863.0/60480 },
411  { +1.0/2, +1.0/12, +1.0/24, +19.0/720, +3.0/160, +863.0/60480 }
412  },
416  INTEGRATOR_ADAMSBASHFORD, INTEGRATOR_ADAMSBASHFORD },
417  { -1.0/2, -5.0/12, -3.0/8, -251.0/720, -95.0/288, -19087.0/60480 },
418  { +1.0/2, +5.0/12, +3.0/8, +251.0/720, +95.0/288, +19087.0/60480 }
419  }
420 };
421 
422 /* The function returns the appropriate integration type for the given
423  corrector integration type and order. */
424 int transient::correctorType (int Method, int order) {
425  return integration_types[Method].integratorType[order - 1];
426 }
427 
428 // Returns the error constant for the given corrector.
429 nr_double_t transient::getCorrectorError (int Method, int order) {
430  return integration_types[Method].corrErrorConstant[order - 1];
431 }
432 
433 // Returns the error constant for the given predictor.
434 nr_double_t transient::getPredictorError (int Method, int order) {
435  return integration_types[Method].predErrorConstant[order - 1];
436 }
437 
438 }
nr_double_t getCorrectorError(int, int)
Definition: transient.cpp:429
int correctorType(char *, int &)
Definition: transient.cpp:331
void calcCorrectorCoeff(int, int, nr_double_t *, nr_double_t *)
Definition: transient.cpp:55
void getConductance(integrator *, nr_double_t, nr_double_t &)
Definition: transient.cpp:247
r
Definition: parse_mdl.y:515
void setIntegrationMethod(circuit *, int)
Definition: transient.cpp:308
nr_double_t getPredictorError(int, int)
Definition: transient.cpp:434
i
Definition: parse_mdl.y:516
void calcPredictorCoeff(int, int, nr_double_t *, nr_double_t *)
Definition: transient.cpp:169
void integrateEuler(integrator *, int, nr_double_t, nr_double_t &, nr_double_t &)
Definition: transient.cpp:254
x
Definition: parse_mdl.y:498
int predictorType(int, int, int &)
Definition: transient.cpp:361
The circuit class header file.
void integrateMoulton(integrator *, int, nr_double_t, nr_double_t &, nr_double_t &)
Definition: transient.cpp:292
nr_double_t corrErrorConstant[6]
Definition: transient.cpp:385
static struct integration_types_t integration_types[]
Definition: transient.cpp:389
#define A(a)
Definition: eqndefined.cpp:72
nr_double_t predErrorConstant[6]
Definition: transient.cpp:386
#define LOG_STATUS
Definition: logging.h:29
void integrateGear(integrator *, int, nr_double_t, nr_double_t &, nr_double_t &)
Definition: transient.cpp:278
void logprint(int level, const char *format,...)
Definition: logging.c:37
void integrateBilinear(integrator *, int, nr_double_t, nr_double_t &, nr_double_t &)
Definition: transient.cpp:266
#define COEFF_G
Definition: transient.cpp:46