Cantera  2.2.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BEulerInt.h
Go to the documentation of this file.
1 /**
2  * @file BEulerInt.h
3  */
4 
5 /*
6  * Copyright 2004 Sandia Corporation. Under the terms of Contract
7  * DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government
8  * retains certain rights in this software.
9  * See file License.txt for licensing information.
10  */
11 #ifndef CT_BEULERINT_H
12 #define CT_BEULERINT_H
13 
17 
18 #define OPT_SIZE 10
19 
20 #define SUCCESS 0
21 #define FAILURE 1
22 
23 #define STEADY 0
24 #define TRANSIENT 1
25 
26 namespace Cantera
27 {
28 
29 enum BEulerMethodType {
30  BEulerFixedStep,
31  BEulerVarStep
32 };
33 
34 /**
35  * Exception class thrown when a BEuler error is encountered.
36  */
37 class BEulerErr : public CanteraError
38 {
39 public:
40  /**
41  * Exception thrown when a BEuler error is encountered. We just call the
42  * Cantera Error handler in the initialization list.
43  */
44  explicit BEulerErr(const std::string& msg);
45 };
46 
47 
48 #define BEULER_JAC_ANAL 2
49 #define BEULER_JAC_NUM 1
50 
51 /*!
52  * Wrapper class for 'beuler' integrator
53  * We derive the class from the class Integrator
54  * @deprecated Unused. To be removed after Cantera 2.2.
55  */
56 class BEulerInt : public Integrator
57 {
58 public:
59  /*!
60  * Constructor. Default settings: dense Jacobian, no user-supplied
61  * Jacobian function, Newton iteration.
62  */
63  BEulerInt();
64 
65  virtual ~BEulerInt();
66 
67  virtual void setTolerances(double reltol, size_t n, double* abstol);
68  virtual void setTolerances(double reltol, double abstol);
69  virtual void setProblemType(int probtype);
70 
71  //! Find the initial conditions for y and ydot.
72  virtual void initializeRJE(double t0, ResidJacEval& func);
73  virtual void reinitializeRJE(double t0, ResidJacEval& func);
74  virtual double integrateRJE(double tout, double tinit = 0.0);
75 
76  // This routine advances the calculations one step using a predictor
77  // corrector approach. We use an implicit algorithm here.
78  virtual doublereal step(double tout);
79 
80  //! Set the solution weights. This is a very important routine as it affects
81  //! quite a few operations involving convergence.
82  virtual void setSolnWeights();
83 
84  virtual double& solution(size_t k) {
85  return m_y_n[k];
86  }
87  double* solution() {
88  return &m_y_n[0];
89  }
90  int nEquations() const {
91  return m_neq;
92  }
93 
94  //! Return the total number of function evaluations
95  virtual int nEvals() const;
96  virtual void setMethodBEMT(BEulerMethodType t);
97  virtual void setIterator(IterType t);
98  virtual void setMaxStep(double hmax);
99  virtual void setMaxNumTimeSteps(int);
100  virtual void setNumInitialConstantDeltaTSteps(int);
101 
102  void print_solnDelta_norm_contrib(const double* const soln0,
103  const char* const s0,
104  const double* const soln1,
105  const char* const s1,
106  const char* const title,
107  const double* const y0,
108  const double* const y1,
109  double damp,
110  int num_entries);
111 
112  //! This routine controls when the solution is printed
113  /*!
114  * @param printSolnStepInterval If greater than 0, then the soln is
115  * printed every printSolnStepInterval steps.
116  * @param printSolnNumberToTout The solution is printed at regular
117  * invervals a total of "printSolnNumberToTout" times.
118  * @param printSolnFirstSteps The solution is printed out the first
119  * "printSolnFirstSteps" steps. After these steps the
120  * other parameters determine the printing. default = 0
121  * @param dumpJacobians Dump Jacobians to disk.
122  */
123  virtual void setPrintSolnOptions(int printSolnStepInterval,
124  int printSolnNumberToTout,
125  int printSolnFirstSteps = 0,
126  bool dumpJacobians = false);
127 
128  //! Set the options for the nonlinear method
129  /*!
130  * Defaults are set in the .h file. These are the defaults:
131  * min_newt_its = 0
132  * matrixConditioning = false
133  * colScaling = false
134  * rowScaling = true
135  */
136  void setNonLinOptions(int min_newt_its = 0,
137  bool matrixConditioning = false,
138  bool colScaling = false,
139  bool rowScaling = true);
140  virtual void setPrintFlag(int print_flag);
141 
142  //! Set the column scaling vector at the current time
143  virtual void setColumnScales();
144 
145  /**
146  * Calculate the solution error norm. if printLargest is true, then a table
147  * of the largest values is printed to standard output.
148  */
149  virtual double soln_error_norm(const double* const,
150  bool printLargest = false);
151  virtual void setInitialTimeStep(double delta_t);
152 
153  /*!
154  * Function called by to evaluate the Jacobian matrix and the current
155  * residual at the current time step.
156  * @param J = Jacobian matrix to be filled in
157  * @param f = Right hand side. This routine returns the current
158  * value of the RHS (output), so that it does
159  * not have to be computed again.
160  */
161  void beuler_jac(GeneralMatrix& J, double* const f,
162  double, double, double* const, double* const, int);
163 
164 protected:
165  //! Internal routine that sets up the fixed length storage based on
166  //! the size of the problem to solve.
167  void internalMalloc();
168 
169  /*!
170  * Function to calculate the predicted solution vector, m_y_pred_n for the
171  * (n+1)th time step. This routine can be used by a first order - forward
172  * Euler / backward Euler predictor / corrector method or for a second order
173  * Adams-Bashforth / Trapezoidal Rule predictor / corrector method. See
174  * Nachos documentation Sand86-1816 and Gresho, Lee, Sani LLNL report UCRL -
175  * 83282 for more information.
176  *
177  * on input:
178  *
179  * N - number of unknowns
180  * order - indicates order of method
181  * = 1 -> first order forward Euler/backward Euler
182  * predictor/corrector
183  * = 2 -> second order Adams-Bashforth/Trapezoidal Rule
184  * predictor/corrector
185  *
186  * delta_t_n - magnitude of time step at time n (i.e., = t_n+1 - t_n)
187  * delta_t_nm1 - magnitude of time step at time n - 1 (i.e., = t_n - t_n-1)
188  * y_n[] - solution vector at time n
189  * y_dot_n[] - acceleration vector from the predictor at time n
190  * y_dot_nm1[] - acceleration vector from the predictor at time n - 1
191  *
192  * on output:
193  *
194  * m_y_pred_n[] - predicted solution vector at time n + 1
195  */
196  void calc_y_pred(int);
197 
198  /*!
199  * Function to calculate the acceleration vector ydot for the first or
200  * second order predictor/corrector time integrator. This routine can be
201  * called by a first order - forward Euler / backward Euler predictor /
202  * corrector or for a second order Adams - Bashforth / Trapezoidal Rule
203  * predictor / corrector. See Nachos documentation Sand86-1816 and Gresho,
204  * Lee, Sani LLNL report UCRL - 83282 for more information.
205  *
206  * on input:
207  *
208  * N - number of local unknowns on the processor
209  * This is equal to internal plus border unknowns.
210  * order - indicates order of method
211  * = 1 -> first order forward Euler/backward Euler
212  * predictor/corrector
213  * = 2 -> second order Adams-Bashforth/Trapezoidal Rule
214  * predictor/corrector
215  *
216  * delta_t_n - Magnitude of the current time step at time n
217  * (i.e., = t_n - t_n-1)
218  * y_curr[] - Current Solution vector at time n
219  * y_nm1[] - Solution vector at time n-1
220  * ydot_nm1[] - Acceleration vector at time n-1
221  *
222  * on output:
223  *
224  * ydot_curr[] - Current acceleration vector at time n
225  *
226  * Note we use the current attribute to denote the possibility that
227  * y_curr[] may not be equal to m_y_n[] during the nonlinear solve
228  * because we may be using a look-ahead scheme.
229  */
230  void calc_ydot(int, double*, double*);
231 
232  /*!
233  * Calculates the time step truncation error estimate from a very simple
234  * formula based on Gresho et al. This routine can be called for a first
235  * order - forward Euler/backward Euler predictor/ corrector and for a
236  * second order Adams- Bashforth/Trapezoidal Rule predictor/corrector. See
237  * Nachos documentation Sand86-1816 and Gresho, Lee, LLNL report UCRL -
238  * 83282 for more information.
239  *
240  * on input:
241  *
242  * abs_error - Generic absolute error tolerance
243  * rel_error - Generic relative error tolerance
244  * x_coor[] - Solution vector from the implicit corrector
245  * x_pred_n[] - Solution vector from the explicit predictor
246  *
247  * on output:
248  *
249  * delta_t_n - Magnitude of next time step at time t_n+1
250  * delta_t_nm1 - Magnitude of previous time step at time t_n
251  */
252  double time_error_norm();
253 
254  /*!
255  * Time step control function for the selection of the time step size based on
256  * a desired accuracy of time integration and on an estimate of the relative
257  * error of the time integration process. This routine can be called for a
258  * first order - forward Euler/backward Euler predictor/ corrector and for a
259  * second order Adams- Bashforth/Trapezoidal Rule predictor/corrector. See
260  * Nachos documentation Sand86-1816 and Gresho, Lee, Sani LLNL report UCRL -
261  * 83282 for more information.
262  *
263  * on input:
264  *
265  * order - indicates order of method
266  * = 1 -> first order forward Euler/backward Euler
267  * predictor/corrector
268  * = 2 -> second order forward Adams-Bashforth/Trapezoidal
269  * rule predictor/corrector
270  *
271  * delta_t_n - Magnitude of time step at time t_n
272  * delta_t_nm1 - Magnitude of time step at time t_n-1
273  * rel_error - Generic relative error tolerance
274  * time_error_factor - Estimated value of the time step truncation error
275  * factor. This value is a ratio of the computed
276  * error norms. The premultiplying constants
277  * and the power are not yet applied to normalize the
278  * predictor/corrector ratio. (see output value)
279  *
280  * on output:
281  *
282  * return - delta_t for the next time step
283  * If delta_t is negative, then the current time step is
284  * rejected because the time-step truncation error is
285  * too large. The return value will contain the negative
286  * of the recommended next time step.
287  *
288  * time_error_factor - This output value is normalized so that
289  * values greater than one indicate the current time
290  * integration error is greater than the user
291  * specified magnitude.
292  */
293  double time_step_control(int m_order, double time_error_factor);
294 
295  //! Solve a nonlinear system
296  /*!
297  * Find the solution to F(X, xprime) = 0 by damped Newton iteration. On
298  * entry, y_comm[] contains an initial estimate of the solution and
299  * ydot_comm[] contains an estimate of the derivative.
300  * On successful return, y_comm[] contains the converged solution
301  * and ydot_comm[] contains the derivative
302  *
303  * @param y_comm[] Contains the input solution. On output y_comm[] contains
304  * the converged solution
305  * @param ydot_comm Contains the input derivative solution. On output
306  * y_comm[] contains the converged derivative solution
307  * @param CJ Inverse of the time step
308  * @param time_curr Current value of the time
309  * @param jac Jacobian
310  * @param num_newt_its number of Newton iterations
311  * @param num_linear_solves number of linear solves
312  * @param num_backtracks number of backtracs
313  * @param loglevel Log level
314  */
315  int solve_nonlinear_problem(double* const y_comm,
316  double* const ydot_comm, double CJ,
317  double time_curr,
318  GeneralMatrix& jac,
319  int& num_newt_its,
320  int& num_linear_solves,
321  int& num_backtracks,
322  int loglevel);
323 
324  /**
325  * Compute the undamped Newton step. The residual function is
326  * evaluated at the current time, t_n, at the current values of the
327  * solution vector, m_y_n, and the solution time derivative, m_ydot_n,
328  * but the Jacobian is not recomputed.
329  */
330  void doNewtonSolve(double, double*, double*, double*,
331  GeneralMatrix&, int);
332 
333 
334  //! Bound the Newton step while relaxing the solution
335  /*!
336  * Return the factor by which the undamped Newton step 'step0'
337  * must be multiplied in order to keep all solution components in
338  * all domains between their specified lower and upper bounds.
339  * Other bounds may be applied here as well.
340  *
341  * Currently the bounds are hard coded into this routine:
342  *
343  * Minimum value for all variables: - 0.01 * m_ewt[i]
344  * Maximum value = none.
345  *
346  * Thus, this means that all solution components are expected
347  * to be numerical greater than zero in the limit of time step
348  * truncation errors going to zero.
349  *
350  * Delta bounds: The idea behind these is that the Jacobian
351  * couldn't possibly be representative if the
352  * variable is changed by a lot. (true for
353  * nonlinear systems, false for linear systems)
354  * Maximum increase in variable in any one Newton iteration: factor of 2
355  * Maximum decrease in variable in any one Newton iteration: factor of 5
356  *
357  * @param y Current value of the solution
358  * @param step0 Current raw step change in y[]
359  * @param loglevel Log level. This routine produces output if loglevel
360  * is greater than one
361  *
362  * @return Returns the damping coefficient
363  */
364  double boundStep(const double* const y, const double* const step0, int loglevel);
365 
366  /*!
367  * On entry, step0 must contain an undamped Newton step for the
368  * solution x0. This method attempts to find a damping coefficient
369  * such that the next undamped step would have a norm smaller than
370  * that of step0. If successful, the new solution after taking the
371  * damped step is returned in y1, and the undamped step at y1 is
372  * returned in step1.
373  */
374  int dampStep(double, const double*, const double*,
375  const double*, double*, double*,
376  double*, double&, GeneralMatrix&, int&, bool, int&);
377 
378  //! Compute Residual Weights
379  void computeResidWts(GeneralMatrix& jac);
380 
381  //! Filter a new step
382  double filterNewStep(double, double*, double*);
383 
384  //! Get the next time to print out
385  double getPrintTime(double time_current);
386 
387  /********************** Member data ***************************/
388  /*********************
389  * METHOD FLAGS
390  *********************/
391 
392  //! IterType is used to specify how the nonlinear equations are
393  //! to be relaxed at each time step.
395  /**
396  * MethodType is used to specify how the time step is to be
397  * chosen. Currently, there are two choices, one is a fixed
398  * step method while the other is based on a predictor-corrector
399  * algorithm and a time-step truncation error tolerance.
400  */
401  BEulerMethodType m_method;
402  /**
403  * m_jacFormMethod determines how a matrix is formed.
404  */
406  /**
407  * m_rowScaling is a boolean. If true then row sum scaling
408  * of the Jacobian matrix is carried out when solving the
409  * linear systems.
410  */
412  /**
413  * m_colScaling is a boolean. If true, then column scaling
414  * is performed on each solution of the linear system.
415  */
417  /**
418  * m_matrixConditioning is a boolean. If true, then the
419  * Jacobian and every RHS is multiplied by the inverse
420  * of a matrix that is suppose to reduce the condition
421  * number of the matrix. This is done before row scaling.
422  */
424  /**
425  * If m_itol =1 then each component has an individual
426  * value of atol. If m_itol = 0, the all atols are equal.
427  */
428  int m_itol;
429 
430  //! Relative time truncation error tolerances
431  double m_reltol;
432 
433  /**
434  * Absolute time truncation error tolerances, when uniform
435  * for all variables.
436  */
437  double m_abstols;
438  /**
439  * Vector of absolute time truncation error tolerance
440  * when not uniform for all variables.
441  */
443 
444  //! Error Weights. This is a surprisingly important quantity.
446 
447  //! Maximum step size
448  double m_hmax;
449 
450  //! Maximum integration order
451  int m_maxord;
452 
453  //! Current integration order
454  int m_order;
455 
456  //! Time step number
458  int m_time_step_attempts;
459 
460  //! Max time steps allowed
462  /**
463  * Number of initial time steps to take where the time truncation error
464  * tolerances are not checked. Instead the delta T is uniform
465  */
467 
468  //! Failure Counter -> keeps track of the number of consecutive failures
470 
471  //! Minimum Number of Newton Iterations per nonlinear step. default = 0
473  /************************
474  * PRINTING OPTIONS
475  ************************/
476  /**
477  * Step Interval at which to print out the solution
478  * default = 1;
479  * If set to zero, there is no printout
480  */
482  /**
483  * Number of evenly spaced printouts of the solution
484  * If zero, there is no printout from this option
485  * default 1
486  * If set to zero there is no printout.
487  */
489 
490  //! Number of initial steps that the solution is printed out. default = 0
492 
493  //! Dump Jacobians to disk. default false
495 
496  /*********************
497  * INTERNAL SOLUTION VALUES
498  *********************/
499 
500  //! Number of equations in the ode integrator
501  int m_neq;
502  vector_fp m_y_n;
503  vector_fp m_y_nm1;
504  vector_fp m_y_pred_n;
505  vector_fp m_ydot_n;
506  vector_fp m_ydot_nm1;
507  /************************
508  * TIME VARIABLES
509  ************************/
510 
511  //! Initial time at the start of the integration
512  double m_t0;
513 
514  //! Final time
515  double m_time_final;
516 
517  double time_n;
518  double time_nm1;
519  double time_nm2;
520  double delta_t_n;
521  double delta_t_nm1;
522  double delta_t_nm2;
523  double delta_t_np1;
524 
525  //! Maximum permissible time step
526  double delta_t_max;
527 
528  vector_fp m_resid;
529  vector_fp m_residWts;
530  vector_fp m_wksp;
531  ResidJacEval* m_func;
532  vector_fp m_rowScales;
533  vector_fp m_colScales;
534 
535  //! Pointer to the Jacobian representing the time dependent problem
537 
538  /**
539  * Determines the level of printing for each time
540  * step.
541  * 0 -> absolutely nothing is printed for
542  * a single time step.
543  * 1 -> One line summary per time step
544  * 2 -> short description, points of interest
545  * 3 -> Lots printed per time step (default)
546  */
548  /***************************************************************************
549  * COUNTERS OF VARIOUS KINDS
550  ***************************************************************************/
551 
552  //! Number of function evaluations
553  int m_nfe;
554 
555  /**
556  * Number of Jacobian Evaluations and
557  * factorization steps (they are the same)
558  */
560 
561  //! Number of total Newton iterations
563 
564  //! Total number of linear iterations
566 
567  //! Total number of convergence failures.
569 
570  //! Total Number of time truncation error failures
572 
573  int num_failures;
574 };
575 
576 } // namespace
577 
578 #endif // CT_BEULER
virtual void setPrintSolnOptions(int printSolnStepInterval, int printSolnNumberToTout, int printSolnFirstSteps=0, bool dumpJacobians=false)
This routine controls when the solution is printed.
Definition: BEulerInt.cpp:124
double m_time_final
Final time.
Definition: BEulerInt.h:515
int m_numTotalTruncFails
Total Number of time truncation error failures.
Definition: BEulerInt.h:571
double time_error_norm()
Definition: BEulerInt.cpp:639
virtual void setProblemType(int probtype)
Set the problem type.
Definition: BEulerInt.cpp:99
int m_printSolnNumberToTout
Number of evenly spaced printouts of the solution If zero, there is no printout from this option defa...
Definition: BEulerInt.h:488
double getPrintTime(double time_current)
Get the next time to print out.
Definition: BEulerInt.cpp:210
virtual void setIterator(IterType t)
Set the linear iterator.
Definition: BEulerInt.cpp:135
int m_printSolnStepInterval
Step Interval at which to print out the solution default = 1; If set to zero, there is no printout...
Definition: BEulerInt.h:481
int m_numTotalNewtIts
Number of total Newton iterations.
Definition: BEulerInt.h:562
int m_neq
Number of equations in the ode integrator.
Definition: BEulerInt.h:501
virtual double soln_error_norm(const double *const, bool printLargest=false)
Calculate the solution error norm.
Definition: BEulerInt.cpp:1295
void doNewtonSolve(double, double *, double *, double *, GeneralMatrix &, int)
Compute the undamped Newton step.
Definition: BEulerInt.cpp:1349
int m_failure_counter
Failure Counter -> keeps track of the number of consecutive failures.
Definition: BEulerInt.h:469
BEulerMethodType m_method
MethodType is used to specify how the time step is to be chosen.
Definition: BEulerInt.h:401
bool m_rowScaling
m_rowScaling is a boolean.
Definition: BEulerInt.h:411
double m_reltol
Relative time truncation error tolerances.
Definition: BEulerInt.h:431
void beuler_jac(GeneralMatrix &J, double *const f, double, double, double *const, double *const, int)
Definition: BEulerInt.cpp:483
double m_t0
Initial time at the start of the integration.
Definition: BEulerInt.h:512
Generic matrix.
Definition: GeneralMatrix.h:24
bool m_matrixConditioning
m_matrixConditioning is a boolean.
Definition: BEulerInt.h:423
int m_min_newt_its
Minimum Number of Newton Iterations per nonlinear step. default = 0.
Definition: BEulerInt.h:472
int m_print_flag
Determines the level of printing for each time step.
Definition: BEulerInt.h:547
virtual double & solution(size_t k)
The current value of the solution of equation k.
Definition: BEulerInt.h:84
vector_fp m_abstol
Vector of absolute time truncation error tolerance when not uniform for all variables.
Definition: BEulerInt.h:442
BEulerErr(const std::string &msg)
Exception thrown when a BEuler error is encountered.
Definition: BEulerInt.cpp:21
IterType m_iter
IterType is used to specify how the nonlinear equations are to be relaxed at each time step...
Definition: BEulerInt.h:394
Wrappers for the function evaluators for Nonlinear solvers and Time steppers.
Definition: ResidJacEval.h:51
double delta_t_max
Maximum permissible time step.
Definition: BEulerInt.h:526
Exception class thrown when a BEuler error is encountered.
Definition: BEulerInt.h:37
int solve_nonlinear_problem(double *const y_comm, double *const ydot_comm, double CJ, double time_curr, GeneralMatrix &jac, int &num_newt_its, int &num_linear_solves, int &num_backtracks, int loglevel)
Solve a nonlinear system.
Definition: BEulerInt.cpp:1731
int m_numTotalLinearSolves
Total number of linear iterations.
Definition: BEulerInt.h:565
GeneralMatrix * tdjac_ptr
Pointer to the Jacobian representing the time dependent problem.
Definition: BEulerInt.h:536
virtual void setSolnWeights()
Set the solution weights.
Definition: BEulerInt.cpp:250
int m_printSolnFirstSteps
Number of initial steps that the solution is printed out. default = 0.
Definition: BEulerInt.h:491
virtual void setColumnScales()
Set the column scaling vector at the current time.
Definition: BEulerInt.cpp:266
int m_numInitialConstantDeltaTSteps
Number of initial time steps to take where the time truncation error tolerances are not checked...
Definition: BEulerInt.h:466
void calc_y_pred(int)
Definition: BEulerInt.cpp:590
double * solution()
The current value of the solution of the system of equations.
Definition: BEulerInt.h:87
int m_order
Current integration order.
Definition: BEulerInt.h:454
int m_nfe
Number of function evaluations.
Definition: BEulerInt.h:553
Dense, Square (not sparse) matrices.
virtual void initializeRJE(double t0, ResidJacEval &func)
Find the initial conditions for y and ydot.
Definition: BEulerInt.cpp:165
virtual void setTolerances(double reltol, size_t n, double *abstol)
Set or reset the number of equations.
Definition: BEulerInt.cpp:78
Abstract base class for ODE system integrators.
Definition: Integrator.h:52
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:99
int m_maxord
Maximum integration order.
Definition: BEulerInt.h:451
double m_abstols
Absolute time truncation error tolerances, when uniform for all variables.
Definition: BEulerInt.h:437
double time_step_control(int m_order, double time_error_factor)
Definition: BEulerInt.cpp:693
int dampStep(double, const double *, const double *, const double *, double *, double *, double *, double &, GeneralMatrix &, int &, bool, int &)
Definition: BEulerInt.cpp:1582
int m_max_time_step_attempts
Max time steps allowed.
Definition: BEulerInt.h:461
int nEquations() const
The number of equations.
Definition: BEulerInt.h:90
double filterNewStep(double, double *, double *)
Filter a new step.
Definition: BEulerInt.cpp:295
int m_jacFormMethod
m_jacFormMethod determines how a matrix is formed.
Definition: BEulerInt.h:405
int m_numTotalConvFails
Total number of convergence failures.
Definition: BEulerInt.h:568
virtual int nEvals() const
Return the total number of function evaluations.
Definition: BEulerInt.cpp:225
double boundStep(const double *const y, const double *const step0, int loglevel)
Bound the Newton step while relaxing the solution.
Definition: BEulerInt.cpp:1515
bool m_colScaling
m_colScaling is a boolean.
Definition: BEulerInt.h:416
IterType
Specifies the method used for iteration.
Definition: Integrator.h:40
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:157
void internalMalloc()
Internal routine that sets up the fixed length storage based on the size of the problem to solve...
Definition: BEulerInt.cpp:230
bool m_dumpJacobians
Dump Jacobians to disk. default false.
Definition: BEulerInt.h:494
virtual doublereal step(double tout)
Integrate the system of equations.
Definition: BEulerInt.cpp:888
Declarations for the class GeneralMatrix which is a virtual base class for matrices handled by solver...
int m_nJacEval
Number of Jacobian Evaluations and factorization steps (they are the same)
Definition: BEulerInt.h:559
void setNonLinOptions(int min_newt_its=0, bool matrixConditioning=false, bool colScaling=false, bool rowScaling=true)
Set the options for the nonlinear method.
Definition: BEulerInt.cpp:140
void calc_ydot(int, double *, double *)
Definition: BEulerInt.cpp:618
double m_hmax
Maximum step size.
Definition: BEulerInt.h:448
int m_itol
If m_itol =1 then each component has an individual value of atol.
Definition: BEulerInt.h:428
vector_fp m_ewt
Error Weights. This is a surprisingly important quantity.
Definition: BEulerInt.h:445
int m_time_step_num
Time step number.
Definition: BEulerInt.h:457
void computeResidWts(GeneralMatrix &jac)
Compute Residual Weights.
Definition: BEulerInt.cpp:271