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 DISCLAIMER: see c http://www.dm.uniba.it/~testset/disclaimer.php c c The most recent version of this source file can be found at c http://www.dm.uniba.it/~testset/src/problems/vdpol.f c c This is revision c $Id: vdpol.F,v 1.2 2006/10/02 10:29:14 testset Exp $ c c----------------------------------------------------------------------- integer function pidate() pidate = 20060828 return end c----------------------------------------------------------------------- subroutine prob(fullnm,problm,type, + neqn,ndisc,t, + numjac,mljac,mujac, + nummas,mlmas,mumas, + ind) character*(*) fullnm, problm, type integer neqn,ndisc,mljac,mujac,mlmas,mumas,ind(*) double precision t(0:*) logical numjac, nummas fullnm = 'Problem VANDERPOL' problm = 'vdpol' type = 'ODE' neqn = 2 ndisc = 0 t(0) = 0d0 t(1) = 2.0d0 numjac = .false. mljac = neqn mujac = neqn return end c----------------------------------------------------------------------- subroutine init(neqn,t,y,yprime,consis) integer neqn double precision t,y(neqn),yprime(neqn) logical consis y(1) = 2d0 y(2) = 0d0 return end c----------------------------------------------------------------------- subroutine settolerances(neqn,rtol,atol,tolvec) integer neqn logical tolvec double precision rtol(neqn), atol(neqn) tolvec = .false. return end c----------------------------------------------------------------------- subroutine setoutput(neqn,solref,printsolout, + nindsol,indsol) logical solref, printsolout integer neqn, nindsol integer indsol(neqn) c the reference solution is available solref = .true. c output file is required printsolout = .true. c default values if printsolout is .true. nindsol = neqn c only nindsol component of indsol are referenced do i=1,nindsol indsol(i) = i end do return end c----------------------------------------------------------------------- subroutine feval(neqn,t,y,yprime,f,ierr,rpar,ipar) integer neqn,ierr,ipar(*) double precision t,y(neqn),yprime(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(ldim,neqn,t,y,yprime,dfdy,ierr,rpar,ipar) integer ldim,neqn,ierr,ipar(*) double precision t,y(neqn),yprime(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 meval(ldim,neqn,t,y,yprime,dfddy,ierr,rpar,ipar) integer ldim,neqn,ierr,ipar(*) double precision t,y(neqn),yprime(neqn),dfddy(ldim,neqn),rpar(*) c c dummy subroutine c return end c----------------------------------------------------------------------- subroutine solut(neqn,t,y) integer neqn double precision t,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