c----------------------------------------------------------------------- c c This file is part of the Test Set for IVP solvers c http://www.dm.uniba.it/~testset/ c c Problem VAN DER POL c ODE of dimension 2 c c----------------------------------------------------------------------- subroutine prob(problm, + neqn,ndisc,t, + ijac,mljac,mujac) character*(*) problm integer neqn,ndisc,mljac,mujac,ijac double precision t(0:*) problm = 'vdpol' neqn = 2 ndisc = 0 t(0) = 0d0 t(1) = 2.0d0 ijac = 1 mljac = neqn mujac = neqn return end c----------------------------------------------------------------------- subroutine init(neqn,y) integer neqn double precision t,y(neqn) y(1) = 2d0 y(2) = 0d0 return end c----------------------------------------------------------------------- subroutine feval(neqn,t,y,f,ierr,rpar,ipar) integer neqn,ierr,ipar(*) double precision t,y(neqn),f(neqn),rpar(*) f(1) = y(2) f(2) = ((1-y(1)**2)*y(2)-y(1))/1.0d-6 return end c----------------------------------------------------------------------- subroutine jeval(neqn,t,y,dfdy,ldim,ierr,rpar,ipar) integer ldim,neqn,ierr,ipar(*) double precision t,y(neqn),dfdy(ldim,neqn),rpar(*) integer i,j dfdy(1,1) = 0d0 dfdy(1,2) = 1d0 dfdy(2,1) = (-2.0d0*y(1)*y(2)-1d0)/1.0d-6 dfdy(2,2) = (1d0-y(1)**2)/1.0d-6 return end c----------------------------------------------------------------------- subroutine solut(neqn,y) integer neqn double precision y(neqn) c c computed using double precision RADAU on an c Alphaserver DS20E, with a 667 MHz EV67 processor. c c uround = 1.01d-19 c rtol = atol = h0 = 1.1d-18 c c y(1) = 0.1706167732170483D+01 y(2) = -0.8928097010247975D+00 return end