Qucs-core  0.0.18
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
trsolver.cpp
Go to the documentation of this file.
1 /*
2  * trsolver.cpp - transient solver class implementation
3  *
4  * Copyright (C) 2004, 2005, 2006, 2007, 2009 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 <string.h>
31 #include <math.h>
32 #include <float.h>
33 
34 #include "compat.h"
35 #include "object.h"
36 #include "logging.h"
37 #include "complex.h"
38 #include "circuit.h"
39 #include "sweep.h"
40 #include "net.h"
41 #include "netdefs.h"
42 #include "analysis.h"
43 #include "nasolver.h"
44 #include "history.h"
45 #include "trsolver.h"
46 #include "transient.h"
47 #include "exception.h"
48 #include "exceptionstack.h"
49 
50 #define STEPDEBUG 0 // set to zero for release
51 #define BREAKPOINTS 0 // exact breakpoint calculation
52 
53 #define dState 0 // delta T state
54 #define sState 1 // solution state
55 
56 // Macro for the n-th state of the solution vector history.
57 #define SOL(state) (solution[(int) getState (sState, (state))])
58 
59 namespace qucs {
60 
61 using namespace transient;
62 
63 // Constructor creates an unnamed instance of the trsolver class.
65  : nasolver<nr_double_t> (), states<nr_double_t> ()
66 {
67  swp = NULL;
69  setDescription ("transient");
70  for (int i = 0; i < 8; i++) solution[i] = NULL;
71  tHistory = NULL;
72  relaxTSR = false;
73  initialDC = true;
74 }
75 
76 // Constructor creates a named instance of the trsolver class.
78  : nasolver<nr_double_t> (n), states<nr_double_t> ()
79 {
80  swp = NULL;
82  setDescription ("transient");
83  for (int i = 0; i < 8; i++) solution[i] = NULL;
84  tHistory = NULL;
85  relaxTSR = false;
86  initialDC = true;
87 }
88 
89 // Destructor deletes the trsolver class object.
91 {
92  if (swp) delete swp;
93  for (int i = 0; i < 8; i++)
94  {
95  if (solution[i] != NULL)
96  {
97  delete solution[i];
98  }
99  }
100  if (tHistory) delete tHistory;
101 }
102 
103 /* The copy constructor creates a new instance of the trsolver class
104  based on the given trsolver object. */
106  : nasolver<nr_double_t> (o), states<nr_double_t> (o)
107 {
108  swp = o.swp ? new sweep (*o.swp) : NULL;
109  for (int i = 0; i < 8; i++) solution[i] = NULL;
110  tHistory = o.tHistory ? new history (*o.tHistory) : NULL;
111  relaxTSR = o.relaxTSR;
112  initialDC = o.initialDC;
113 }
114 
115 // This function creates the time sweep if necessary.
117 {
118  if (swp != NULL) delete swp;
119  swp = createSweep ("time");
120 }
121 
122 // Performs the initial DC analysis.
124 {
125  int error = 0;
126 
127  // First calculate a initial state using the non-linear DC analysis.
128  setDescription ("initial DC");
129  initDC ();
131  solve_pre ();
132  applyNodeset ();
133 
134  // Run the DC solver once.
135  try_running ()
136  {
137  error = solve_nonlinear ();
138  }
139  // Appropriate exception handling.
140  catch_exception ()
141  {
143  pop_exception ();
145  logprint (LOG_ERROR, "WARNING: %s: %s analysis failed, using line search "
146  "fallback\n", getName (), getDescription ());
147  applyNodeset ();
148  restart ();
149  error = solve_nonlinear ();
150  break;
151  default:
152  // Otherwise return.
153  estack.print ();
154  error++;
155  break;
156  }
157 
158  // Save the DC solution.
159  storeSolution ();
160 
161  // Cleanup nodal analysis solver.
162  solve_post ();
163 
164  // Really failed to find initial DC solution?
165  if (error)
166  {
167  logprint (LOG_ERROR, "ERROR: %s: %s analysis failed\n",
168  getName (), getDescription ());
169  }
170  return error;
171 }
172 
173 /* This is the transient netlist solver. It prepares the circuit list
174  for each requested time and solves it then. */
175 int trsolver::solve (void)
176 {
177  nr_double_t time, saveCurrent;
178  int error = 0, convError = 0;
179  char * solver = getPropertyString ("Solver");
180  relaxTSR = !strcmp (getPropertyString ("relaxTSR"), "yes") ? true : false;
181  initialDC = !strcmp (getPropertyString ("initialDC"), "yes") ? true : false;
182 
183  runs++;
184  saveCurrent = current = 0;
185  stepDelta = -1;
186  converged = 0;
187  fixpoint = 0;
189 
190  // Choose a solver.
191  if (!strcmp (solver, "CroutLU"))
193  else if (!strcmp (solver, "DoolittleLU"))
195  else if (!strcmp (solver, "HouseholderQR"))
197  else if (!strcmp (solver, "HouseholderLQ"))
199  else if (!strcmp (solver, "GolubSVD"))
201 
202  // Perform initial DC analysis.
203  if (initialDC)
204  {
205  error = dcAnalysis ();
206  if (error)
207  return -1;
208  }
209 
210  // Initialize transient analysis.
211  setDescription ("transient");
212  initTR ();
214  solve_pre ();
215 
216  // Create time sweep if necessary.
217  initSteps ();
218  swp->reset ();
219 
220  // Recall the DC solution.
221  recallSolution ();
222 
223  // Apply the nodesets and adjust previous solutions.
224  applyNodeset (false);
225  fillSolution (x);
226 
227  // Tell integrators to be initialized.
228  setMode (MODE_INIT);
229 
230  int running = 0;
231  rejected = 0;
232  delta /= 10;
234  adjustOrder (1);
235 
236  // Start to sweep through time.
237  for (int i = 0; i < swp->getSize (); i++)
238  {
239  time = swp->next ();
240  if (progress) logprogressbar (i, swp->getSize (), 40);
241 
242 #if DEBUG && 0
243  logprint (LOG_STATUS, "NOTIFY: %s: solving netlist for t = %e\n",
244  getName (), (double) time);
245 #endif
246 
247  do // while (saveCurrent < time), i.e. until a requested breakpoint is hit
248  {
249 #if STEPDEBUG
250  if (delta == deltaMin)
251  {
252  // the integrator step size has become smaller than the
253  // specified allowed minimum, Qucs is unable to solve the circuit
254  // while meeting the tolerance conditions
256  "WARNING: %s: minimum delta h = %.3e at t = %.3e\n",
257  getName (), (double) delta, (double) current);
258  }
259 #endif
260  // updates the integrator coefficients, and updates the array of prev
261  // 8 deltas with the new delta for this step
263 
264  // Run predictor to get a start value for the solution vector for
265  // the successive iterative corrector process
266  error += predictor ();
267 
268  // restart Newton iteration
269  if (rejected)
270  {
271  restart (); // restart non-linear devices
272  rejected = 0;
273  }
274 
275  // Run corrector process with appropriate exception handling.
276  // The corrector iterates through the solutions of the integration
277  // process until a certain error tolerance has been reached.
278  try_running () // #defined as: do {
279  {
280  error += corrector ();
281  }
282  catch_exception () // #defined as: } while (0); if (estack.top ()) switch (estack.top()->getCode ())
283  {
285  pop_exception ();
286 
287  // step back from the current time value to the previous time
288  if (current > 0) current -= delta;
289  // Reduce step-size (by half) if failed to converge.
290  delta /= 2;
291  if (delta <= deltaMin)
292  {
293  // but do not reduce the step size below a specified minimum
294  delta = deltaMin;
295  // instead reduce the order of the integration
296  adjustOrder (1);
297  }
298  // step forward to the new current time value
299  if (current > 0) current += delta;
300 
301  // Update statistics.
302  statRejected++;
303  statConvergence++;
304  rejected++; // mark the previous step size choice as rejected
305  converged = 0;
306  error = 0;
307 
308  // Start using damped Newton-Raphson.
310  convError = 2;
311 #if DEBUG
312  logprint (LOG_ERROR, "WARNING: delta rejected at t = %.3e, h = %.3e "
313  "(no convergence)\n", (double) saveCurrent, (double) delta);
314 #endif
315  break;
316  default:
317  // Otherwise return.
318  estack.print ();
319  error++;
320  break;
321  }
322  // return if any errors occured other than convergence failure
323  if (error) return -1;
324 
325  // if the step was rejected, the solution loop is restarted here
326  if (rejected) continue;
327 
328  // check whether Jacobian matrix is still non-singular
329  if (!A->isFinite ())
330  {
331  logprint (LOG_ERROR, "ERROR: %s: Jacobian singular at t = %.3e, "
332  "aborting %s analysis\n", getName (), (double) current,
333  getDescription ());
334  return -1;
335  }
336 
337  // Update statistics and no more damped Newton-Raphson.
339  if (--convError < 0) convHelper = 0;
340 
341  // Now advance in time or not...
342  if (running > 1)
343  {
344  adjustDelta (time);
345  adjustOrder ();
346  }
347  else
348  {
349  fillStates ();
350  nextStates ();
351  rejected = 0;
352  }
353 
354  saveCurrent = current;
355  current += delta;
356  running++;
357  converged++;
358 
359  // Tell integrators to be running.
360  setMode (MODE_NONE);
361 
362  // Initialize or update history.
363  if (running > 1)
364  {
365  updateHistory (saveCurrent);
366  }
367  else
368  {
369  initHistory (saveCurrent);
370  }
371  }
372  while (saveCurrent < time); // Hit a requested time point?
373 
374  // Save results.
375 #if STEPDEBUG
376  logprint (LOG_STATUS, "DEBUG: save point at t = %.3e, h = %.3e\n",
377  (double) saveCurrent, (double) delta);
378 #endif
379 
380 #if BREAKPOINTS
381  saveAllResults (saveCurrent);
382 #else
383  saveAllResults (time);
384 #endif
385  } // for (int i = 0; i < swp->getSize (); i++)
386 
387  solve_post ();
388  if (progress) logprogressclear (40);
389  logprint (LOG_STATUS, "NOTIFY: %s: average time-step %g, %d rejections\n",
390  getName (), (double) (saveCurrent / statSteps), statRejected);
391  logprint (LOG_STATUS, "NOTIFY: %s: average NR-iterations %g, "
392  "%d non-convergences\n", getName (),
394 
395  // cleanup
396  deinitTR ();
397  return 0;
398 }
399 
400 // The function initializes the history.
401 void trsolver::initHistory (nr_double_t t)
402 {
403  // initialize time vector
404  tHistory = new history ();
405  tHistory->append (t);
406  tHistory->self ();
407  // initialize circuit histories
408  nr_double_t age = 0.0;
409  circuit * root = subnet->getRoot ();
410  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ())
411  {
412  if (c->hasHistory ())
413  {
414  c->applyHistory (tHistory);
415  saveHistory (c);
416  if (c->getHistoryAge () > age)
417  {
418  age = c->getHistoryAge ();
419  }
420  }
421  }
422  // set maximum required age for all circuits
423  tHistory->setAge (age);
424 }
425 
426 /* The following function updates the histories for the circuits which
427  requested them. */
428 void trsolver::updateHistory (nr_double_t t)
429 {
430  if (t > tHistory->last ())
431  {
432  // update time vector
433  tHistory->append (t);
434  // update circuit histories
435  circuit * root = subnet->getRoot ();
436  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ())
437  {
438  if (c->hasHistory ()) saveHistory (c);
439  }
440  tHistory->drop ();
441  }
442 }
443 
444 // Stores node voltages and branch currents in the given circuits history.
446 {
447 
448  int N = countNodes ();
449  int r, i, s = c->getSize ();
450 
451  for (i = 0; i < s; i++)
452  {
453  // save node voltages
454  r = findAssignedNode (c, i);
455  if (r < 0)
456  // the node was not found, append a zero to the history
457  // matching this index
458  c->appendHistory (i, 0.0);
459  else
460  // the node was found, append the voltage value to
461  // that node's history
462  c->appendHistory (i, x->get (r));
463  }
464 
465  for (i = 0; i < c->getVoltageSources (); i++)
466  {
467  // save branch currents
468  r = c->getVoltageSource () + i;
469  c->appendHistory (i + s, x->get (r + N));
470  }
471 
472 }
473 
474 /* This function predicts a start value for the solution vector for
475  the successive iterative corrector process. */
477 {
478  int error = 0;
479  switch (predType)
480  {
481  case INTEGRATOR_GEAR: // explicit GEAR
482  predictGear ();
483  break;
484  case INTEGRATOR_ADAMSBASHFORD: // ADAMS-BASHFORD
485  predictBashford ();
486  break;
487  case INTEGRATOR_EULER: // FORWARD EULER
488  predictEuler ();
489  break;
490  default:
491  *x = *SOL (1); // This is too a simple predictor...
492  break;
493  }
494  saveSolution ();
495  *SOL (0) = *x;
496  return error;
497 }
498 
499 // Stores the given vector into all the solution vectors.
501 {
502  for (int i = 0; i < 8; i++) *SOL (i) = *s;
503 }
504 
505 /* The function predicts the successive solution vector using the
506  explicit Adams-Bashford integration formula. */
508 {
509  int N = countNodes ();
510  int M = countVoltageSources ();
511  nr_double_t xn, dd, hn;
512 
513  // go through each solution
514  for (int r = 0; r < N + M; r++)
515  {
516  xn = predCoeff[0] * SOL(1)->get (r); // a0 coefficient
517  for (int o = 1; o <= predOrder; o++)
518  {
519  hn = getState (dState, o); // previous time-step
520  // divided differences
521  dd = (SOL(o)->get (r) - SOL(o + 1)->get (r)) / hn;
522  xn += predCoeff[o] * dd; // b0, b1, ... coefficients
523  }
524  x->set (r, xn); // save prediction
525  }
526 }
527 
528 /* The function predicts the successive solution vector using the
529  explicit forward Euler integration formula. Actually this is
530  Adams-Bashford order 1. */
532 {
533  int N = countNodes ();
534  int M = countVoltageSources ();
535  nr_double_t xn, dd, hn;
536 
537  for (int r = 0; r < N + M; r++)
538  {
539  xn = predCoeff[0] * SOL(1)->get (r);
540  hn = getState (dState, 1);
541  dd = (SOL(1)->get (r) - SOL(2)->get (r)) / hn;
542  xn += predCoeff[1] * dd;
543  x->set (r, xn);
544  }
545 }
546 
547 /* The function predicts the successive solution vector using the
548  explicit Gear integration formula. */
550 {
551  int N = countNodes ();
552  int M = countVoltageSources ();
553  nr_double_t xn;
554 
555  // go through each solution
556  for (int r = 0; r < N + M; r++)
557  {
558  xn = 0;
559  for (int o = 0; o <= predOrder; o++)
560  {
561  // a0, a1, ... coefficients
562  xn += predCoeff[o] * SOL(o + 1)->get (r);
563  }
564  x->set (r, xn); // save prediction
565  }
566 }
567 
568 /* The function iterates through the solutions of the integration
569  process until a certain error tolerance has been reached. */
571 {
572  int error = 0;
573  error += solve_nonlinear ();
574  return error;
575 }
576 
577 // The function advances one more time-step.
579 {
580  circuit * root = subnet->getRoot ();
581  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ())
582  {
583  // for each circuit get the next state
584  c->nextState ();
585  }
586 
587  *SOL (0) = *x; // save current solution
588  nextState ();
589  statSteps++;
590 }
591 
592 /* This function stores the current state of each circuit into all
593  other states as well. It is useful for higher order integration
594  methods in order to initialize the states after the initial
595  transient solution. */
597 {
598  circuit * root = subnet->getRoot ();
599  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ())
600  {
601  for (int s = 0; s < c->getStates (); s++)
602  c->fillState (s, c->getState (s));
603  }
604 }
605 
606 // The function modifies the circuit lists integrator mode.
607 void trsolver::setMode (int state)
608 {
609  circuit * root = subnet->getRoot ();
610  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ())
611  c->setMode (state);
612 }
613 
614 // The function passes the time delta array to the circuit list.
616 {
617  circuit * root = subnet->getRoot ();
618  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ())
619  c->setDelta (deltas);
620 }
621 
622 /* This function tries to adapt the current time-step according to the
623  global truncation error. */
624 void trsolver::adjustDelta (nr_double_t t)
625 {
626  deltaOld = delta;
627  delta = checkDelta ();
628  if (delta > deltaMax) delta = deltaMax;
629  if (delta < deltaMin) delta = deltaMin;
630 
631  // delta correction in order to hit exact breakpoint
632  int good = 0;
633  if (!relaxTSR) // relaxed step raster?
634  {
635  if (!statConvergence || converged > 64) /* Is this a good guess? */
636  {
637  // check next breakpoint
638  if (stepDelta > 0.0)
639  {
640  // restore last valid delta
641  delta = stepDelta;
642  stepDelta = -1.0;
643  }
644  else
645  {
646  if (delta > (t - current) && t > current)
647  {
648  // save last valid delta and set exact step
650  delta = t - current;
651  good = 1;
652  }
653  else
654  {
655  stepDelta = -1.0;
656  }
657  }
658  if (delta > deltaMax) delta = deltaMax;
659  if (delta < deltaMin) delta = deltaMin;
660  }
661  }
662 
663  // usual delta correction
664  if (delta > 0.9 * deltaOld || good) // accept current delta
665  {
666  nextStates ();
667  rejected = 0;
668 #if STEPDEBUG
670  "DEBUG: delta accepted at t = %.3e, h = %.3e\n",
671  (double) current, (double) delta);
672 #endif
673  }
674  else if (deltaOld > delta) // reject current delta
675  {
676  rejected++;
677  statRejected++;
678 #if STEPDEBUG
680  "DEBUG: delta rejected at t = %.3e, h = %.3e\n",
681  (double) current, (double) delta);
682 #endif
683  if (current > 0) current -= deltaOld;
684  }
685  else
686  {
687  nextStates ();
688  rejected = 0;
689  }
690 }
691 
692 /* The function can be used to increase the current order of the
693  integration method or to reduce it. */
694 void trsolver::adjustOrder (int reduce)
695 {
696  if ((corrOrder < corrMaxOrder && !rejected) || reduce)
697  {
698  if (reduce)
699  {
700  corrOrder = 1;
701  }
702  else if (!rejected)
703  {
704  corrOrder++;
705  }
706 
707  // adjust type and order of corrector and predictor
710 
711  // apply new corrector method and order to each circuit
712  circuit * root = subnet->getRoot ();
713  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ())
714  {
715  c->setOrder (corrOrder);
717  }
718  }
719 }
720 
721 /* Goes through the list of circuit objects and runs its calcDC()
722  function. */
724 {
725  circuit * root = self->getNet()->getRoot ();
726  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ())
727  {
728  c->calcDC ();
729  }
730 }
731 
732 /* Goes through the list of circuit objects and runs its calcTR()
733  function. */
735 {
736  circuit * root = self->getNet()->getRoot ();
737  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ())
738  {
739  c->calcTR (self->current);
740  }
741 }
742 
743 /* Goes through the list of non-linear circuit objects and runs its
744  restartDC() function. */
745 void trsolver::restart (void)
746 {
747  circuit * root = subnet->getRoot ();
748  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ())
749  {
750  if (c->isNonLinear ()) c->restartDC ();
751  }
752 }
753 
754 /* Goes through the list of circuit objects and runs its initDC()
755  function. */
756 void trsolver::initDC (void)
757 {
758  circuit * root = subnet->getRoot ();
759  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ())
760  {
761  c->initDC ();
762  }
763 }
764 
765 /* Goes through the list of circuit objects and runs its initTR()
766  function. */
767 void trsolver::initTR (void)
768 {
769  char * IMethod = getPropertyString ("IntegrationMethod");
770  nr_double_t start = getPropertyDouble ("Start");
771  nr_double_t stop = getPropertyDouble ("Stop");
772  nr_double_t points = getPropertyDouble ("Points");
773 
774  // fetch corrector integration method and determine predicor method
775  corrMaxOrder = getPropertyInteger ("Order");
780 
781  // initialize step values
782  delta = getPropertyDouble ("InitialStep");
783  deltaMin = getPropertyDouble ("MinStep");
784  deltaMax = getPropertyDouble ("MaxStep");
785  if (deltaMax == 0.0)
786  deltaMax = MIN ((stop - start) / (points - 1), stop / 200);
787  if (deltaMin == 0.0)
788  deltaMin = NR_TINY * 10 * deltaMax;
789  if (delta == 0.0)
790  delta = MIN (stop / 200, deltaMax) / 10;
791  if (delta < deltaMin) delta = deltaMin;
792  if (delta > deltaMax) delta = deltaMax;
793 
794  // initialize step history
795  setStates (2);
796  initStates ();
797  // initialise the history of states, setting them all to 'delta'
799 
800  // copy the initialised states to the 'deltas' array
802  // copy the deltas to all the circuits
803  setDelta ();
804  // set the initial corrector and predictor coefficients
807 
808  // initialize history of solution vectors (solutions)
809  for (int i = 0; i < 8; i++)
810  {
811  // solution contains the last sets of node voltages and branch
812  // currents at each of the last 8 'deltas'.
813  // Note for convenience the definition:
814  // #define SOL(state) (solution[(int) getState (sState, (state))])
815  // is provided and used elsewhere to update the solutions
817  setState (sState, (nr_double_t) i, i);
818  }
819 
820  // tell circuits about the transient analysis
821  circuit *c, * root = subnet->getRoot ();
822  for (c = root; c != NULL; c = (circuit *) c->getNext ())
823  initCircuitTR (c);
824  // also initialize created circuits
825  for (c = root; c != NULL; c = (circuit *) c->getPrev ())
826  initCircuitTR (c);
827 }
828 
829 // This function cleans up some memory used by the transient analysis.
831 {
832  // cleanup solutions
833  for (int i = 0; i < 8; i++)
834  {
835  delete solution[i];
836  solution[i] = NULL;
837  }
838  // cleanup history
839  if (tHistory)
840  {
841  delete tHistory;
842  tHistory = NULL;
843  }
844 }
845 
846 // The function initialize a single circuit.
848 {
849  c->initTR ();
850  c->initStates ();
852  c->setOrder (corrOrder);
854 }
855 
856 /* This function saves the results of a single solve() functionality
857  (for the given timestamp) into the output dataset. */
858 void trsolver::saveAllResults (nr_double_t time)
859 {
860  qucs::vector * t;
861  // add current frequency to the dependency of the output dataset
862  if ((t = data->findDependency ("time")) == NULL)
863  {
864  t = new qucs::vector ("time");
865  data->addDependency (t);
866  }
867  if (runs == 1) t->add (time);
868  saveResults ("Vt", "It", 0, t);
869 }
870 
871 /* This function is meant to adapt the current time-step the transient
872  analysis advanced. For the computation of the new time-step the
873  truncation error depending on the integration method is used. */
874 nr_double_t trsolver::checkDelta (void)
875 {
876  nr_double_t LTEreltol = getPropertyDouble ("LTEreltol");
877  nr_double_t LTEabstol = getPropertyDouble ("LTEabstol");
878  nr_double_t LTEfactor = getPropertyDouble ("LTEfactor");
879  nr_double_t dif, rel, tol, lte, q, n = NR_MAX;
880  int N = countNodes ();
881  int M = countVoltageSources ();
882 
883  // cec = corrector error constant
884  nr_double_t cec = getCorrectorError (corrType, corrOrder);
885  // pec = predictor error constant
886  nr_double_t pec = getPredictorError (predType, predOrder);
887 
888  // go through each solution
889  for (int r = 0; r < N + M; r++)
890  {
891 
892  // skip real voltage sources
893  if (r >= N)
894  {
895  if (findVoltageSource(r - N)->isVSource ())
896  continue;
897  }
898 
899  dif = x->get (r) - SOL(0)->get (r);
900  if (std::isfinite (dif) && dif != 0)
901  {
902  // use Milne' estimate for the local truncation error
903  rel = MAX (fabs (x->get (r)), fabs (SOL(0)->get (r)));
904  tol = LTEreltol * rel + LTEabstol;
905  lte = LTEfactor * (cec / (pec - cec)) * dif;
906  q = delta * exp (log (fabs (tol / lte)) / (corrOrder + 1));
907  n = MIN (n, q);
908  }
909  }
910 #if STEPDEBUG
911  logprint (LOG_STATUS, "DEBUG: delta according to local truncation "
912  "error h = %.3e\n", (double) n);
913 #endif
914  delta = MIN ((n > 1.9 * delta) ? 2 * delta : delta, n);
915  return delta;
916 }
917 
918 // The function updates the integration coefficients.
919 void trsolver::updateCoefficients (nr_double_t delta)
920 {
921  setState (dState, delta);
925 }
926 
927 // properties
928 PROP_REQ [] =
929 {
930  {
931  "Type", PROP_STR, { PROP_NO_VAL, "lin" },
932  PROP_RNG_STR2 ("lin", "log")
933  },
934  { "Start", PROP_REAL, { 0, PROP_NO_STR }, PROP_POS_RANGE },
935  { "Stop", PROP_REAL, { 1e-3, PROP_NO_STR }, PROP_POS_RANGE },
936  { "Points", PROP_INT, { 10, PROP_NO_STR }, PROP_MIN_VAL (2) },
938 };
939 PROP_OPT [] =
940 {
941  {
942  "IntegrationMethod", PROP_STR, { PROP_NO_VAL, "Trapezoidal" },
943  PROP_RNG_STR4 ("Euler", "Trapezoidal", "Gear", "AdamsMoulton")
944  },
945  { "Order", PROP_INT, { 2, PROP_NO_STR }, PROP_RNGII (1, 6) },
946  { "InitialStep", PROP_REAL, { 1e-9, PROP_NO_STR }, PROP_POS_RANGE },
947  { "MinStep", PROP_REAL, { 1e-16, PROP_NO_STR }, PROP_POS_RANGE },
948  { "MaxStep", PROP_REAL, { 0, PROP_NO_STR }, PROP_POS_RANGE },
949  { "MaxIter", PROP_INT, { 150, PROP_NO_STR }, PROP_RNGII (2, 10000) },
950  { "abstol", PROP_REAL, { 1e-12, PROP_NO_STR }, PROP_RNG_X01I },
951  { "vntol", PROP_REAL, { 1e-6, PROP_NO_STR }, PROP_RNG_X01I },
952  { "reltol", PROP_REAL, { 1e-3, PROP_NO_STR }, PROP_RNG_X01I },
953  { "LTEabstol", PROP_REAL, { 1e-6, PROP_NO_STR }, PROP_RNG_X01I },
954  { "LTEreltol", PROP_REAL, { 1e-3, PROP_NO_STR }, PROP_RNG_X01I },
955  { "LTEfactor", PROP_REAL, { 1, PROP_NO_STR }, PROP_RNGII (1, 16) },
956  { "Temp", PROP_REAL, { 26.85, PROP_NO_STR }, PROP_MIN_VAL (K) },
957  { "Solver", PROP_STR, { PROP_NO_VAL, "CroutLU" }, PROP_RNG_SOL },
958  { "relaxTSR", PROP_STR, { PROP_NO_VAL, "no" }, PROP_RNG_YESNO },
959  { "initialDC", PROP_STR, { PROP_NO_VAL, "yes" }, PROP_RNG_YESNO },
961 };
962 struct define_t trsolver::anadef =
964 
965 } // namespace qucs
start
Definition: parse_zvr.y:126
int statIterations
Definition: trsolver.h:97
#define PROP_POS_RANGE
Definition: netdefs.h:129
nr_double_t getCorrectorError(int, int)
Definition: transient.cpp:429
exceptionstack estack
circuit * getRoot(void)
Definition: net.h:46
net * subnet
Definition: analysis.h:273
#define PROP_RNGII(f, t)
Definition: netdefs.h:138
#define N(n)
Definition: equation.cpp:58
int corrector(void)
Definition: trsolver.cpp:570
#define dState
Definition: trsolver.cpp:53
static void calcTR(trsolver *)
Definition: trsolver.cpp:734
void predictGear(void)
Definition: trsolver.cpp:549
#define PROP_DEF
Definition: netdefs.h:189
#define PROP_RNG_STR2(s1, s2)
Definition: netdefs.h:145
trsolver(char *)
Definition: trsolver.cpp:77
#define MODE_NONE
Definition: integrator.h:30
void fillState(int, state_type_t)
Definition: states.cpp:132
PROP_OPT[]
Definition: acsolver.cpp:232
int correctorType(char *, int &)
Definition: transient.cpp:331
void recallSolution(void)
Definition: nasolver.cpp:1282
nr_double_t getPropertyDouble(const char *)
Definition: object.cpp:176
nr_double_t delta
Definition: trsolver.h:78
void applyNodeset(bool nokeep=true)
Definition: nasolver.cpp:247
#define PROP_REAL
Definition: netdefs.h:174
void setStates(int n)
Definition: states.h:52
void appendHistory(int, nr_double_t)
Definition: circuit.cpp:906
t
Definition: parse_vcd.y:290
#define PROP_NO_PROP
Definition: netdefs.h:122
void calcCorrectorCoeff(int, int, nr_double_t *, nr_double_t *)
Definition: transient.cpp:55
void logprogressclear(int points)
Definition: logging.c:90
void setDelta(void)
Definition: trsolver.cpp:615
nr_double_t last(void)
Definition: history.cpp:113
#define K
Absolute 0 in centigrade.
Definition: constants.h:59
int getPropertyInteger(const char *)
Definition: object.cpp:198
int countNodes(void)
Definition: nasolver.cpp:912
int predMaxOrder
Definition: trsolver.h:86
#define PROP_NO_STR
Definition: netdefs.h:125
int getVoltageSource(void)
Definition: circuit.h:187
int corrMaxOrder
Definition: trsolver.h:85
#define catch_exception()
n
Definition: parse_citi.y:147
#define PROP_LINEAR
Definition: netdefs.h:120
tvector< nr_double_t > * solution[8]
Definition: trsolver.h:93
void append(nr_double_t)
Definition: history.cpp:62
r
Definition: parse_mdl.y:515
tmatrix< nr_type_t > * A
Definition: nasolver.h:125
object * getNext(void)
Definition: object.h:59
#define PROP_INT
Definition: netdefs.h:173
void restart(void)
Definition: trsolver.cpp:745
void setIntegrationMethod(circuit *, int)
Definition: transient.cpp:308
void initStates(void)
Definition: states.cpp:75
int solve(void)
placehoder for solution function
Definition: trsolver.cpp:175
void setDescription(const char *n)
Definition: nasolver.h:65
int statConvergence
Definition: trsolver.h:98
#define PROP_RNG_X01I
Definition: netdefs.h:137
int getSize(void)
Get the number of ports the circuit element has.
Definition: circuit.h:143
void initDC(void)
Definition: trsolver.cpp:756
nr_double_t getPredictorError(int, int)
Definition: transient.cpp:434
#define PROP_ACTION
Definition: netdefs.h:115
void deinitTR(void)
Definition: trsolver.cpp:830
#define try_running()
void storeSolution(void)
Definition: nasolver.cpp:1258
#define PROP_RNG_SOL
Definition: netdefs.h:163
i
Definition: parse_mdl.y:516
void adjustDelta(nr_double_t)
Definition: trsolver.cpp:624
template class for storing state variables.
Definition: states.h:38
int getVoltageSources(void)
Definition: circuit.cpp:602
int dcAnalysis(void)
Definition: trsolver.cpp:123
void calcPredictorCoeff(int, int, nr_double_t *, nr_double_t *)
Definition: transient.cpp:169
void predictBashford(void)
Definition: trsolver.cpp:507
#define pop_exception()
void saveAllResults(nr_double_t)
Definition: trsolver.cpp:858
stop
Definition: parse_zvr.y:127
void initHistory(nr_double_t)
Definition: trsolver.cpp:401
points
Definition: parse_zvr.y:129
void fillStates(void)
Definition: trsolver.cpp:596
void initSteps(void)
Definition: trsolver.cpp:116
void reset(void)
Definition: sweep.h:56
base class for qucs circuit elements.
Definition: circuit.h:92
double tol
void solve_pre(void)
Definition: nasolver.cpp:213
void nextState(void)
Definition: states.cpp:117
nr_double_t deltas[8]
Definition: trsolver.h:77
void fillSolution(tvector< nr_double_t > *)
Definition: trsolver.cpp:500
void saveResults(const char *, const char *, int, qucs::vector *f=NULL)
Definition: nasolver.cpp:1308
nr_double_t deltaOld
Definition: trsolver.h:81
int getSize(void)
Definition: sweep.h:47
void add(nr_complex_t)
Definition: vector.cpp:151
net * getNet(void)
Definition: circuit.h:173
void setCalculation(calculate_func_t f)
Definition: nasolver.h:69
void setMode(int)
Definition: trsolver.cpp:607
void adjustOrder(int reduce=0)
Definition: trsolver.cpp:694
nr_double_t checkDelta(void)
Definition: trsolver.cpp:874
The analysis class header file.
circuit * findVoltageSource(int)
Definition: nasolver.cpp:953
void nextStates(void)
Definition: trsolver.cpp:578
dataset * data
Definition: analysis.h:274
void updateCoefficients(nr_double_t)
Definition: trsolver.cpp:919
void print(const char *prefix=NULL)
int predictor(void)
Definition: trsolver.cpp:476
void setAge(nr_double_t a)
Definition: history.h:42
nr_double_t deltaMax
Definition: trsolver.h:79
#define PROP_RNG_STR4(s1, s2, s3, s4)
Definition: netdefs.h:149
int predictorType(int, int, int &)
Definition: transient.cpp:361
state_type_t getState(int, int n=0)
Definition: states.cpp:99
void setOrder(int o)
Definition: integrator.h:52
The circuit class header file.
int solve_nonlinear(void)
Definition: nasolver.cpp:449
bool initialDC
Definition: trsolver.h:101
virtual void initTR(void)
Definition: circuit.h:122
void saveHistory(circuit *)
Definition: trsolver.cpp:445
type
Definition: parse_vcd.y:164
#define PROP_MIN_VAL(k)
Definition: netdefs.h:133
int countVoltageSources(void)
Definition: nasolver.cpp:944
#define MIN(x, y)
Minimum of x and y.
Definition: constants.h:132
void logprogressbar(nr_double_t current, nr_double_t final, int points)
Definition: logging.c:63
#define CONV_LineSearch
Definition: nasolver.h:38
void predictEuler(void)
Definition: trsolver.cpp:531
void setCoefficients(nr_double_t *c)
Definition: integrator.h:56
void saveState(int, state_type_t *)
Definition: states.cpp:146
nr_double_t corrCoeff[8]
Definition: trsolver.h:76
nr_double_t deltaMin
Definition: trsolver.h:80
#define PROP_STR
Definition: netdefs.h:175
const char * getDescription(void)
Definition: nasolver.h:66
object * getPrev(void)
Definition: object.h:61
sweep * createSweep(const char *)
create a named sweep object
Definition: analysis.cpp:107
nr_double_t current
Definition: trsolver.h:94
void initCircuitTR(circuit *)
Definition: trsolver.cpp:847
nr_complex_t exp(const nr_complex_t z)
Compute complex exponential.
Definition: complex.cpp:205
#define LOG_ERROR
Definition: logging.h:28
char * getName(void)
Definition: object.cpp:84
void self(void)
Definition: history.h:54
nr_double_t stepDelta
Definition: trsolver.h:82
void initTR(void)
Definition: trsolver.cpp:767
void drop(void)
Definition: history.cpp:124
#define PROP_NO_VAL
Definition: netdefs.h:124
void solve_post(void)
Definition: nasolver.cpp:205
#define LOG_STATUS
Definition: logging.h:29
#define SOL(state)
Definition: trsolver.cpp:57
nr_double_t next(void)
Definition: sweep.cpp:140
PROP_REQ[]
Definition: acsolver.cpp:229
void saveSolution(void)
Definition: nasolver.cpp:1250
#define CONV_SteepestDescent
Definition: nasolver.h:39
#define PROP_RNG_YESNO
Definition: netdefs.h:158
char * getPropertyString(const char *)
Definition: object.cpp:159
void setState(int, state_type_t, int n=0)
Definition: states.cpp:109
sweep * swp
Definition: trsolver.h:74
#define M(con)
Definition: evaluate.cpp:64
void logprint(int level, const char *format,...)
Definition: logging.c:37
nr_complex_t log(const nr_complex_t z)
Compute principal value of natural logarithm of z.
Definition: complex.cpp:215
tvector< nr_type_t > * x
Definition: nasolver.h:122
int statRejected
Definition: trsolver.h:96
int findAssignedNode(circuit *, int)
Definition: nasolver.cpp:927
void(* calculate_func_t)(nasolver< nr_type_t > *)
Definition: nasolver.h:68
#define sState
Definition: trsolver.cpp:54
static void calcDC(trsolver *)
Definition: trsolver.cpp:723
nr_double_t predCoeff[8]
Definition: trsolver.h:75
history * tHistory
Definition: trsolver.h:99
#define MODE_INIT
Definition: integrator.h:31
#define PROP_NO_SUBSTRATE
Definition: netdefs.h:118
void updateHistory(nr_double_t)
Definition: trsolver.cpp:428
#define MAX(x, y)
Maximum of x and y.
Definition: constants.h:127