Purpose
Interface for using a common entry point, DSblock compatible for defining Differential Algebraic Equations using several packages.
The equations follow the form (CASE A): F(dx(t)/dt, x(t), u(t), p, t) = 0 y(t) = g(dx(t)/dx, x(t), u(t), p, t) for the most general model which can only be solved by DASSL and DASSPK. A restricted case can be solved with RADAU5, LSODI, LSOIBT, if the system is expressed as (CASE B): F(x(t), u(t), p, t)*dx(t)/dt = A(x(t), u(t), p, t) y(t) = g(dx(t)/dx, x(t), u(t), p, t) And finally, the GELDA package is able to solve DAEs with the expression (CASE C): F(u(t), p, t)*dx(t)/dt = A(u(t), p, t)*x(t) + E(u(t), p, t) The user must define the subroutines: DAEDF: F(dx(t)/dt, x(t), u(t), p, t) for CASES: A, B and C DAEDA: A(x(t), u(t), p, t) for CASES: B and C DAEDE: E(u(t), p, t) for CASES: C DAEOUT: g(dx(t)/dx, x(t), u(t), p, t) and the Jacobians (JACFX, JACFU, JACFP) if used. The interface adapts the structure to fit all the codesSpecification
SUBROUTINE DAESolver(ISOLVER,CDAEDF_,CDAEDA_,CDAEDE_,CDAEOUT_, $ CJACFX_,CJACFU_,CJACFP_,CJACFXDOT_, $ NX, NY, NU, NP, TINI, TOUT, $ X, XDOTI, Y, U, P, $ IPAR, RPAR, RTOL, ATOL, $ IWORK, LIWORK, DWORK, LDWORK, $ IWARN, INFO) .. Scalar Arguments .. DOUBLE PRECISION TINI, TOUT INTEGER ISOLVER, IWARN, INFO, $ NX, NY, NU, NP, $ LDWORK, LIWORK CHARACTER*9 CDAEDF_, CDAEDA_,CDAEDE_, CDAEOUT_, $ CJACFX_, CJACFU_, CJACFP_, CJACFXDOT_, $ CDAEDF, CDAEDA,CDAEDE, CDAEOUT, $ CJACFX, CJACFU, CJACFP, CJACFXDOT .. Array Arguments .. DOUBLE PRECISION DWORK(LDWORK), RPAR(*), ATOL(*), RTOL(*) $ X(NX), XDOTI(NX), Y(NY), U(NU), P(NP) INTEGER IWORK(LIWORK), IPAR(*)Arguments
Mode Parameters
ISOLVER INTEGER Indicates the nonlinear solver packages to be used = 1: LSODI, = 2: LSOIBT, = 3: RADAU5, = 4: DASSL, = 5: DASPK, = 6: DGELDA.Input/Output Parameters
DAEDF (input) EXTERNAL Evaluates the F(dx(t)/dt, x(t), u(t), p, t). DAEDA (input) EXTERNAL Evaluates the A(x(t), u(t), p, t). DAEDE (input) EXTERNAL Evaluates the E(u(t), p, t). DAEOUT (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. NP (input) INTEGER Dimension of the parameter 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. XDOTI (input) DOUBLE PRECISION array, dimension (NX) Array containing dx(t)/dt at initial point. 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. U (input) DOUBLE PRECISION array, dimension (NU) Array containing the input initial values. P (input) DOUBLE PRECISION array, dimension (NP) Array containing the parameter variables. IPAR (input/output) INTEGER array, dimension (201) 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 a scalar and atol is a vector 2 : both rtol and atol are vectors IPAR(2), Compute Output Values only at TOUT (and not at the intermediate step). (1:Yes, 0:No) IPAR(3), mfjac, Method flag for jacobian 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(6), ml, lower half-bandwithds of the banded jacobian, excluding tne main diagonal. IPAR(7), mu, upper half-bandwithds of the banded jacobian, excluding the main diagonal. (Note: IPAR(6) and IPAR(7) are obligatories only if the jacobian matrix is banded) IPAR(101) = Number of steps taken for the problem. IPAR(102) = Number of residual evaluations. IPAR(103) = Number of jacobian evaluations. Common parameters for RADAU5, ODEPACK and DGELDA: IPAR(9), mfmass, Method flag for mass-matrix 0 : No mass-matrix used (non-stiff method). 1 : User supplied full mass-matrix (stiff). 2 : User supplied banded mass-matrix (stiff). 10 : Identity mass-matrix is used (stiff). IPAR(10), mlmass, lower half-bandwithds of the banded mass matrix, excluding the main diagonal. IPAR(11), mumass, upper half-bandwithds of the banded mass matrix, excluding the main diagonal. IPAR(12), Maximum number of steps allowed during one call to the solver. Common parameters for ODEPACK, DASSL, DASPK and DGELDA: IPAR(13), Maximum order to be allowed. default values : 12 if meth = 1 5 if meth = 2 If exceds the default value, it will be reduced to the default value. In DASSL, DASPK and DGELDA : (1 .LE. MAXORD .LE. 5) IPAR(111) = The method order last used(successfully). IPAR(112) = The order to be attempted on the next step. Common parameters for ODEPACK package: IPAR(16), Status Flag IPAR(17), Optional inputs, must be 0 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. - LSOIBIT IPAR(24), mb, block size. (mb .GE. 1) and mb*IPAR(28) = NX IPAR(25), nb, number of blocks in the main diagonal. (nb .ge. 4) and nb*IPAR(27) = NX - RADAU5 IPAR(26) Transforms the Jacobian matrix to Hessenberg form.(Only if IPAR(9)=1 and IPAR(3)=1 or 10) IPAR(27) Maximum number of Newton iterations in each step. IPAR(28) Starting values for Newton's method .EQ. 0 -> is taken the extrapolated collocation solution .NE. 0 -> zero values are used. IPAR(29) Dimension of the index 1 variables( >0 ). 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 (default 0). IPAR(34) Value of M2 (default(M2=M1). 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, DASPK and DGELDA 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(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. Common parameters for DASSL and DASPK solvers: IPAR(37), code solve the problem without invoking any special non negativity constraints: 0: Yes 1: To have constraint checking only in the initial condition calculation. 2: To enforze nonnegativity in X during the integration. 3: To enforce both options 1 and 2. IPAR(137), Total number of convergence test failures. - DASPK IPAR(39), DASPK use: 0: direct methods (dense or band) 1: Krylov method (iterative) 2: Krylov method + Jac (iterative) IPAR(41), Proceed to the integration after the initial condition calculation is done. Used when IPAR(38) > 0: 0: Yes 1: No IPAR(42), Errors are controled localy on all the variables: 0:Yes 1: No IPAR(8), Extra printing 0, no printing 1, for minimal printing 2, for full printing IPAR(44), maximum number of iterations in the SPIGMR algorithm. (.LE. NX) IPAR(45), number of vectors on which orthogonalization is done in the SPIGMR algorithm. (.LE. IPAR(44)) IPAR(46), maximum number of restarts of the SPIGMR algorithm per nonlinear iteration. (.GE. 0) IPAR(47), maximum number of Newton iterations per Jacobian or preconditioner evaluation. (> 0) IPAR(48), maximum number of Jacobian or preconditioner evaluations. (> 0) IPAR(49), maximum number of values of the artificial stepsize parameter H to be tried if IPAR(38) = 1. (> 0). IPAR(50), flag to turn off the linesearch algorithm. 0 : ON 1 : OFF (default) 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. - DGELDA IPAR(51), contains the strangeness index. IPAR(52), number of differential components IPAR(53), number of algebraic components IPAR(54), number of undetermined components IPAR(55), method used: if 1 then uses the BDF solver 2 then uses the Runge-Kutta solver IPAR(56), E(t) and A(t) are: 1 time dependent 0 constants IPAR(57), Maximum index of the problem. ( .GE. 0 ) IPAR(58), Step size strategy: 0, Mod. predictive controlled of Gustafsson(safer) 1, classical step size control(faster) RPAR (input/output) DOUBLE PRECISION array, dimension (201) 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 parameters for solvers: RPAR(1), Initial step size guess.Obligatory in RADAU5. RPAR(2), Maximum absolute step size allowed. Common parameters for ODEPACK, DASSL, DASPK and DGELDA: RPAR(111), Step size in t last used (successfully). RPAR(112), Step size to be attempted on the next step. RPAR(113), Current value of the independent variable which the solver has actually reached Common parameters for ODEPACK solver: RPAR(16), Critical value of t which the solver is not overshoot. RPAR(17), Minimum absolute step size allowed. RPAR(18), Tolerance scale factor, greater than 1.0. Parameters for RADAU5 solver: RPAR(26), The rounding unit, default 1E-16. RPAR(27), The safety factor in step size prediction, default 0.9D0. RPAR(28), Decides whether the jacobian should be recomputed, default 0.001D0. Increase when jacobian evaluations are costly For small systems should be smaller. RPAR(29), Stopping criterion for Newton's method, default MIN(0.03D0, RTOL(1)**0.5D0). RPAR(30), RPAR(31): This saves, together with a large RPAR(28), LU-decompositions and computing time for large systems. Small systems: RPAR(30)=1.D0, RPAR(31)=1.2D0 Large full systems: RPAR(30)=0.99D0, RPAR(31)=2.D0 might be good. RPAR(32), RPAR(33), Parameters for step size selection.Condition: RPAR(32)<=HNEW/HOLD<=RPAR(33) Parameters for DASSL, DASPK and DGELDA solvers: RPAR(36), Stopping point (Tstop) - DASPK RPAR(37), convergence test constant in SPIGMR algorithm. (0 .LT. RPAR(37) .LT. 1.0) RPAR(38), minimum scaled step in linesearch algorithm. The default is = (unit roundoff)**(2/3). (> 0) RPAR(39), swing factor in the Newton iteration convergence test. (default 0.1) (> 0) - DASPK RPAR(40), safety factor used in step size prediction. RPAR(41) and RPAR(42) restric the relation between the new and old stepsize in step size selection. 1/RPAR(41) .LE. Hnew/Hold .LE. 1/RPAR(42) RPAR(43), RPAR(44) QUOT1 and QUOT2 repectively. If QUOT1 < Hnew/Hold < QUOT2 and A and E are constants, the work can be saved by setting Hnew=Hold and using the system matrix of the previous step.Tolerances
RTOL DOUBLE PRECISION Relative Tolerance. ATOL DOUBLE PRECISION Absolute Tolerance.Workspace
IWORK INTEGER array, dimension (LIWORK) LIWORK INTEGER Minimum size of DWORK, depending on solver: - LSODI, LSOIBT, DASSL 20 + NX - RADAU5 3*N+20 DWORK DOUBLE PRECISION array, dimension (LDWORK) LDWORK INTEGER Size of DWORK, depending on solver: - LSODI 22 + 9*NX + NX**2 , IPAR(3) = 1 or 10 22 + 10*NX + (2*ML + MU)*NX , IPAR(3) = 2 or 11 - LSOIBT 20 + nyh*(maxord + 1) + 3*NX + lenw where nyh = Initial value of NX maxord = Maximum order allowed(default or IPAR(13) lenw = 3*mb*mb*nb + 2 - RADAU5 N*(LJAC+LMAS+3*LE+12)+20 where LJAC=N if (full jacobian) LJAC=MLJAC+MUJAC+1 if (banded jacobian) and LMAS=0 if (IPAR(9) = 10 or 11) LMAS=N if (IPAR(9) = 1) LMAS=MLMAS+MUMAS+1 if (IPAR(9) = 2) and LE=N if (IPAR(9) = 1 or 10) LE=2*MLJAC+MUJAC+1 if (IPAR(9) = 2 or 11) - 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: LSODI/LSOIBT/RADAU5 do not use the input vector as argument; = 2: LSODI/LSOIBT do not use the param vector as argument; = 3: RTOL and ATOL are used as scalars;Error Indicator
INFO INTEGER = 0: Successful exit; < 0: If INFO = -i, the i-th argument had an illegal value; = 1: Wrong tolerance mode; = 2: Method (IPAR(9)) is not allowed for ODEPACK/RADAU5; = 3: Method (IPAR(3)) is not allowed for LSODE/RADAU5/DASSL; = 4: Option not allowed for IPAR(37); = 5: Option not allowed for IPAR(38); = 100+ERROR: RADAU5 returned -ERROR; = 200+ERROR: DASSL returned -ERROR; = 300+ERROR: DASPK returned -ERROR; = 400+ERROR: DGELDA returned -ERROR.Method
Since the package integrates 8 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
* DAESOLVER EXAMPLE PROGRAM TEXT FOR LSODIX
PROBLEM
* Copyright (c) 2002-2017 NICONET e.V.
*
* .. Parameters ..
INTEGER
NIN, NOUT
PARAMETER
( NIN = 5, NOUT = 6 )
INTEGER LSODI_, LSOIBT_, RADAU5_,
DASSL_, DASPK_, GELDA_
PARAMETER (LSODI_ = 1, LSOIBT_
= 2)
PARAMETER (RADAU5_ = 3, DASSL_
= 4, DASPK_ = 5)
PARAMETER (GELDA_ = 6)
* .. Executable Statements ..
*
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
: LSODI, LSOIBT, RADAU5, DASSL, DASPK, GELD
&A'
ELSE
*
CALL GETARG_(1, SOLVER)
*
WRITE (*,*) 'Problem:
LSODIX Solver: ',SOLVER(1:7)
*
IF (SOLVER(1:5) .EQ.
'LSODI') THEN
CALL TEST(LSODI_)
ELSEIF (SOLVER(1:6)
.EQ. 'LSOIBT') THEN
CALL TEST(LSOIBT_)
ELSEIF (SOLVER(1:6)
.EQ. 'RADAU5') THEN
CALL TEST(RADAU5_)
ELSEIF (SOLVER(1:5)
.EQ. 'GELDA') THEN
CALL TEST(GELDA_)
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 (' DAESOLVER EXAMPLE PROGRAM RESULTS FOR LSODIX PROBLEM'
.
,/1X)
END
SUBROUTINE TEST( ISOLVER )
*
* Copyright (c) 2002-2017 NICONET e.V.
*
* PURPOSE
*
* Testing subroutine DAESolver
*
* ARGUMENTS
*
* Input/Output Parameters
*
* ISOLVER (input) INTEGER
*
Indicates the nonlinear solver package to be used:
*
= 1: LSODI,
*
= 2: LSOIBT,
*
= 3: RADAU5,
*
= 4: DASSL,
*
= 5: DASPK,
*
= 6: DGELDA.
*
* METHOD
*
* REFERENCES
*
* CONTRIBUTORS
*
* REVISIONS
*
* -
*
* KEYWORDS
*
*
* ******************************************************************
* .. Parameters ..
INTEGER LSODI_, LSOIBT_, RADAU5_,
DASSL_, DASPK_, GELDA_
PARAMETER (LSODI_ = 1, LSOIBT_
= 2)
PARAMETER (RADAU5_ = 3, DASSL_
= 4, DASPK_ = 5)
PARAMETER (GELDA_ = 6)
INTEGER
NIN, NOUT
PARAMETER
( NIN = 5, NOUT = 6 )
INTEGER
MD, ND, LPAR, LWORK
PARAMETER
( MD = 400, ND = 100, LPAR = 201,
$
LWORK = 10000 )
* .. Common variables ..
COMMON /TESTING/ ISOLVER2
INTEGER ISOLVER2
* .. Scalar Arguments ..
INTEGER ISOLVER
* .. Local Scalars ..
INTEGER
NEQN, NDISC, MLJAC, MUJAC, MLMAS, MUMAS
INTEGER
IWARN, INFO
DOUBLE PRECISION ATOL, RTOL, NORM
LOGICAL
NUMJAC, NUMMAS, CONSIS
* .. Local Arrays ..
CHARACTER FULLNM*40, PROBLM*8, TYPE*3
CHARACTER*9 CDAEDF,CDAEDA,CDAEDE,CDAEOUT,
$
CJACFX,CJACFU,CJACFP,CJACFXDOT
INTEGER
IND(MD), IPAR(LPAR), IWORK(LWORK)
DOUBLE PRECISION T(0:ND), RPAR(LPAR),
DWORK(LWORK)
DOUBLE PRECISION X(MD), XPRIME(MD),
Y(MD), U(MD), P(MD), SOLU(MD)
* .. External Functions ..
DOUBLE PRECISION DNRM2
EXTERNAL
DNRM2
* .. External Subroutines ..
EXTERNAL
PLSODIX, ILSODIX, SLSODIX
EXTERNAL
DAXPY
* .. Executable Statements ..
*
ISOLVER2 = ISOLVER
DO 20 I=1,NEQN
Y(I)=0D0
U(I)=0D0
P(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
* Get the problem dependent parameters.
RTOL=1D-4
ATOL=1D-6
IPAR(1)=0
IPAR(2)=1
IPAR(3)=1
IPAR(12)= 10000
IF (ISOLVER .EQ. LSODI_ .OR. ISOLVER
.EQ. RADAU5_) THEN
IPAR(9)=1
IPAR(16)=1
C IPAR(17)=0
RPAR(1)=1D-3
ELSE
C (ISOLVER
.EQ. DASSL_ .OR. ISOLVER .EQ. DASPK_)
C IPAR(36)=0
C IPAR(37)=0
C IPAR(38)=0
IPAR(39)=1
END IF
CALL PLSODIX(FULLNM,PROBLM,TYPE,NEQN,NDISC,T,NUMJAC,MLJAC,
$
MUJAC,NUMMAS,MLMAS,MUMAS,IND)
CALL ILSODIX(NEQN,T(0),X,XPRIME,CONSIS)
x(1) = 1.0d0
x(2) = 0.0d0
x(3) = 0.0d0
xprime(1) = -0.04D0
xprime(2) = 0.04D0
xprime(3) = 0.0D0
CALL SLSODIX(NEQN,T(1),SOLU)
IF ( TYPE.NE.'DAE' ) THEN
WRITE ( NOUT,
FMT = 99998 )
ELSE
WRITE ( NOUT,
FMT = 99997 ) FULLNM, PROBLM, TYPE, ISOLVER
CDAEDF=''
CDAEDA=''
CDAEDE=''
CDAEOUT=''
CJACFX=''
CJACFU=''
CJACFP=''
CJACFXDOT=''
CALL DAESolver(
ISOLVER, CDAEDF, CDAEDA, CDAEDE, CDAEOUT,
$
CJACFX, CJACFU, CJACFP, CJACFXDOT,
$
NEQN, NEQN, NEQN, NEQN, T(0), T(1),
$
X, XPRIME, Y, U, P,
$
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. 30 ) THEN
WRITE ( NOUT, FMT = 99994 )
DO 80 I=1,NEQN
WRITE ( NOUT, FMT = 99993 ) I, X(I), SOLU(I)
80
CONTINUE
END IF
NORM=DNRM2(NEQN,SOLU,1)
IF ( NORM.EQ.0D0 ) THEN
NORM=1D0
END IF
CALL DAXPY(NEQN,-1D0,X,1,SOLU,1)
NORM=DNRM2(NEQN,SOLU,1)/NORM
WRITE ( NOUT, FMT = 99992 ) NORM
END IF
END IF
*
99998 FORMAT (' ERROR: This test is only intended for DAE problems')
99997 FORMAT (' ',A,' (',A,' , ',A,') with SOLVER ',I2)
99996 FORMAT (' INFO on exit from DAESolver = ',I3)
99995 FORMAT (' IWARN on exit from DAESolver = ',I3)
99994 FORMAT (' Solution: (calculated) (reference)')
99993 FORMAT (I,F,F)
99992 FORMAT (' Relative error comparing with the reference solution:'
$
,E,/1X)
* *** Last line of TEST ***
END
SUBROUTINE DAEDA_( RPAR, NRP, IPAR,
NIP, X, NX, U, NU, P, NP,
$
F, LDF, T, INFO )
*
* Copyright (c) 2002-2017 NICONET e.V.
*
* PURPOSE
*
* Interface routine between DAESolver and
the problem function FEVAL
*
* ARGUMENTS
*
* Input/Output Parameters
*
* RPAR (input/output)
DOUBLE PRECISION array, dimension (NRP)
*
Array for communication between the driver and FEVAL.
*
* NRP (input)
INTEGER
*
Dimension of RPAR array.
*
* IPAR (input/output)
INTEGER array, dimension (NIP)
*
Array for communication between the driver and FEVAL.
*
* NIP (input)
INTEGER
*
Dimension of IPAR array.
*
* X
(input) DOUBLE PRECISION array, dimension (NX)
*
Array containing the state variables.
*
* NX
(input) INTEGER
*
Dimension of the state vector.
*
* U
(input) DOUBLE PRECISION array, dimension (NU)
*
Array containing the input values.
*
* NU
(input) INTEGER
*
Dimension of the input vector.
*
* P
(input) DOUBLE PRECISION array, dimension (NP)
*
Array containing the parameter values.
*
* NP
(input) INTEGER
*
Dimension of the parameter vector.
*
* F
(output) DOUBLE PRECISION array, dimension (LDF,NX)
*
The resulting function value f(T,X).
*
* LDF (input)
INTEGER
*
The leading dimension of F.
*
* T
(input) INTEGER
*
The time point where the function is evaluated.
*
* Error Indicator
*
* INFO INTEGER
*
Returns values of error from FEVAL or 100 in case
*
a bad problem was choosen.
*
* METHOD
*
* REFERENCES
*
* CONTRIBUTORS
*
* REVISIONS
*
* -
*
* KEYWORDS
*
*
* ******************************************************************
*
* .. Common variables ..
COMMON /TESTING/ ISOLVER
INTEGER LSODI_, LSOIBT_, RADAU5_,
DASSL_, DASPK_, GELDA_
PARAMETER (LSODI_ = 1, LSOIBT_
= 2)
PARAMETER (RADAU5_ = 3, DASSL_
= 4, DASPK_ = 5)
PARAMETER (GELDA_ = 6)
* .. Scalar Arguments ..
INTEGER
NRP, NIP, NX, NU, NP, LDF, INFO
DOUBLE PRECISION T
* .. Array Arguments ..
INTEGER
IPAR(NIP)
DOUBLE PRECISION RPAR(NRP), X(NX),
U(NU), P(NP),
$ F(LDF,NX)
* .. External Subroutines ..
EXTERNAL
FLSODIX
* .. Executable Statements ..
CALL FLSODIX(NX,T,X,X,F,INFO,RPAR,IPAR)
* *** Last line of DAEDA_ ***
END
SUBROUTINE DAEDF_( RPAR, NRP, IPAR,
NIP, X, XPRIME, NX,
$
U, NU, P, NP, T, F, LDF, INFO )
*
* Copyright (c) 2002-2017 NICONET e.V.
*
* PURPOSE
*
* Interface routine between DAESolver and
the problem function MEVAL
*
* ARGUMENTS
*
* Input/Output Parameters
*
* RPAR (input/output)
DOUBLE PRECISION array, dimension (NRP)
*
Array for communication between the driver and MEVAL.
*
* NRP (input)
INTEGER
*
Dimension of RPAR array.
*
* IPAR (input/output)
INTEGER array, dimension (NIP)
*
Array for communication between the driver and MEVAL.
*
* NIP (input)
INTEGER
*
Dimension of IPAR array.
*
* X
(input) DOUBLE PRECISION array, dimension (NX)
*
Array containing the state variables.
*
* XPRIME (input) DOUBLE PRECISION
array, dimension (NX)
*
Array containing the state variables derivative.
*
* NX
(input) INTEGER
*
Dimension of the state vector.
*
* U
(input) DOUBLE PRECISION array, dimension (NU)
*
Array containing the input values.
*
* NU
(input) INTEGER
*
Dimension of the input vector.
*
* P
(input) DOUBLE PRECISION array, dimension (NP)
*
Array containing the parameter values.
*
* NP
(input) INTEGER
*
Dimension of the parameter vector.
*
* T
(input) INTEGER
*
The time point where the function is evaluated.
*
* F
(output) DOUBLE PRECISION array, dimension (LDF,NX)
*
The resulting function value f(T,X).
*
* LDF (input)
INTEGER
*
The leading dimension of F.
*
* Error Indicator
*
* INFO INTEGER
*
Returns values of error from MEVAL or 100 in case
*
a bad problem was choosen.
*
* METHOD
*
* REFERENCES
*
* CONTRIBUTORS
*
* REVISIONS
*
* -
*
* KEYWORDS
*
*
* ******************************************************************
*
* .. Common variables ..
COMMON /TESTING/ ISOLVER
INTEGER ISOLVER
INTEGER LSODI_, LSOIBT_, RADAU5_,
DASSL_, DASPK_, GELDA_
PARAMETER (LSODI_ = 1, LSOIBT_
= 2)
PARAMETER (RADAU5_ = 3, DASSL_
= 4, DASPK_ = 5)
PARAMETER (GELDA_ = 6)
* .. Scalar Arguments ..
INTEGER
NRP, NIP, NX, NU, NP, LDF, INFO
DOUBLE PRECISION T
* .. Array Arguments ..
INTEGER
IPAR(NIP)
DOUBLE PRECISION RPAR(NRP), X(NX),
XPRIME(NX), U(NU), P(NP),
$ F(LDF,NX)
* .. Local Scalars ..
INTEGER
I
* .. External Subroutines ..
EXTERNAL
MLSODIX, RLSODIX
* .. Executable Statements ..
IF (ISOLVER .EQ. DASSL_ .OR. ISOLVER
.EQ. DASPK_) THEN
CALL RLSODIX(LDF,NX,T,X,XPRIME,F,INFO,RPAR,IPAR)
ELSE
CALL MLSODIX(LDF,NX,T,X,XPRIME,F,INFO,RPAR,IPAR)
ENDIF
* *** Last line of DAEDF_ ***
END
SUBROUTINE JACFX_( NRP, NIP,
RPAR, IPAR, NX, NU,
$
NP, X, U, P, T, FX, LDFX,
$
INFO )
*
* Copyright (c) 2002-2017 NICONET e.V.
*
* PURPOSE
*
* Interface routine between DAESolver and
the problem function JEVAL
*
* ARGUMENTS
*
* Input/Output Parameters
*
* NRP (input)
INTEGER
*
Dimension of RPAR array.
*
* NIP (input)
INTEGER
*
Dimension of IPAR array.
*
* RPAR (input/output)
DOUBLE PRECISION array
*
Array for communication between the driver and JEVAL.
*
* IPAR (input/output)
INTEGER array
*
Array for communication between the driver and JEVAL.
*
* NX
(input) INTEGER
*
Dimension of the state vector.
*
* NU
(input) INTEGER
*
Dimension of the input vector.
*
* NP
(input) INTEGER
*
Dimension of the parameter vector.
*
* X
(input) DOUBLE PRECISION array, dimension (NX)
*
Array containing the state variables.
*
* U
(input) DOUBLE PRECISION array, dimension (NU)
*
Array containing the input values.
*
* P
(input) DOUBLE PRECISION array, dimension (NP)
*
Array containing the parameter values.
*
* T
(input) INTEGER
*
The time point where the derivative is evaluated.
*
* FX
(output) DOUBLE PRECISION array, dimension (LDFX,NX)
*
The array with the resulting Jacobian matrix.
*
* LDFX (input)
INTEGER
*
The leading dimension of the array FX.
*
* Error Indicator
*
* INFO INTEGER
*
Returns values of error from JEVAL or 100 in case
*
a bad problem was choosen.
*
* METHOD
*
* REFERENCES
*
* CONTRIBUTORS
*
* REVISIONS
*
* -
*
* KEYWORDS
*
*
* ******************************************************************
*
* .. Common variables ..
COMMON /TESTING/ ISOLVER
INTEGER LSODI_, LSOIBT_, RADAU5_,
DASSL_, DASPK_, GELDA_
PARAMETER (LSODI_ = 1, LSOIBT_
= 2)
PARAMETER (RADAU5_ = 3, DASSL_
= 4, DASPK_ = 5)
PARAMETER (GELDA_ = 6)
* .. Scalar Arguments ..
INTEGER
NRP, NIP, NX, NU, NP, LDFX, INFO
DOUBLE PRECISION T
* .. Array Arguments ..
INTEGER
IPAR(NIP)
DOUBLE PRECISION X(NX), U(NU), P(NP),
RPAR(NRP), FX(LDFX,NX)
* .. External Subroutines ..
EXTERNAL
JLSODIX
* .. Executable Statements ..
CALL JLSODIX(LDFX,NX,T,X,X,FX,INFO,RPAR,IPAR)
* *** Last line of JACFX_ ***
END
SUBROUTINE JACFXDOT_( NRP, NIP, RPAR,
IPAR,
$
NX, NU, NP, XPRIME, U, P, T, J, LDJ, INFO )
*
* Copyright (c) 2002-2017 NICONET e.V.
*
* PURPOSE
*
* MATJACFXDOT routine for TRANSAMP problem
*
* ARGUMENTS
*
* Input/Output Parameters
*
* NRP (input)
INTEGER
*
Dimension of RPAR array.
*
* NIP (input)
INTEGER
*
Dimension of IPAR array.
*
* RPAR (input/output)
DOUBLE PRECISION array
*
Array for communication with the driver.
*
* IPAR (input/output)
INTEGER array
*
Array for communication with the driver.
*
* NX
(input) INTEGER
*
Dimension of the state vector.
*
*
* NU
(input) INTEGER
*
Dimension of the input vector.
*
* NP
(input) INTEGER
*
Dimension of the parameter vector.
*
* XPRIME (input) DOUBLE PRECISION
array, dimension (NX)
*
Array containing the derivative of the state variables.
*
* U
(input) DOUBLE PRECISION array, dimension (NU)
*
Array containing the input values.
*
* P
(input) DOUBLE PRECISION array, dimension (NP)
*
Array containing the parameter values.
*
* T
(input) INTEGER
*
The time point where the derivative is evaluated.
*
* J
(output) DOUBLE PRECISION array, dimension (LDJ,NX)
*
The array with the resulting derivative matrix.
*
* LDJ (input)
INTEGER
*
The leading dimension of the array J.
*
* Error Indicator
*
* INFO INTEGER
*
Returns 1 in case a bad problem was choosen.
*
* METHOD
*
* REFERENCES
*
* CONTRIBUTORS
*
* REVISIONS
*
* -
*
* KEYWORDS
*
*
* ******************************************************************
*
* .. Common variables ..
COMMON /TESTING/ ISOLVER
INTEGER LSODI_, LSOIBT_, RADAU5_,
DASSL_, DASPK_, GELDA_
PARAMETER (LSODI_ = 1, LSOIBT_
= 2)
PARAMETER (RADAU5_ = 3, DASSL_
= 4, DASPK_ = 5)
PARAMETER (GELDA_ = 6)
* .. Scalar Arguments ..
INTEGER
NRP, NIP, NX, NU, NP, LDJ, INFO
DOUBLE PRECISION T
* .. Array Arguments ..
INTEGER
IPAR(NIP)
DOUBLE PRECISION XPRIME(NX), U(NU),
P(NP), RPAR(NRP), J(LDJ,NX)
* .. Executable Statements ..
*
CALL JDOTLSODIX(LDJ,NX,T,XPRIME,XPRIME,J,INFO,RPAR,IPAR)
ENDIF
* *** Last line of JACFXDOT_ ***
END
Program Data
No data requiredProgram Results
DAESOLVER EXAMPLE PROGRAM RESULTS Problem: LSODIX Solver: LSODI lsodix (lsodix , DAE) with SOLVER 1 IWARN on exit from DAESolver = 2 Solution: (calculated) (reference) 6.462112224297606E-07 1.255974374648338E-10 6.117680951077711E-07 Relative error comparing with the reference solution: .8898590685949503E-06