Purpose
Interface for using a common entry point, DSblock compatible for defining Differential Algebraic Equations using several packages. The equations follow the form: dx/dt = f(x(t), (t), p, t) y(t) = g(x(t), (t), p, t) The user must define only the subroutines ODEDER and ODEOUT and the Jacobians (JACFX, JACFU, JACFP) if used, and the interface adapts the structure to fit all the codesSpecification
SUBROUTINE ODESolver(ISOLVER, CODEDER_, CODEOUT_, CJACFX_, $ CJACFU_, CJACFP_, $ NX, NY, NU, TINI, TOUT, X, U, Y, $ IPAR, DPAR, RTOL, ATOL, $ IWORK, LIWORK, DWORK, LDWORK, $ IWARN, INFO) .. Scalar Arguments .. DOUBLE PRECISION RTOL, TINI, TOUT INTEGER ISOLVER, IWARN, INFO, NX, NY, NU, $ LDWORK, LIWORK CHARACTER*9 CODEDER_, CODEOUT_, CJACFX_, CJACFU_, CJACFP_ .. Array Arguments .. DOUBLE PRECISION ATOL(*), DWORK(LDWORK), DPAR(*), $ X(NX), Y(NY), U(NU) INTEGER IWORK(LIWORK), IPAR(*)Arguments
Mode Parameters
ISOLVER (input) INTEGER Indicates the nonlinear solver package to be used: = 1: LSODE, = 2: LSODA, = 3: LSODES, = 4: RADAU5, = 5: DASSL, = 6: DASPK, = 7: DGELDA.Input/Output Parameters
ODEDER (input) EXTERNAL Evaluates the right hand side f of the ODE. ODEOUT (input) EXTERNAL Evaluates the output signals function g. JACFX (input) EXTERNAL Evaluates the jacobian matrix with respect to X. JACFU (input) EXTERNAL Evaluates the jacobian matrix with respect to U. JACFP (input) EXTERNAL Evaluates the jacobian matrix with respect to P. NX (input) INTEGER Dimension of the state vector. NY (input) INTEGER Dimension of the output vector. NU (input) INTEGER Dimension of the input vector.
TINI (input) DOUBLE PRECISION Initial value of time. TOUT (input) DOUBLE PRECISION Final value of time. X (input/output) DOUBLE PRECISION array, dimension (NX) On entry, array containing the initial state variables. On exit, it has the last value of the state variables. U (input) DOUBLE PRECISION array, dimension (NU) Array containing the input initial values. Y (input/output) DOUBLE PRECISION array, dimension (NY) On entry, array containing the initial values of Y. On exit, it has the results of the system. IPAR (input/output) INTEGER array, dimension (230) INPUT: 1..15 General 16..25 ODEPACK 26..35 RADAU5 36..50 DASSL/PK 51..60 GELDA 61..100 Reserved OUTPUT: 101..110 General 111..125 ODEPACK 126..135 RADAU5 136..145 DASSL/PK 146..155 GELDA 156..200 Reserved Any Mode: 201.. User Available Common integer parameters for SOLVERS: IPAR(1), Tolerance mode 0 : both RTOL and ATOL are scalars 1 : RTOL is scalar and ATOL is vector 2 : both RTOL and ATOL are vectors IPAR(2), Compute Output Values, must be 1 IPAR(3), mf, Method flag 0 : No jacobian used (non-stiff method). 1 : User supplied full jacobian (stiff). 2 : User supplied banded jacobian (stiff). 3 : User supplied sparse jacobian (stiff). 10 : internally generated full jacobian (stiff). 11 : internally generated banded jacobian (stiff). 12 : internally generated sparse jacobian (stiff). IPAR(5), Maximum number of steps allowed during one call to the solver. IPAR(6), ml, lower half-bandwithds of the banded jacobian, excluding tne main diagonal. IPAR(7), mu, upper half-bandwithds of the banded jacobian, excluding tne main diagonal. IPAR(8), Flag to generate extra printing at method switches: 0 means no extra printing 1 for minimal printing 2 for full printing IPAR(101) = Number of steps taken for the problem. IPAR(102) = Number of f evaluations. IPAR(103) = Number of jacobian evaluations. Common parameters for ODEPACK, DASSL and DASPK solvers: IPAR(111) = The method order last used(successfully). IPAR(112) = The order to be attempted on the next step. Common parameters for ODEPACK solver: IPAR(16), Status Flag IPAR(17), Optional inputs IPAR(18), Maximum number of messages printed, default value is 10. IPAR(113) = Index of the component of largest in the weighted local error vector ( e(i)/ewt(i) ). IPAR(114) = Length of rwork actually required. IPAR(115) = Length of iwork actually required. - LSODE and LSODES IPAR(19), Maximum order to be allowed. 12 if meth = 1 5 if meth = 2 If exceds the default value, it will be reduced to the default value. - LSODES IPAR(118), Number of nonzero elements in the jacobian matrix, including the diagonal (miter = 1 or 2). IPAR(119), Number of groups of column indices, used in difference quotient jacobian aproximations if miter = 2. IPAR(120), Number of sparse LU decompositions. IPAR(121), Base address in rwork of the history array. IPAR(122), Base address of the structure descriptor array ian. IPAR(123), Base address of the structure descriptor array jan. IPAR(124), Number of nonzero elements in the strict lower triangle of the LU factorization. IPAR(125), Number of nonzero elements in the strict upper triangle of the LU factorization. - LSODA IPAR(22), Maximum order to be allowed for the nonstiff method, default value is 12. If exceds the default value, it will be reduced to the default value. IPAR(23), Maximum order to be allowed for the stiff method, default value is 5. If exceds the default value, it will be reduced to the default value. IPAR(116), Method indicator for the last successful step 1 adams (nonstiff) 2 bdf (stiff) IPAR(117), Current method indicator 1 adams (nonstiff) 2 bdf (stiff) Parameters for RADAU5 solver: IPAR(26) Transforms the Jacobian matrix to Hessenberg form. IPAR(27) Maximum number of Newton iterations. IPAR(28) Starting values for Newton's method if 0 then is taken the extrapolated collocation solution if not equal 0 zero values are used. IPAR(29) Dimension of the index 1 variables. IPAR(30) Dimension of the index 2 variables. IPAR(31) Dimension of the index 3 variables. IPAR(32) Switch for step size strategy 0,1 Mod. Predictive controller(Gustafsson) 2 Classical step size control IPAR(33) Value of M1. IPAR(34) Value of M2. IPAR(126), Number of accepted steps. IPAR(127), Number of rejected steps. IPAR(128), Number of LU-Decompositions of both matrices IPAR(129), Number of forward-backward substitutions, of both systems. Common parameters for DASSL and DASPK solvers: IPAR(36), this parameter enables the code to initialize itself. Must set to 0 to indicate the start of every new problem. 0: Yes. (On each new problem) 1: No. (Allows 500 new steps) IPAR(37), code solve the problem without invoking any special non negativity contraints: 0: Yes 1: To have constraint checking only in the initial condition calculation. 2: To enforce nonnegativity in X during the integration. 3: To enforce both options 1 and 2. IPAR(38), Solver try to compute the initial T, X and XPRIME: 0: The initial T, X and XPRIME are consistent. 1: Given X_d calculate X_a and X'_d 2: Given X' calculate X. ( X_d differential variables in X X_a algebrac variables in X ) IPAR(136), Total number of error test failures so far. IPAR(137), Total number of convergence test failures. -Parameters for DASPK IPAR(39), DASPK use: 0: direct methods (compatible with DASSL) 1: Krylov method 2: Krylov method + Jac IPAR(40), DASPK uses scalars MAXLm KMP, NRMAX and EPLI when uses Krylov method. 0: uses default values. 1: uses user values. IPAR(41), Proceed to the integration after the initial condition calculation is done. Used when INFOV(11)>0 0: Yes 1: No IPAR(42), Errors are controled localy on all the variables. 0: Yes 1: No IPAR(43), Use default values for initial condition heuristic controls. 0: Yes 1: No and provide MXNIT, MXNJ, MXNH, LSOFF, STPTOL, EPINIT. IPAR(138), number of convergence failures of the linear iteration IPAR(139), length of IWORK actually required. IPAR(140), length of RWORK actually required. IPAR(141), total number of nonlinear iterations. IPAR(142), total number of linear (Krylov) iterations IPAR(143), number of PSOL calls. DPAR (input/output) DOUBLE PRECISION array, dimension (202) INPUT: 1..15 General 16..25 ODEPACK 26..35 RADAU5 36..50 DASSL/PK 51..60 GELDA 61..100 Reserved OUTPUT: 101..110 General 111..125 ODEPACK 126..135 RADAU5 136..145 DASSL/PK 146..155 GELDA 156..200 Reserved Any Mode: 201.. User Available Common real parameters for SOLVERS: DPAR(1), Initial step size guess. Optional in: ODEPACK, DASSL, .. DPAR(2), Maximum absolute step size allowed. Common parameters for ODEPACK and DASSL: DPAR(111), Step size in t last used (successfully). DPAR(113), Current value of the independent variable which the solver has actually reached Common parameters for ODEPACK solvers: DPAR(16), Critical value of t which the solver is not overshoot. DPAR(17), Minimum absolute step size allowed. DPAR(112), Step size to be attempted on the next step. DPAR(18), Tolerance scale factor, greater than 1.0. - LSODA DPAR(115) Value of t at the time of the last method switch, if any. - LSODA DPAR(115) Value of t at the time of the last method switch, if any. - LSODES DPAR(19), The element threshhold for sparsity determination when moss = 1 or 2. Parameters for RADAU5 solver: DPAR(26), The rounding unit, default 1E-16. DPAR(27), The safety factor in step size prediction, default 0.9D0. DPAR(28), Decides whether the jacobian should be recomputed, default 0.001D0. DPAR(29), Stopping criterion for Newton's method, default MIN(0.03D0, RTOL(1)**0.5D0) DPAR(30), DPAR(31): This saves, together with a large DPAR(28), LU-decompositions and computing time for large systems. DPAR(32), DPAR(33), Parameters for step size selection. Parameters for DASSL and DASPK solvers: DPAR(36), Stopping point (Tstop)Tolerances
RTOL DOUBLE PREISION Relative Tolerance. ATOL DOUBLE PREISION Absolute Tolerance.Workspace
IWORK INTEGER array, dimension (LIWORK) LIWORK INTEGER Size of IWORK, depending on solver: - LSODE 20 for mf = 10, 20 + neq for mf = 21, 22, 24, or 25. if mf = 24, 25, input in iwork(1),iwork(2) the lower and upper half-bandwidths ml,mu. -LSODA 20+NX -LSODES 30 -DASSL 20+NEQ DWORK DOUBLE PREISION array, dimension (LDWORK) LDWORK INTEGER Size of DWORK, depending on solver: - LSODE 20+16*NX , IPAR(3) = 10, 22+ 9*NX+NX**2 , IPAR(3) = 21 or 22, 22+10*NX+(2*IPAR(4)+IPAR(9))*NX, IPAR(3) = 24 or 25. - LSODA 22+NX*max(16,NX+9) - LSODES 20+16*NX , mf=10 20+(2+1./lenrat)*nnz + (11+9./lenrat)*NX, mf=121,222 - DASSL >= 40 LRW .GE. 40+(MAXORD+4)*NEQ+NEQ**2, IPAR(3) = 1 or 10 >= 40+(MAXORD+4)*NEQ+(2*ML+MU+1)*NEQ, IPAR(3) = 2 >= 40+(MAXORD+4)*NEQ+(2*ML+MU+1)*NEQ +2*(NEQ/(ML+MU+1)+1), IPAR(3) = 11Warning Indicator
IWARN INTEGER = 0: no warning; = 1: LSODE/LSODA/LSODES/RADAU5 do not use the input vector as argument; = 2: Only the 1st element of RTOL is used; = 3: Method (IPAR(3)) not allowed with LSODE/LSODA/LSODES/RADAU5/DASSL/DASPK; = 4: Only the 1st element of ATOL is used; = 5: Option not allowed for IPAR(37); = 6: Option not allowed for IPAR(38).Error Indicator
INFO INTEGER = 0: Successful exit; < 0: If INFO = -i, the i-th argument had an illegal value; = 1: Wrong tolerance mode; = 2: Sparse storage (IPAR(4)=1) incompatible with LSODE/LSODA/RADAU5; = 3: Dense storage (IPAR(4)=0) incompatible with LSODES = 100+ERROR: ODEDER returned ERROR = 200+ERROR: RADAU5 returned -ERROR = 300+ERROR: DDASSL returned -ERROR = 400+ERROR: DDASPK returned -ERROR = 500+ERROR: DGELDA returned -ERROR
Since the package integrates 9 different solvers, it is possible to solve differential equations by means of Backward Differential Formulas, Runge-Kutta, using direct or iterative methods (including preconditioning) for the linear system associated, differential equations with time-varying coefficients or of order higher than one. The interface facilitates the user the work of changing the integrator and testing the results, thus leading a more robust and efficient integrated package.References
[1] A.C. Hindmarsh, Brief Description of ODEPACK: A Systematized Collection of ODE Solvers, http://www.netlib.org/odepack/doc [2] L.R. Petzold DASSL Library Documentation, http://www.netlib.org/ode/ [3] P.N. Brown, A.C. Hindmarsh, L.R. Petzold, DASPK Package 1995 Revision [4] R.S. Maier, Using DASPK on the TMC CM5. Experiences with Two Programming Models, Minesota Supercomputer Center, Technical Report. [5] E. Hairer, G. Wanner, Solving Ordinary Dirential Equations II. Stiánd Dirential- Algebraic Problems., Springer Seried in Computational Mathermatics 14, Springer-Verlag 1991, Second Edition 1996. [6] P. Kunkel, V. Mehrmann, W. Rath und J. Weickert, `GELDA: A Software Package for the Solution of General Linear Dirential Algebraic equations', SIAM Journal Scienti^Lc Computing, Vol. 18, 1997, pp. 115 - 138. [7] M. Otter, DSblock: A neutral description of dynamic systems. Version 3.3, http://www.netlib.org/odepack/doc [8] M. Otter, H. Elmqvist, The DSblock model interface for exchanging model components, Proceedings of EUROSIM 95, ed. F.Brenenecker, Vienna, Sep. 11-15, 1995 [9] M. Otter, The DSblock model interface, version 4.0, Incomplete Draft, http://dv.op.dlr.de/~otter7dsblock/dsblock4.0a.html [10] Ch. Lubich, U. Novak, U. Pohle, Ch. Engstler, MEXX - Numerical Software for the Integration of Constrained Mechanical Multibody Systems, http://www.netlib.org/odepack/doc [11] Working Group on Software (WGS), SLICOT Implementation and Documentation Standards (version 1.0), WGS-Report 90-1, Eindhoven University of Technology, May 1990. [12] P. Kunkel and V. Mehrmann, Canonical forms for linear differential- algebraic equations with variable coeÆcients., J. Comput. Appl. Math., 56:225{259, 1994. [13] Working Group on Software (WGS), SLICOT Implementation and Documentation Standards, WGS-Report 96-1, Eindhoven University of Technology, updated: Feb. 1998, ../../REPORTS/rep96-1.ps.Z. [14] A. Varga, Standarization of Interface for Nonlinear Systems Software in SLICOT, Deutsches Zentrum ur Luft un Raumfahrt, DLR. SLICOT-Working Note 1998-4, 1998, Available at ../../REPORTS/SLWN1998-4.ps.Z. [15] D. Kirk, Optimal Control Theory: An Introduction, Prentice-Hall. Englewood Cli, NJ, 1970. [16] F.L. Lewis and V.L. Syrmos, Optimal Control, Addison-Wesley. New York, 1995. [17] W.M.Lioen, J.J.B de Swart, Test Set for Initial Value Problem Solvers, Technical Report NM-R9615, CWI, Amsterdam, 1996. http://www.cwi.nl/cwi/projects/IVPTestset/. [18] V.Hernandez, I.Blanquer, E.Arias, and P.Ruiz, Definition and Implementation of a SLICOT Standard Interface and the associated MATLAB Gateway for the Solution of Nonlinear Control Systems by using ODE and DAE Packages}, Universidad Politecnica de Valencia, DSIC. SLICOT Working Note 2000-3: July 2000. Available at ../../REPORTS/SLWN2000-3.ps.Z. [19] J.J.B. de Swart, W.M. Lioen, W.A. van der Veen, SIDE, November 25, 1998. Available at http://www.cwi.nl/cwi/projects/PSIDE/. [20] Kim, H.Young, F.L.Lewis, D.M.Dawson, Intelligent optimal control of robotic manipulators using neural networks. [21] J.C.Fernandez, E.Arias, V.Hernandez, L.Penalver, High Performance Algorithm for Tracking Trajectories of Robot Manipulators, Preprints of the Proceedings of the 6th IFAC International Workshop on Algorithms and Architectures for Real-Time Control (AARTC-2000), pages 127-134.Numerical Aspects
The numerical aspects of the routine lie on the features of the different packages integrated. Several packages are more robust than others, and other packages simply cannot deal with problems that others do. For a detailed description of the numerical aspects of each method is recommended to check the references above.Further Comments
Several packages (LSODES, LSOIBT) deal only with sparse matrices. The interface checks the suitability of the methods to the parameters and show a warning message if problems could arise.Example
Program Text
* ODESOLVER EXAMPLE PROGRAM TEXT FOR LSODEX
PROBLEM
* Copyright (c) 2002-2017 NICONET e.V.
*
* .. Parameters ..
INTEGER
NIN, NOUT
PARAMETER
( NIN = 5, NOUT = 6 )
INTEGER LSODE_, LSODA_, LSODES_,
RADAU5_, DASSL_, DASPK_, DGELDA_
PARAMETER (LSODE_ = 1, LSODA_
= 2, LSODES_ = 3)
PARAMETER (RADAU5_ = 4, DASSL_ =
5, DASPK_ = 6)
PARAMETER (DGELDA_ = 7)
*
EXTERNAL IARGC_
INTEGER IARGC_
INTEGER NUMARGS
CHARACTER*80 NAME
CHARACTER*80 SOLVER
*
* .. Executable Statements ..
*
WRITE ( NOUT, FMT = 99999 )
*
NUMARGS = IARGC_()
*
CALL GETARG_(0, NAME)
IF (NUMARGS .NE. 1) THEN
WRITE (*,*) 'Syntax
Error: ',NAME(1:8),' <solver>'
WRITE (*,*) 'Solvers
: LSODE, LSODA, LSODES, RADAU5, DASSL, DASP
&K'
ELSE
*
CALL GETARG_(1, SOLVER)
*
WRITE (*,*) 'Problem:
LSODEX Solver: ',SOLVER(1:7)
*
IF (SOLVER(1:5) .EQ.
'LSODE') THEN
CALL TEST(LSODE_)
ELSEIF (SOLVER(1:5)
.EQ. 'LSODA') THEN
CALL TEST(LSODA_)
ELSEIF (SOLVER(1:6)
.EQ. 'LSODES') THEN
CALL TEST(LSODES_)
ELSEIF (SOLVER(1:6)
.EQ. 'RADAU5') THEN
CALL TEST(RADAU5_)
ELSEIF (SOLVER(1:5)
.EQ. 'DASSL') THEN
CALL TEST(DASSL_)
ELSEIF (SOLVER(1:5)
.EQ. 'DASPK') THEN
CALL TEST(DASPK_)
ELSE
WRITE (*,*)
'Error: Solver: ', SOLVER,' unknown'
ENDIF
ENDIF
*
99999 FORMAT (' ODESOLVER EXAMPLE PROGRAM RESULTS FOR LSODEX PROBLEM'
.
,/1X)
END
*
*
*
*
SUBROUTINE TEST( ISOLVER )
*
* Copyright (c) 2002-2017 NICONET e.V.
*
* PURPOSE
*
* Testing subroutine ODESolver
*
* ARGUMENTS
*
* Input/Output Parameters
*
* ISOLVER (input) INTEGER
*
Indicates the nonlinear solver package to be used:
*
= 1: LSODE,
*
= 2: LSODA,
*
= 3: LSODES,
*
= 4: RADAU5,
*
= 5: DASSL,
*
= 6: DASPK,
*
= 7: DGELDA.
*
* METHOD
*
* REFERENCES
*
* CONTRIBUTORS
*
* REVISIONS
*
* -
*
* KEYWORDS
*
*
* ******************************************************************
* .. Parameters ..
INTEGER LSODE_, LSODA_, LSODES_,
RADAU5_, DASSL_, DASPK_, DGELDA_
PARAMETER (LSODE_ = 1, LSODA_
= 2, LSODES_ = 3)
PARAMETER (RADAU5_ = 4, DASSL_ =
5, DASPK_ = 6)
PARAMETER (DGELDA_ = 7)
INTEGER
NIN, NOUT
PARAMETER
( NIN = 5, NOUT = 6 )
INTEGER
MD, ND, LPAR, LWORK
PARAMETER
( MD = 400, ND = 100, LPAR = 250,
$
LWORK = 650000 )
* .. Scalar Arguments ..
INTEGER ISOLVER
* .. Local Scalars ..
INTEGER
NEQN, NDISC, MLJAC, MUJAC, MLMAS, MUMAS
INTEGER
IWARN, INFO
DOUBLE PRECISION ATOL(MD), RTOL,
NORM
LOGICAL
NUMJAC, NUMMAS, CONSIS
* .. Local Arrays ..
CHARACTER FULLNM*40, PROBLM*8, TYPE*3
CHARACTER*9 ODEDER, ODEOUT, JACFX,
JACFU, JACFP
INTEGER
IND(MD), IPAR(LPAR), IWORK(LWORK)
DOUBLE PRECISION T(0:ND), RPAR(LPAR),
DWORK(LWORK)
DOUBLE PRECISION X(MD), XPRIME(MD),
U(MD), Y(MD)
* .. External Functions ..
DOUBLE PRECISION DNRM2
EXTERNAL
DNRM2
* .. External Subroutines ..
EXTERNAL
PLSODEX,ILSODEX,SLSODEX
EXTERNAL
DAXPY
* .. Executable Statements ..
*
DO 20 I=1,NEQN
U(I)=0D0
Y(I)=0D0
20 CONTINUE
DO 40 I=1,LPAR
IPAR(I)=0
RPAR(I)=0D0
40 CONTINUE
DO 60 I=1,LWORK
IWORK(I)=0
DWORK(I)=0D0
60 CONTINUE
IPAR(2)=1
* Get the problem dependent parameters.
RPAR(1)=1D-3
IPAR(1)=0
ATOL(1)=1D-6
ATOL(2)=1D-10
ATOL(3)=1D-6
RTOL=1D-4
CALL PLSODEX(FULLNM,PROBLM,TYPE,NEQN,NDISC,T,NUMJAC,MLJAC,
$
MUJAC,NUMMAS,MLMAS,MUMAS,IND)
CALL ILSODEX(NEQN,T(0),X,XPRIME,CONSIS)
CALL SLSODEX(NEQN,T(1),XPRIME)
IF ( TYPE.NE.'ODE' ) THEN
WRITE ( NOUT,
FMT = 99998 )
ELSE
WRITE ( NOUT,
FMT = 99997 ) FULLNM, PROBLM, TYPE, ISOLVER
IF ( NUMJAC )
THEN
IPAR(3)=0
ELSE
IPAR(3)=1
END IF
IPAR(6)=MLJAC
IPAR(7)=MUJAC
ODEDER=''
ODEOUT=''
JACFX=''
JACFU=''
JACFP=''
CALL ODESolver(
ISOLVER, ODEDER, ODEOUT, JACFX, JACFU, JACFP,
$
NEQN, NEQN, NEQN, T(0), T(1), X, U, Y,
$
IPAR, RPAR, RTOL, ATOL,
$
IWORK, LWORK, DWORK, LWORK, IWARN, INFO )
IF ( INFO.NE.0
) THEN
WRITE ( NOUT, FMT = 99996 ) INFO
ELSE
IF ( IWARN.NE.0 ) THEN
WRITE ( NOUT, FMT = 99995 ) IWARN
ENDIF
IF ( NEQN .LE. 10 ) THEN
WRITE ( NOUT, FMT = 99994 )
DO 80 I=1,NEQN
WRITE ( NOUT, FMT = 99993 ) X(I), XPRIME(I)
80
CONTINUE
END IF
NORM=DNRM2(NEQN,XPRIME,1)
IF ( NORM.EQ.0D0 ) THEN
NORM=1D0
END IF
CALL DAXPY(NEQN,-1D0,X,1,XPRIME,1)
NORM=DNRM2(NEQN,XPRIME,1)/NORM
WRITE ( NOUT, FMT = 99992 ) NORM
END IF
END IF
*
99998 FORMAT (' ERROR: This test is only intended for ODE problems')
99997 FORMAT (' ',A,' (',A,' , ',A,') with SOLVER ',I2)
99996 FORMAT (' INFO on exit from ODESolver = ',I3)
99995 FORMAT (' IWARN on exit from ODESolver = ',I3)
99994 FORMAT (' Solution: (calculated) (reference)')
99993 FORMAT (F,F)
99992 FORMAT (' Relative error comparing with the reference solution:'
$
,E,/1X)
* *** Last line of TEST ***
END
SUBROUTINE ODEDER_( NX, NU, T, X,
U, RPAR, IPAR, F, INFO )
*
* Copyright (c) 2002-2017 NICONET e.V.
*
* PURPOSE
*
* Interface routine between ODESolver and
the problem function FEVAL
*
* ARGUMENTS
*
* Input/Output Parameters
*
* NX
(input) INTEGER
*
Dimension of the state vector.
*
* NU
(input) INTEGER
*
Dimension of the input vector.
*
* T
(input) INTEGER
*
The time point where the function is evaluated.
*
* X
(input) DOUBLE PRECISION array, dimension (NX)
*
Array containing the state variables.
*
* U
(input) DOUBLE PRECISION array, dimension (NU)
*
Array containing the input values.
*
* RPAR (input/output)
DOUBLE PRECISION array
*
Array for communication between the driver and FEVAL.
*
* IPAR (input/output)
INTEGER array
*
Array for communication between the driver and FEVAL.
*
* F
(output) DOUBLE PRECISION array, dimension (NX)
*
The resulting function value f(T,X).
*
* Error Indicator
*
* INFO INTEGER
*
Return values of error from FEVAL or 100 in case
*
a bad problem was choosen.
*
* METHOD
*
* REFERENCES
*
* CONTRIBUTORS
*
* REVISIONS
*
* -
*
* KEYWORDS
*
*
* ******************************************************************
*
* .. Scalar Arguments ..
INTEGER
NX, NU, INFO
DOUBLE PRECISION T
* .. Array Arguments ..
INTEGER
IPAR(*)
DOUBLE PRECISION X(NX), U(NU), RPAR(*),
F(NX)
* .. External Subroutines ..
EXTERNAL
FLSODEX
* .. Executable Statements ..
CALL FLSODEX(NX,T,X,X,F,INFO,RPAR,IPAR)
* *** Last line of ODEDER_ ***
END
SUBROUTINE JACFX_( NX, DUMMY,
LDFX, T, X, DUMMY2, RPAR, IPAR, FX,
$
INFO )
*
* Copyright (c) 2002-2017 NICONET e.V.
*
* PURPOSE
*
* Interface routine between ODESolver and
the problem function JEVAL
*
* ARGUMENTS
*
* Input/Output Parameters
*
* NX
(input) INTEGER
*
Dimension of the state vector.
*
* DUMMY (input) INTEGER
*
* LDFX (input)
INTEGER
*
The leading dimension of the array FX.
*
* T
(input) INTEGER
*
The time point where the derivative is evaluated.
*
* X
(input) DOUBLE PRECISION array, dimension (NX)
*
Array containing the state variables.
*
* DUMMY2 (input) DOUBLE PRECISION
*
* RPAR (input/output)
DOUBLE PRECISION array
*
Array for communication between the driver and FEVAL.
*
* IPAR (input/output)
INTEGER array
*
Array for communication between the driver and FEVAL.
*
* FX
(output) DOUBLE PRECISION array, dimension (LDFX,NX)
*
The array with the resulting Jacobian matrix.
*
* Error Indicator
*
* INFO INTEGER
*
Return values of error from JEVAL or 100 in case
*
a bad problem was choosen.
*
* METHOD
*
* REFERENCES
*
* CONTRIBUTORS
*
* REVISIONS
*
* -
*
* KEYWORDS
*
*
* ******************************************************************
*
* .. Scalar Arguments ..
INTEGER
NX, DUMMY, LDFX, INFO
DOUBLE PRECISION T
* .. Array Arguments ..
INTEGER
IPAR(*)
DOUBLE PRECISION X(NX), DUMMY2(*),
RPAR(*), FX(NX)
* .. External Subroutines ..
EXTERNAL
JLSODEX
* .. Executable Statements ..
CALL JLSODEX(LDFX,NX,T,X,X,FX,INFO,RPAR,IPAR)
* *** Last line of JACFX_ ***
END
Program Data
No data requiredProgram Results
ODESOLVER EXAMPLE PROGRAM RESULTS Problem: LSODEX Solver: LSODE IWARN on exit from ODESolver = 1 Solution: (calculated) 8.287534436182735E-08 3.329129749822125E-13 1.118553835127275E-07