######################################################################

ODE :=
`See ?directionfield, ?phaseplot, ?impeuler, ?rungekutta, ?rungekuttahf`:
print(`[ directionfield, phaseplot,impeuler,rungekutta,rungekuttahf ]`);
print(`directionfield is the new name for fieldplot in Release 2.`);
#
# The user routines in this package are
#
#  directionfield:  create a direction field with integral curves
#                   for a first order ODE
#  phaseplot:       create a phase plot with flow lines for two first
#                   order ODE's
#  impeuler:	    solves a first order ODE using the "improved Euler method"
#  rungekutta:	    numerical solution to a system of first order 
#                   ODEs using a 4'th order Runge-Kutta method.  
#  rungekuttahf     is a hardware floating point version of rungekutta.
#
# There is a Maple V help file included for each of these routines.
# After reading in the procedures, to get help for phaseplot type ?phaseplot
#
# April 1993.
# Renamed fieldplot routine to directionfield. Added simple step size 
# adjustments for Runge/Kutta solver built into directionfield.  Jumpy 
# plots are thus less likely.  
# Rewrote rungekuttahf to handle up to 50 equations.  Routine is now 
# much faster for complicated equations.
#
# January 1993.
# Fixed a bug which caused display to break when grid =[0,0] was used as
# an option.  Thanks to Robert Lopez for pointing this out and 
# helping me track it down.  Fixed rungekuttahf to handle more
# types of functions. In particular rungekuttahf now handles
# functions with if statements.
#
# June 1992.
# Fixed arrow routines which would not scale arrows in extreme cases.
# Changed (improved?) the shape of the arrows.
# Added type checking to improve the error messages the user would
# see when the parameters were incorrectly input, thanks to Michael
# Monagan for his help on this.
#
# November 1991.
# Added help files for rungekutta, rungekuttahf, firsteuler, impeuler.
# Fixed a bug in rungekuttahf which was restricting accuracy to current
# Digits setting, rather than the accuracy of evalhf.
#
# April 1991.
# This file contains Maple procedures which work with Maple V. They are
# used to analyze systems of first order differential equations.
#
# Hard copy documentation for the procedures is available as chapter 5
# of a book published by Brooks/Cole Publishing Company called Maple V
# Flight Manual.
#
# Thanks to George Labahn of the Waterloo group for sharing his code,
# parts of which helped to speed the following code. Also, the technique
# of putting arrows on vectors is his. Michael Monagan has been very
# helpful testing the routines and making suggestions for general
# improvements. Others in the Waterloo group have also been helpful for
# testing and suggestions.
#
# Daniel Schwalbe
# 6635 County Road 101
# Corcoran, MN 55340
#
# e-mail dschwalb@daisy.uwaterloo.ca
#
#

`help/text/directionfield` := TEXT(
`FUNCTION: directionfield - create a direction field with flow lines `,
`for first order differential equation.`,
`   `,
`CALLING SEQUENCE:`,
` directionfield(f, h, v)`,
` directionfield(f, h, v,{inits},...)`,
`   `,
`PARAMETERS:`,
` f - function of two variables where D(y)(t)=f(t,y).`,
` h - horizontal range`,
` v - vertical range`,
` {inits} - set of initial values for flow lines.`,
`   `,
`SYNOPSIS: - A typical call to the directionfield function is`,
`directionfield(f,a..b,c..d), where f is a real function of two`,
`variables and a..b,c..d specifies the horizontal and vertical real`,
`range. For each point, (tp,yp), of a grid a short line segment is`,
`drawn with slope f(tp,yp).`,
`   `,
`- Each initial point is specified as a list, [t0,y0]. A fourth`,
`order Runge/Kutta scheme is used to calculate the solution to`,
`D(y)(t)=f(t,y) through the point (t0,y0) over the entire`,
`horizontal range. The routine will adjust the step size`,
`down for very steep solutions.`,
`   `,
`- Remaining arguments are interpreted as options which are`,
`specified as equations of the form option = value. There are`,
`currently 3 options -- grid, stepsize, iterations`,
`The value for grid is a list, [m,n], where the current default is`,
`[16,12]. The value for stepsize controls the step size in the`,
`Runge/Kutta numerical scheme and determines how many points will`,
`be plotted, the current default is (b-a)/20 and 20 points are`,
`plotted. It is possible to decrease the stepsize without`,
`increasing the number of points plotted (which slows down the`,
`routine and leads to unreasonably large plot structures) by`,
`setting the value for iterations. iterations = 10 decreases the`,
`stepsize by a factor of 10, and every 10th point will be stored and`,
`plotted. `,
`   `,
`   `,
`EXAMPLES:`,
`directionfield((t,y)->t-y,-4..4,-3..3);`,
`directionfield((t,y)->y,0..1,0..3,{[0,1]});`,
`directionfield((t,y)->-t/y,-4..4,-4..4,{[0,1],[0,2],[0,3],[0,4]});`,
`directionfield((t,y)->sin(t+y),-4..4,-4..4,grid=[16,14]);`,
`   `,
`eq:=(t,y)->sin(t+y):`,
`inits:={seq(seq([2*i,2*j],i=-2..2),j=-2..2)}:`,
`directionfield(eq,-4..4,-4..4,inits,grid=[0,0]);`,
`   `,
`directionfield((t,y)->-y^2+t,0..20,0..10,[0,1]);`,
`directionfield((t,y)->-y^2+t,0..20,0..10,[0,1],iterations=10);`,
`   `,
`eq4:=(x,y)->(1-x^2-y^2)/(y-x+2):`,
`directionfield(eq4,-2..2,-2..2,grid=[15,15]);`,
`   `,
`eq4:=(x,y)->(1-x^2-y^2)/(y-x+2):`,
`init4:=seq(seq([j,i],i=-2..2),j=-2..2):`,
`directionfield(eq4,-2..2,-2..2,{init4});`,
`   `,
`   `,
`eq5:=(t,y)->y*tan(t);`,
`init5:=seq(seq([i,j],i=-2..2),j=-2..2):`,
`directionfield(eq5,-6..6,-6..6,grid=[0,0],{init5},`,
` iterations=5,stepsize=0.1);`,
`   `,
`eq6:=(t,y)->-y*tan(t);`,
`init6:=seq(seq([2*i,2*j],i=-2..2),j=-2..2):`,
`directionfield(eq6,-6..6,-6..6,grid=[0,0],{init6},`,
` iterations=5,stepsize=.1);`,
` `,
`field:=directionfield((t,x)->sin(t*x),0..10,0..10,`,
`{seq([0,i],i=0..10)},stepsize=0.1,grid=[25,25]):`,
`isoclines:=plot({seq(i*Pi/t,i=1..10)},t=0..10,0..10):`,
`plots[display]({field,isoclines});`,
` `,
`see also: firsteuler, impeuler, rungekutta, rungekuttahf, phaseplot`):

`help/text/phaseplot` := TEXT(
`FUNCTION: phaseplot - draw phase space with flow lines for second`,
` order autonomous differential equation.`,
` `,
`CALLING SEQUENCE:`,
` phaseplot([f1,f2], h, v)`,
` phaseplot([f1,f2], h, v,{inits},...)`,
`   `,
`PARAMETERS: f1,f2 - two functions of three variables where`,
`D(x)(t)=f1(t,x,y) and D(y)(t)=f2(t,x,y) where normally f1 and f2`,
`do not depend on t. t=0 is substituted to draw the phase space`,
`otherwise. Note it will still be interesting to draw flow`,
`lines in the non-autonomous case, i.e. where f1 and f2 depend`,
`explicity on time. h - horizontal range v - vertical range {inits}`,
`- set of initial values for flow lines.`,
`   `,
`SYNOPSIS: - A typical call to the phaseplot function is`,
`phaseplot([f1,f2],a..b,c..d), where f1 and f2 are real functions of`,
`three variables and a..b,c..d specifies the horizontal and`,
`vertical real range. For each point, (xp,yp), of a grid, a scaled`,
`vector with basepoint (xp,yp) and endpoint`,
`(f1(0,xp,yp),f2(0,xp,yp)) is drawn.`,
`   `,
`- Each initial point is specified as a list, [t0,x0,y0]. A`,
`fourth order Runge/Kutta scheme is used to calculate the solution`,
`to the system of equations, D(x)(t)=f1(t,x,y),`,
`D(y)(t)=f2(t,x,y),x(t0)=x0,y(t0)=y0. The default is to use a step`,
`size of 0.1 and do 50 iterations.`,
`   `,
`- Remaining arguments are interpreted as options which are`,
`specified as equations of the form option = value. There are`,
`currently 4 options -- grid, stepsize, numsteps, and iterations.`,
`The value for grid is a list, [m,n], where the current default is`,
`[16,14]. The value for stepsize controls the step size in the`,
`Runge/Kutta numerical scheme, the current default is 0.1. The`,
`value for numsteps specifies how many steps forward and backward`,
` will be calculated and stored by the Runge/Kutta routine. The`,
`current default is 50 forward and back. It is possible to`,
`decrease the stepsize without increasing the number of points`,
`plotted (which slows down the routine and leads to unreasonably`,
`large plot structures) by setting the value for iterations.`,
`iterations = 10, would decrease the stepsize by a factor of 10`,
`and every 10th point would be stored and plotted.`,
` `,
`   `,
`EXAMPLES: `,
`phaseplot([(t,x,y)->-x,(t,x,y)->y],-4..4,-4..4);`,
`   `,
`eq1:=(t,x,y)-> 5*x-x^2-x*y: `,
`eq2:=(t,x,y)->x*y-2*y:`,
`inits:={seq([0,5,i],i=1..9),seq([0,i,0.1],i=1..9)}:`,
`phaseplot([eq1,eq2],0..10,0..10,inits,grid=[15,15]);`,
`   `,
`readlib(Heaviside): `,
`eq3:=(t,y,v)->v: `,
`eq4:=(t,y,v)->-Heaviside(y):`,
`phaseplot([eq3,eq4],-4..4,-4..4); `,
`   `,
`readlib(Heaviside): `,
`eq5:=(t,y,v)->v:`,
`eq6:=(t,y,v)->-Heaviside(y-5)*v/2-(1-Heaviside(y-5))*v/10+1:`,
`phaseplot([eq5,eq6],0..10,0..5,[0,0,0],iterations=10,`,
` numsteps=100,grid=[15,15]);`,
`   `,
`eq7:=(t,x,y)-> -x+x*y: `,
`eq8:=(t,x,y)->-y+2*x*y:`,
`inits2:={seq(seq([0,i/2,j/2],i=0..2),j=0..2)}:`,
`phaseplot([eq7,eq8],-.5..1.5,-.5..1.5,inits2,grid=[15,15]);`,
`   `,
`pend3:=(t,y,v) -> v;`,
`pend4:=(t,y,v) -> -sin(y) -0.2*v;`,
`pendinits:=seq([0,j/2],j=2..8):`,
`phaseplot([pend3,pend4],-10..10,-4..4,{pendinits},`,
` numsteps=200);`,
`   `,
`ropedv:= proc(t,v,y) if v > 0 then (9800 -3*980*y - 3*v^2)/(3*y)`,
`   else (9800 -3*980*y)/(3*y) fi end;`,
`ropedy:= (t,v,y) -> v;`,
`phaseplot([ropedv,ropedy],-59..59,0.1..5);`,
`   `,
`eq11:=(t,y,v)->v^2-y^2;`,
`eq12:=(t,y,v)->y-2*v;`,
`eq1init:=seq(seq([i,j],j=-2..2),i=-2..2):`,
`phaseplot([eq11,eq12],-5..5,-5..5,grid=[12,12],{eq1init});`,
`   `,
`eq21:=(t,x,y)->cos(x+y);`,
`eq22:=(t,x,y)->-sin(y);`,
`phaseplot([eq21,eq22],-6..6,-6..6,grid=[25,25]);`,
`   `,
`see also: firsteuler, impeuler, rungekutta, rungekuttahf, directionfield`):

`help/text/rungekutta` := TEXT(
`FUNCTION: rungekutta - approximate a solution to a system of first`,
`order ordinary differential equations using a fourth order Runge/Kutta`,
`numerical method. Solution is returned as an array of lists.`,
`rungekuttahf is an alternative form of this procedure which uses`,
`hardware floating point operations and is not affected by the Digits`,
`settings. rungekuttahf is limited to double precision on most systems`,
`and will be about as accurate as rungekutta using a Digits setting `,
`of 15. It is however, MUCH faster.`,
`   `,
`CALLING SEQUENCE:`,
` rungekutta(fcnlist,initvalue,stepsize,numsteps)`,
` rungekutta(fcnlist,initvalue,stepsize,numsteps,iterations)`,
`   `,
`PARAMETERS:`,
` fcnlist - list of functions, [f1,f2,...,fn], of n+1 real variables.  `,
`           The differential equation to be solved is of the form:`,
`				dx1/dt = f1(t,x1,x2,...,xn)`,
`				dx2/dt = f2(t,x1,x2,...,xn)`,
`				   .   .          .`,
`				   .   .          .`,
`				   .   .          .`,
`				dxn/dt = fn(t,x1,x2,...,xn)`,
` initvalue - list of n+1 numbers, [i0,i1,i2,...,in],  which represents `,
`             the initial value of the differential equation. `,
` stepsize - step size used by the numerical method.`,
` numsteps - number of steps taken by numerical method.`,
` iterations - optional value which decreases the step size without `,
`              storing more points.  New step size will be `,
`			  stepsize/iterations.`,
`   `,
`SYNOPSIS: - A typical call to the rungekutta function is`,
`rungekutta([f1,f2,...,fn],[i0,i1,i2,...,in],stepsize,numsteps), where`,
`f1,f2,...,fn are real valued functions of n+1 variables, stepsize is`,
`the step size used by the numerical method, and numsteps specifies the`,
`number of steps to be taken by the numerical method. An array of (n+1)`,
`lists of (numsteps+1) numbers is returned, where each list represents`,
`a point calculated by the numerical method.`,
`   `,
`- An optional fifth argument is used to reduce the step size without `,
`storing more points.  This is especially useful when graphing `,
`solutions which require a small step size and a large plot structure `,
`is undesirable.  Also useful when only the final value is wanted.`,
`   `,
`- makelist is a procedure in the ODE file to extract a list of values`,
`from an array of lists. This is useful for making two dimensional`,
`plots of solutions calculated by rungekutta.  makelist(A,m,n) will `,
`make a list of values from the mth and nth components from each `,
`entry of A. m defaults to 1 and n defaults to 2 if no values are `,
`entered.`,
`   `,
`EXAMPLES: `,
`rkpts1:=rungekutta((t,y)->t-y,[0,1],.1,10);`,
`plot(makelist(rkpts1));`,
`   `,
`eq1:=(t,x,y)-> y;`,
`eq2:=(t,x,y)->-x;`,
`rkpts2:=rungekutta([eq1,eq2],[0,1,1],.1,60);`,
`plot(makelist(rkpts2,2,3));`,
`   `,
`dp1:=(t,x,xp,y,yp)->xp;`,
`dp2:=(t,x,xp,y,yp)->-2*sin(x)+sin(y);`,
`dp3:=(t,x,xp,y,yp)->yp;`,
`dp4:=(t,x,xp,y,yp)->-2*sin(y)+2*sin(x);`,
`dpinit:=[0,0,0,0,2];`,
`dppts:=rungekutta([dp1,dp2,dp3,dp4],dpinit,1,10,1):`,
`dppts2:=rungekutta([dp1,dp2,dp3,dp4],dpinit,1,10,2):`,
`dppts3:=rungekutta([dp1,dp2,dp3,dp4],dpinit,1,10,4):`,
`plot({makelist(dppts,2,4),makelist(dppts2,2,4),`,
`  makelist(dppts3,2,4)});`,
`   `,
`see also: firsteuler, impeuler, rungekuttahf, directionfield, phaseplot`):

`help/text/rungekuttahf` := TEXT(
`FUNCTION: rungekuttahf - approximate a solution to a system of up to`,
`fifty first order ordinary differential equations using a fourth order`,
`Runge/Kutta numerical method. Solution is returned as an array of`,
`lists. rungekutta is an alternative form of this procedure which uses`,
`Maple's arbritrary precision floating point operations and is thus`,
`affected by the Digits settings. rungekutta is useful for studying the`,
`effect of the round off error that occurs with numerical methods or`,
`with differential equations for which rungekuttahf fails. rungekuttahf`,
`is limited to double precision on most systems and will be about as`,
`accurate as rungekutta using a Digits setting of 15. rungekuttahf is`,
`however, MUCH faster.`,
`   `,
`CALLING SEQUENCE:`,
` rungekuttahf(fcnlist,initvalue,stepsize,numsteps)`,
` rungekuttahf(fcnlist,initvalue,stepsize,numsteps,iterations)`,
`   `,
`PARAMETERS:`,
` fcnlist - list of functions, [f1,f2,...,fn], of n+1 real variables.  `,
`           The differential equation to be solved is of the form:`,
`				dx1/dt = f1(t,x1,x2,...,xn)`,
`				dx2/dt = f2(t,x1,x2,...,xn)`,
`				   .   .          .`,
`				   .   .          .`,
`				   .   .          .`,
`				dxn/dt = fn(t,x1,x2,...,xn)`,
` initvalue - list of n+1 numbers, [i0,i1,i2,...,in],  which represents `,
`             the initial value of the differential equation. `,
` stepsize - step size used by the numerical method.`,
` numsteps - number of steps taken by numerical method.`,
` iterations - optional value which decreases the step size without `,
`              storing more points.  New step size will be `,
`			  stepsize/iterations.`,
`   `,
`SYNOPSIS: - A typical call to the rungekuttahf function is`,
`rungekuttahf([f1,f2,...,fn],[i0,i1,i2,...,in],stepsize,numsteps),`,
`where f1,f2,...,fn are real valued functions of n+1 variables,`,
`stepsize is the step size used by the numerical method, and numsteps`,
`specifies the number of steps to be taken by the numerical method.`,
`Alternatively, rungekuttahf(f,[i0,i1,i2,...,in],stepsize,numsteps)`,
`where f is a vector valued function also works. An array of`,
`(n+1) lists of (numsteps+1) numbers is returned, where each list`,
`represents a point calculated by the numerical method.`,
`   `,
`- An optional fifth argument is used to reduce the step size without `,
`storing more points.  This is especially useful when graphing `,
`solutions which require a small step size and a large plot structure `,
`is undesirable.  Also useful when only the final value is wanted.`,
`   `,
`- makelist is a procedure in the ODE file to extract a list of values`,
`from an array of lists. This is useful for making two dimensional`,
`plots of solutions calculated by rungekuttahf.  makelist(A,m,n) will `,
`make a list of values from the mth and nth components from each `,
`entry of A. m defaults to 1 and n defaults to 2 if no values are `,
`entered.`,
`   `,
`EXAMPLES: `,
`rkpts1:=rungekuttahf((t,y)->t-y,[0,1],.1,10);`,
`plot(makelist(rkpts1));`,
`   `,
`eq1:=(t,x,y)-> y;`,
`eq2:=(t,x,y)->-x;`,
`rkpts2:=rungekuttahf([eq1,eq2],[0,1,1],.1,60):`,
`plot(makelist(rkpts2,2,3));`,
`   `,
`dp:=(t,x,xp,y,yp)->[xp,-2*sin(x)+sin(y),yp,-2*sin(y)+2*sin(x)];`,
`dpinit:=[0,0,0,0,2];`,
`dppts:=rungekuttahf(dp,dpinit,.1,100):`,
`plot({makelist(dppts,2,4)});`,
`   `,
`eq3:=(t,x,y)->x*sin(t)+5*y;`,
`eq4:=(t,x,y)->-10*x+y*sin(t);`,
`rkpts3:=rungekuttahf([eq3,eq4],[0,2,0],.05,200,5):`,
`plot(makelist(rkpts3,2,3));`,
`plot({makelist(rkpts3,1,2),makelist(rkpts3,1,3)});`,
`#for the Maple expert and the brave alike, construct a PLOT3D structure`,
`plot1:=PLOT3D(CURVELIST([seq([rkpts3[i][1],rkpts3[i][2],rkpts3[i][3]],`,
`    i=0..100)]),AXES(BOXED)):";`,
`	`,
`pulse:=proc(t) if not type(t,numeric) then 'pulse'(t)`,
`   elif trunc(trunc(t)/2)=trunc(t)/2 then 1 else 0 fi end;`,
`pulsepts:=rungekuttahf((t,x)->pulse(t)*x,[0,1],.1,100):`,
`plot(makelist(pulsepts));`,
`	`,
`   `,
`see also: firsteuler, impeuler, rungekutta, directionfield, phaseplot`):

`help/text/firsteuler` := TEXT(
`FUNCTION: firsteuler - approximate a solution to a first order`,
`ordinary differential equations using Euler's method.`,
`Solution is returned as an array of lists. `,
`   `,
`CALLING SEQUENCE:`,
` firsteuler(f,i,h,n)`,
`   `,
`PARAMETERS:`,
` f - function of 2 real variables.  The differential equation to be `,
`     solved is dy/dt = f(t,y)  `,
` i - list of 2 numbers. `,
` h - step size used by the numerical method.`,
` n - number of steps taken by numerical method.`,
`   `,
`SYNOPSIS: - A typical call to the firsteuler function is`,
`firsteuler(f,i,h,n) where f is a real valued functions of 2 variables,`,
`h is the step size used by the numerical method, and n specifies the`,
`number of steps to be taken by the numerical method. An array of (n+1)`,
`lists of 2 numbers is returned, where each list represents a point`,
`calculated by the numerical method.`,
`   `,
`- makelist is a procedure in the ODE file to extract a list of values`,
`from an array of lists. This is useful for making a plot of the`,
`approximate solution calculated by firsteuler. makelist(A) will make a`,
`list of values from the entries of A. `,
`   `,
`EXAMPLES: `,
`eulerpts:=firsteuler((t,y)->t-y,[0,1],.1,10);`,
`plot(makelist(eulerpts));`,
`   `,
`eq:=(t,y)->y;`,
`pts1:=firsteuler(eq,[0,1],.5,2);`,
`pts2:=firsteuler(eq,[0,1],.25,4);`,
`pts3:=firsteuler(eq,[0,1],.1,10);`,
`plot({makelist(pts1),makelist(pts2),makelist(pts3),`,
`      exp(t)},0..1,style=LINE);`,
`	  `,
`eq2:=(t,y)->-t/y;`,
`eq2pts:=firsteuler(eq2,[0,4],.2,25):`,
`plot(makelist(eq2pts));`,
`   `,
`eq3:=(t,y)->5*y-6*exp(-t);`,
`eq3pts:=firsteuler(eq3,[0,1],.1,10);`,
`eq3sol:=dsolve({diff(y(t),t)=eq3(t,y),y(0)=1},y(t));`,
`plot({makelist(eq3pts),rhs(eq3sol)},0..1);`,
`   `,
`eq4:=(t,y)->25*y*(1-y);`,
`eq4pts:=firsteuler(eq4,[0,1.3],.1,10):`,
`plot(makelist(eq4pts));`,
`   `,
`eq4:=(t,y)->25*y*(1-y);`,
`eq4pts:=firsteuler(eq4,[0,1.3],.1,10):`,
`eq4pts2:=firsteuler(eq4,[0,1.3],.05,20):`,
`eq4pts3:=firsteuler(eq4,[0,1.3],.025,40):`,
`plot({makelist(eq4pts),makelist(eq4pts2),makelist(eq4pts3)});`,
`   `,
`   `,
`see also: impeuler, rungekutta, rungekuttahf, directionfield, phaseplot`):

`help/text/impeuler` := TEXT(
`FUNCTION: impeuler - approximate a solution to a first order ordinary`,
`differential equations using the improved Euler method. Solution is`,
`returned as an array of lists. `,
`   `,
`CALLING SEQUENCE:`,
` impeuler(f,i,h,n)`,
`   `,
`PARAMETERS:`,
` f - function of 2 real variables.  The differential equation to be `,
`     solved is dy/dt = f(t,y)  `,
` i - list of 2 numbers. `,
` h - step size used by the numerical method.`,
` n - number of steps taken by numerical method.`,
`   `,
`SYNOPSIS: - A typical call to the impeuler function is`,
`impeuler(f,i,h,n) where f is a real valued functions of 2 variables,`,
`h is the step size used by the numerical method, and n specifies the`,
`number of steps to be taken by the numerical method. An array of (n+1)`,
`lists of 2 numbers is returned, where each list represents a point`,
`calculated by the numerical method.`,
`   `,
`- makelist is a procedure in the ODE file to extract a list of values`,
`from an array of lists. This is useful for making a plot of the`,
`approximate solution calculated by impeuler. makelist(A) will make a`,
`list of values from the entries of A. `,
`   `,
`EXAMPLES: `,
`impts:=impeuler((t,y)->t-y,[0,1],.1,10);`,
`plot(makelist(impts));`,
`   `,
`eq:=(t,y)->y;`,
`impts1:=impeuler(eq,[0,1],.5,2);`,
`impts2:=impeuler(eq,[0,1],.25,4);`,
`impts3:=impeuler(eq,[0,1],.1,10);`,
`plot({makelist(impts1),makelist(impts2),makelist(impts3),`,
`      exp(t)},0..1,style=LINE);`,
`	  `,
`eq2:=(t,y)->-t/y;`,
`eq2impts:=impeuler(eq2,[0,4],.2,25):`,
`plot(makelist(eq2impts));`,
`   `,
`eq3:=(t,y)->5*y-6*exp(-t);`,
`eq3impts:=impeuler(eq3,[0,1],.1,10);`,
`eq3sol:=dsolve({diff(y(t),t)=eq3(t,y),y(0)=1},y(t));`,
`plot({makelist(eq3impts),rhs(eq3sol)},0..1);`,
`   `,
`eq4:=(t,y)->25*y*(1-y);`,
`eq4impts:=impeuler(eq4,[0,1.3],.1,10):`,
`plot(makelist(eq4impts));`,
`   `,
`eq4:=(t,y)->25*y*(1-y);`,
`eq4impts:=impeuler(eq4,[0,1.3],.1,10):`,
`eq4impts2:=impeuler(eq4,[0,1.3],.05,20):`,
`eq4impts3:=impeuler(eq4,[0,1.3],.025,40):`,
`eq4implot:=plot({makelist(eq4impts),makelist(eq4impts2),`,
`              makelist(eq4impts3)}):`,
`eq4field:=directionfield(eq4,0..1,0.8..1.5):`,
`with(plots):`,
`display({eq4implot,eq4field});`,
`   `,
`   `,
`see also: firsteuler, rungekutta, rungekuttahf, directionfield, phaseplot`):


#################################################

firsteuler:= proc(f,initlist,h0,n,iterations)
local i,xk,yk,pts,h,npts;
options `Copyright 1991 by Daniel Schwalbe`;
pts:=array(0..n);
if nargs>4 then
   npts:=iterations;
else npts:=1 fi;
h:=evalf(h0/npts);
xk:=evalf(initlist[1]); yk:=evalf(initlist[2]);
pts[0] := [xk,yk];
for i to n while type(yk,numeric) do
  to npts while type(yk,numeric) do
    yk:=yk + h*f(xk,yk);
    xk:=xk+h od;
  pts[i]:=[xk,yk] od:
if i <= n then ERROR(`undefined procedure`)
else op(pts) fi:
end:

##################################################

impeuler:=proc(f,initlist,h0,n,iterations)
local h,xk,yk,pts,j,npts:
options `Copyright 1991 by Daniel Schwalbe`;
pts:=array(0..n);
if nargs>4 then
  npts:=iterations
  else npts:=1 fi;
h:=evalf(h0/npts):
xk:=evalf(initlist[1]); yk:=evalf(initlist[2]);
pts[0] := [xk,yk];
for j to n while type(yk,numeric) do
  to npts while type(yk,numeric) do
    yk:=yk+h/2*(f(xk,yk)+f(xk+h,yk+h*f(xk,yk))):
    xk:=xk+h od:
  pts[j]:=[xk,yk] od:
if j <= n then ERROR(`undefined procedure`)
 else op(pts) fi:
end:
##################################################
rungekutta:=proc(f,initlist,h0,n,iterations)
local k1,k2,k3,k4,fl,h,pts,j,i,k,m,npts,temp:
options `Copyright 1991 by Daniel Schwalbe`;
m:=nops([op(initlist)])-1;
pts:=array(0..n);

if type(f,list) then fl:=f
  elif nops([op(1,eval(f))]) > 2 then
    fl:=[seq(unapply(f(op(1,eval(f)))[i],op(1,eval(f))),i=1..m)]
  else fl:=[eval(f)] fi;
if nops([op(1,eval(fl[1]))]) <> m+1 then 
  ERROR(`Initial condition does not match number of unknowns`) fi;

if nargs > 4 then
  npts:=iterations
  else npts:=1 fi;
h:=evalf(h0/npts):
pts[0] := map(evalf,[op(initlist)]);
temp:=pts[0];
for j from 0 to n-1 while
  not member(false,map(type,(temp,numeric))) do
    to npts while not member(false,map(type,(temp,numeric))) do
      for i to m do k1[i]:=
         eval(evalf(fl[i](op(temp)))) od:
      for i to m do k2[i]:= eval(evalf(fl[i](temp[1]+h/2,
         seq(temp[k+1] + h*k1[k]/2,k=1..m)))) od:
      for i to m do k3[i]:=eval(evalf(fl[i](temp[1]+h/2,
         seq(temp[k+1] + h*k2[k]/2,k=1..m)))) od:
      for i to m do k4[i]:=eval(evalf(fl[i](pts[j][1]+h,
         seq(temp[k+1] + h*k3[k],k=1..m)))) od:
    temp:=[temp[1]+h, seq(temp[i+1]+
      h*(k1[i]+2*k2[i]+2*k3[i]+k4[i])/6,i=1..m)] od:
  pts[j+1]:= temp od:
if j<=n-1 then ERROR(`undefined procedure`)
 else op(pts) fi:
end:


##################################################

rkhf:=proc(f,B,CC)
local j,i,h,n,npts,t,k1,k2,k3,k4,m,y:
options `Copyright 1993 by Daniel Schwalbe`;
h:=CC[1]; n:=CC[2]; npts:=CC[3]; m:=CC[4];
if m<=20 then t:=array(1..m+1) else t:=array(1..51) fi;
if m<=20 then y:=array(1..m+1) else y:=array(1..51) fi;
for i to m+1 do t[i]:=B[0,i] od;
if m<=20 then
  k1:=array(1..m);k2:=array(1..m);
  k3:=array(1..m);k4:=array(1..m);
else
  k1:=array(1..50);k2:=array(1..50);
  k3:=array(1..50);k4:=array(1..50);
fi;

for j from 0 to n-1  do
  to npts  do
    f(t,k1):
	y[1]:=t[1]+h/2;
	for i from 2 to m+1 do y[i]:=t[i]+h*k1[i-1]/2 od;
    f(y,k2):
	for i from 2 to m+1 do y[i]:=t[i]+h*k2[i-1]/2 od;
    f(y,k3):
	y[1]:=t[1]+h;
	for i from 2 to m+1 do y[i]:=t[i]+h*k3[i-1] od;
    f(y,k4):
    t[1]:=y[1];
    for i from 2 to m+1 do
      t[i]:= t[i]+h*(k1[i-1]+2*k2[i-1]+2*k3[i-1]+k4[i-1])/6 
  od;od:
  for i to m+1 do B[j+1,i]:= t[i] od:od:
end:


macro(_AA=`plots/ODE/AA`);
macro(_body=`plots/ODE/body`);
rungekuttahf:=proc(f,init,h0,n,iterations)
local h,m,A,B,k,precision,npts,CC,fl,i,j,g,eq;
options `Copyright 1993 by Daniel Schwalbe`;

precision := max(Digits,trunc(evalhf(Digits)));
m:=nops([op(init)])-1;
if m>50 then 
  RETURN(`rungekuttahf is limited to 50 equations.`); fi;
if m<=20 then A:=array(0..n,1..m+1)
  else A:=array(0..n,1..51) fi;
if nargs > 4 then npts:=iterations else npts:=1 fi;
h:=evalf(h0/npts,precision):
for i to m+1 do A[0,i]:=evalf(init[i],precision) od;
CC:=array(1..4);
CC[1]:=h; CC[2]:=n; CC[3]:=npts; CC[4]:=m;

if type(f,list) then 
  if nops([op(1,eval(f[1]))]) <> m+1 then 
    ERROR(`Initial condition does not match number of unknowns`) fi;
  fl:=traperror([seq(eval(f[i])(seq(_AA[j],j=1..m+1)),i=1..m)]);
  if fl = lasterror then
  fl:=[seq(subs(g=eval(f[i]),g(seq(_AA[j],j=1..m+1))),i=1..m)];fi;
else 
  if nops([op(1,eval(f))]) <> m+1 then 
    ERROR(`Initial condition does not match number of unknowns`) fi;
  if m=1 then 
    fl:=traperror(f(_AA[1],_AA[2]) );
    if fl = lasterror then
       fl:= subs(g=eval(f),g(_AA[1],_AA[2])) fi;
  else
    fl:=f(seq(_AA[j],j=1..m+1)) fi;
fi;

if m=1 then
  fl:=subs('_arg'=_AA,'_body'=fl, 
      proc(_arg,eq) 
        eq[1]:=_body;
      end)
elif m=2 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2] 
      end); 
elif m=3 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3] 
      end) 
elif m=4 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
      end)
elif m=5 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
      end) 
elif m=6 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
      end) 
elif m=7 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
      end) 
elif m=8 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
      end) 
elif m=9 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
		eq[9]:=_body[9]; 
      end) 
elif m=10 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
		eq[9]:=_body[9]; 
		eq[10]:=_body[10]; 
      end) 
elif m=11 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
		eq[9]:=_body[9]; 
		eq[10]:=_body[10]; 
		eq[11]:=_body[11]; 
      end) 
elif m=12 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
		eq[9]:=_body[9]; 
		eq[10]:=_body[10]; 
		eq[11]:=_body[11]; 
		eq[12]:=_body[12]; 
      end) 
elif m=13 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
		eq[9]:=_body[9]; 
		eq[10]:=_body[10]; 
		eq[11]:=_body[11]; 
		eq[12]:=_body[12]; 
		eq[13]:=_body[13]; 
      end) 
elif m=14 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
		eq[9]:=_body[9]; 
		eq[10]:=_body[10]; 
		eq[11]:=_body[11]; 
		eq[12]:=_body[12]; 
		eq[13]:=_body[13]; 
		eq[14]:=_body[14]; 
      end) 
elif m=15 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
		eq[9]:=_body[9]; 
		eq[10]:=_body[10]; 
		eq[11]:=_body[11]; 
		eq[12]:=_body[12]; 
		eq[13]:=_body[13]; 
		eq[14]:=_body[14]; 
		eq[15]:=_body[15]; 
      end) 
elif m=16 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
		eq[9]:=_body[9]; 
		eq[10]:=_body[10]; 
		eq[11]:=_body[11]; 
		eq[12]:=_body[12]; 
		eq[13]:=_body[13]; 
		eq[14]:=_body[14]; 
		eq[15]:=_body[15]; 
		eq[16]:=_body[16]; 
      end) 
elif m=17 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
		eq[9]:=_body[9]; 
		eq[10]:=_body[10]; 
		eq[11]:=_body[11]; 
		eq[12]:=_body[12]; 
		eq[13]:=_body[13]; 
		eq[14]:=_body[14]; 
		eq[15]:=_body[15]; 
		eq[16]:=_body[16]; 
		eq[17]:=_body[17]; 
      end) 
elif m=18 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
		eq[9]:=_body[9]; 
		eq[10]:=_body[10]; 
		eq[11]:=_body[11]; 
		eq[12]:=_body[12]; 
		eq[13]:=_body[13]; 
		eq[14]:=_body[14]; 
		eq[15]:=_body[15]; 
		eq[16]:=_body[16]; 
		eq[17]:=_body[17]; 
		eq[18]:=_body[18]; 
      end) 
elif m=19 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
		eq[9]:=_body[9]; 
		eq[10]:=_body[10]; 
		eq[11]:=_body[11]; 
		eq[12]:=_body[12]; 
		eq[13]:=_body[13]; 
		eq[14]:=_body[14]; 
		eq[15]:=_body[15]; 
		eq[16]:=_body[16]; 
		eq[17]:=_body[17]; 
		eq[18]:=_body[18]; 
		eq[19]:=_body[19]; 
      end) 
elif m=20 then
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m), 
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
		eq[9]:=_body[9]; 
		eq[10]:=_body[10]; 
		eq[11]:=_body[11]; 
		eq[12]:=_body[12]; 
		eq[13]:=_body[13]; 
		eq[14]:=_body[14]; 
		eq[15]:=_body[15]; 
		eq[16]:=_body[16]; 
		eq[17]:=_body[17]; 
		eq[18]:=_body[18]; 
		eq[19]:=_body[19]; 
		eq[20]:=_body[20]; 
      end) 
else
  fl:=subs('_arg'=_AA,
   seq(_body[j]=fl[j],j=1..m),
   seq(_body[j]=0,j=m+1..50),
      proc(_arg,eq) 
        eq[1]:=_body[1];
        eq[2]:=_body[2];
		eq[3]:=_body[3]; 
		eq[4]:=_body[4]; 
		eq[5]:=_body[5]; 
		eq[6]:=_body[6]; 
		eq[7]:=_body[7]; 
		eq[8]:=_body[8]; 
		eq[9]:=_body[9]; 
		eq[10]:=_body[10]; 
		eq[11]:=_body[11]; 
		eq[12]:=_body[12]; 
		eq[13]:=_body[13]; 
		eq[14]:=_body[14]; 
		eq[15]:=_body[15]; 
		eq[16]:=_body[16]; 
		eq[17]:=_body[17]; 
		eq[18]:=_body[18]; 
		eq[19]:=_body[19]; 
		eq[20]:=_body[20]; 
		eq[21]:=_body[21]; 
		eq[22]:=_body[22]; 
		eq[23]:=_body[23]; 
		eq[24]:=_body[24]; 
		eq[25]:=_body[25]; 
		eq[26]:=_body[26]; 
		eq[27]:=_body[27]; 
		eq[28]:=_body[28]; 
		eq[29]:=_body[29]; 
		eq[30]:=_body[30]; 
		eq[31]:=_body[31]; 
		eq[32]:=_body[32]; 
		eq[33]:=_body[33]; 
		eq[34]:=_body[34]; 
		eq[35]:=_body[35]; 
		eq[36]:=_body[36]; 
		eq[37]:=_body[37]; 
		eq[38]:=_body[38]; 
		eq[39]:=_body[39]; 
		eq[40]:=_body[40]; 
		eq[41]:=_body[41]; 
		eq[42]:=_body[42]; 
		eq[43]:=_body[43]; 
		eq[44]:=_body[44]; 
		eq[45]:=_body[45]; 
		eq[46]:=_body[46]; 
		eq[47]:=_body[47]; 
		eq[48]:=_body[48]; 
		eq[49]:=_body[49]; 
		eq[50]:=_body[50]; 
      end) 
fi;
 
evalhf(rkhf(fl,var(A),var(CC)));
B:=array(0..n);
for i from 0 to n do B[i]:=[seq(A[i,k],k=1..m+1)] od;
op(B)
end:

macro(_AA=_AA);
macro(_body=_body);


##################################################

makelist:=proc(A,c10,c20)
local c1,c2,ri,rf,i;
options `Copyright 1991 by Daniel Schwalbe`;
if nargs = 1 then
   ri:=op(1,op(2,eval(A)));
   rf:=op(2,op(2,eval(A)));
   c1:=1;
   c2:=2;
else
   c1:=c10;
   c2:=c20;
   ri:=op(1,op(2,eval(A)));
   rf:=op(2,op(2,eval(A))) fi;
[seq([A[i][c1],A[i][c2]],i=ri..rf)];
end:


##################################################
dadaptiverk:=proc(f,CC,h0)
local F,j,tt,Fval,yy,A,p,q,h,t,y;
options `Copyright 1991 by Daniel Schwalbe`;
p:=CC[3]: q:=CC[4]:
tt:=CC[1]:yy:=CC[2]:h:=h0:
A[0]:=[tt,yy]:
F:=dsolve({D(y)(t)=f(t,y(t)),y(tt)=yy},y(t),numeric);
j:=0; Fval:=traperror(F(tt+h));
while Fval<>lasterror and tt+h<=q do
  tt:=tt+h; j:=j+1;
  A[j]:=subs(Fval,[t,y(t)]):
  Fval:=traperror(F(tt+h)) od;
CC[8]:=j:
j:=0;tt:=CC[1]: Fval:=traperror(F(tt-h));
while Fval<>lasterror and tt-h>=p do
  tt:=tt-h; j:=j-1;
  A[j]:=subs(Fval,[t,y(t)]):
  Fval:=traperror(F(tt-h)) od;
CC[7]:=j:
op(A)
end:



fieldrkhf:=proc(A,f,CC,h0,npts,nsteps)
local pt,j,h,p,q,r,s,k1,k2,k3,k4,yold,ynew,oldpt2,oldpt1;
options `Copyright 1993 by Daniel Schwalbe`;
pt:=array(1..2):
h:=h0/npts;
r:=CC[5];s:=CC[6];
p:=CC[3];q:=CC[4];
pt[1]:=CC[1]: pt[2]:=CC[2]:
A[1,0]:=CC[1]; A[2,0]:=CC[2];
j:=0;CC[7]:=1:
yold:=f(pt[1],pt[2]);oldpt2:=pt[2];oldpt1:=pt[1];
while pt[1] <= q and pt[1] >= p and pt[2] < s + 10*(s-r) 
  and pt[2] > r - 10*(s-r) and j<nsteps 
  do
  CC[7]:=j+1:
  j:=j+1;
  to npts  do
    k1:=f(pt[1],pt[2]):
    k2:=f(pt[1]+h/2,pt[2]+h*k1/2):
    k3:=f(pt[1]+h/2,pt[2]+h*k2/2):
    k4:=f(pt[1]+h,pt[2]+h*k3):
    ynew:=(k1+2*k2+2*k3+k4)/6;
	if ((abs(yold)>1) and yold*ynew<0)
          or abs(h*ynew)>(s-r)/20 then 
	  h:=h/2;
	  pt[1]:=oldpt1;
	  pt[2]:=oldpt2;
	else
	 oldpt1:=pt[1];
     pt[1]:=pt[1]+h;
	 oldpt2:=pt[2];
     pt[2]:= pt[2]+h*ynew fi;
  od;
  yold:=ynew;
  A[1,j]:=pt[1]:
  A[2,j]:=pt[2]:
od:

end:


ccomputelineshf:=proc(f,Ax,Ay,CC)
local m,n,dx,dy,k,slp,xinc,yinc,a,b,
 dydx,p,q,r,s,i,j,ldx,ldy;
m:=CC[1]; n:=CC[2];
p:=CC[3]; q:=CC[4];
r:=CC[5]; s:=CC[6];
dx:=(q-p)/m; dy:=(s-r)/n;
dydx:=dy/dx;
ldx:=dx/4; ldy:=dy/4;
a:=p;
for i from 0 to m do
  b:=r;
  for j from 0 to n do
    k:=2*j;
    slp:=f(a,b);
    if slp > dydx or -slp > dydx then
      xinc:=ldy/slp; yinc:=ldy
    else xinc:=ldx; yinc:=ldx*slp fi;
    Ax[i,k]:=a-xinc;
    Ay[i,k]:=b-yinc;
    Ax[i,k+1]:=a+xinc;
    Ay[i,k+1]:=b+yinc;
    b:=b+dy od;
  a:=a+dx od;
end:


ccomputelines:=proc(f,Ax,Ay,p,q,r,s,m,n)
local dx,dy,k,ldx,slp,ldy,xinc,yinc,a,b,i,j;
dx:=(q-p)/m; dy:=(s-r)/n;
ldx:=dx/4; ldy:=dy/4;
a:=p;
for i from 0 to m do
  b:=r;
  for j from 0 to n do
    k:=2*j;
    slp:=traperror(evalf(f(a,b)));
    if lasterror=slp then
      xinc:=0; yinc:=0;
    elif not type(slp,'numeric') then
      xinc:=0; yinc:=0;
    elif slp > dy/dx or -slp > dy/dx then
      xinc:=ldy/slp; yinc:=ldy;
    else xinc:=ldx; yinc:=ldx*slp fi;
    Ax[i,k]:=a-xinc;
    Ay[i,k]:=b-yinc;
    Ax[i,k+1]:=a+xinc;
    Ay[i,k+1]:=b+yinc;
    b:=b+dy od;
  a:=a+dx od;
end:

macro( INITCOND1 = {[constant,constant]} );
macro(_body=`plots/ODE/body`);
directionfield:=proc(F,H,V)
local npts,init,flow,A,Ax,Ay,i,j,nsteps,f,
 CC,h,m,n,p,q,r,s,cargs,flows,precision,step,numstyle,
 fl, flines,carg,initlist,Aer,parm:
options `Copyright 1991 by Daniel Schwalbe`;
precision := max(Digits,trunc(evalhf(Digits)));
npts:=2; m:=16; n:=12; nsteps:=29:

if not type(F,{procedure,numeric}) then
   ERROR(`1st argument (functions) should be a procedure`) fi;
if not type([H,V],[constant..constant,constant..constant]) then
   ERROR(`2nd and 3rd arguments (domains) must be real ranges`) fi;

step:= proc(x) if x<0 then 0 else 1 fi end;
parm:=op(1,eval(F));

f:=traperror(subs('_body'=F(parm),Heaviside=step,
       '_arg'=parm, proc(_arg) _body end)):
if f=lasterror then
  f:=subs('_body'=g(parm),g=eval(F),Heaviside=step,
       '_arg'=parm, proc(_arg) _body end) fi;

p:=evalf(lhs(H),precision); q:=evalf(rhs(H),precision);
r:=evalf(lhs(V),precision); s:=evalf(rhs(V),precision);
h:=evalf((q-p)/nsteps,precision):
initlist:={};
cargs:={args[4..nargs]};
for carg in cargs do
  if type(carg,'set') then
    if type(carg,'set(INITCOND1)') then initlist := carg
   else ERROR(`initial conditions must be a set of points of type`,INITCOND1)
    fi;
  elif type(carg,'list') then 
    if type(carg,INITCOND1) then initlist:={carg}
	else ERROR(`initial condition must be of type`,INITCOND1)
    fi;
  elif not type(carg,equation) then
	ERROR(`optional arguments must be equations`)
  elif op(1,carg)='stepsize' then
     h:=evalf(op(2,carg),precision);
	 nsteps:=trunc((q-p)/h + 1);
     if not type(h,numeric) then
	ERROR(`stepsize must be a numerical value`)
     fi;
  elif op(1,carg)='numsteps' then
     nsteps:=op(2,carg);
     if not type(nsteps,posint) then
	ERROR(`numsteps must be a positive integer`)
     fi;
  elif op(1,carg)='iterations' then
     npts:=op(2,carg);
     if not type(npts,posint) then
	ERROR(`number of iterations must be a positive integer`)
     fi;
  elif op(1,carg)='scheme' then numstyle := 'adaptive'
  elif op(1,carg)='grid' then
     if not type(rhs(carg),[integer,integer]) then
	ERROR(`grid must be a list of two integers`) fi;
     m:=op(2,carg)[1];
     n:=op(2,carg)[2]
  else ERROR(`unknown option`,lhs(carg))
  fi
od;

CC:=array(1..8);
CC[1]:=m: CC[2]:=n: 
CC[3]:=p: CC[4]:=q: CC[5]:=r: CC[6]:=s:
if m>0 and n>0 then
  Ax:=array(0..m,0..2*n+1);
  Ay:=array(0..m,0..2*n+1);
  fl:=traperror(evalhf(ccomputelineshf(
            f,var(Ax),var(Ay),var(CC)))):
  if lasterror=fl then #print(`ccomputelineshf failed`);
    ccomputelines(f,Ax,Ay,p,q,r,s,m,n) fi;
  flines:=seq(seq(CURVE(BLACK,LINE,
      [Ax[i,2*j],Ay[i,2*j],Ax[i,2*j+1],
           Ay[i,2*j+1]]),i=0..m),j=0..n):
else flines:= NULL fi:
i:=1:
for init in initlist do
if init[1]<p or init[1]>q or init[2]<r or init[2]>s then
     print(`initial point `,init,` is not in range`);
  else
  CC[1]:=evalf(init[1],precision);
  CC[2]:=evalf(init[2],precision);
 if numstyle<>'adaptive' then
  A:=array(1..2,0..nsteps):
  Aer:=traperror(evalhf(fieldrkhf(var(A),f,var(CC),h,npts,nsteps)));
  if Aer=lasterror then  
    Aer:= traperror(fieldrkhf(A,f,CC,h,npts,nsteps)) fi;
  if Aer<>lasterror then
    flow[i]:=seq(op([A[1,j],A[2,j]]),j=0..round(CC[7])):
    i:=i+1 
  else print(`solution not found for`,init, Aer);	
	fi;
  A:=array(1..2,0..nsteps):
  Aer:=traperror(evalhf(fieldrkhf(var(A),f,var(CC),-h,npts,nsteps)));
  if Aer=lasterror then  
    Aer:= traperror(fieldrkhf(A,f,CC,-h,npts,nsteps)) fi;
  if Aer<>lasterror then
    flow[i]:=seq(op([A[1,j],A[2,j]]),j=0..round(CC[7])):
    i:=i+1 
  else print(`solution not found for`,init, Aer);	
	fi;
 else
    A:=dadaptiverk(f,CC,h):
    flow[i]:=seq(op(A[j]),j=round(CC[7])..round(CC[8]));
    i:=i+1 fi 
fi;
od:

flows:=seq(CURVE(RED,LINE,[flow[j]]),j=1..i-1);
m:=max(m,16);n:=max(n,12);
PLOT(AXIS(HORIZONTAL,SCALE(4,LINEAR,NUMBER)),
 AXIS(VERTICAL,SCALE(4,LINEAR,NUMBER)),
 RANGE(HORIZONTAL,p-(q-p)/m..q+(q-p)/m),
 RANGE(VERTICAL,r-(s-r)/n..s+(s-r)/n),
 flines,flows):
end:

macro( INITCOND1 = INITCOND1 );
macro(_body=_body);

##################################################

computevectshf2:=proc(Ax,Ay,f1,f2,p,q,r,s,m,n)
local i,j,k,maxlen,ap,bp,bqp,bsr,t,x,y,
 scale,wx,wy,zx,zy,blocks,dx,dy,X,srqp,qpsr;
options `Copyright 1991 by Daniel Schwalbe`;
dx:=(q-p)/m: dy:=(s-r)/n:
if m > n then blocks:= m else blocks:= n fi;
maxlen:=0; 
bqp:=blocks/(q-p);bsr:=blocks/(s-r);
X:=array(0..n,0..m,1..2);
x:=p;
for j from 0 to m do
y:=r;
  for i from 0 to n do
    Ax[j,6*i]:=x;
	Ay[j,6*i]:=y;
	X[i,j,1]:=f1(t,x,y);
	X[i,j,2]:=f2(t,x,y);
    maxlen:=max(sqrt((bqp*X[i,j,1])^2+(bsr*X[i,j,2])^2),maxlen);
      y:=y+dy od; x:=x+dx od;
if blocks < 20 then blocks:=20 fi;
bqp:=blocks/(q-p);bsr:=blocks/(s-r);
srqp:=(q-p)/(s-r);
qpsr:=(s-r)/(q-p);
for j from 0 to m do
  for i from 0 to n do
    k:=6*i;
	x:=X[i,j,1]/maxlen;y:=X[i,j,2]/maxlen;
    if x=0 and y=0 then
	  wx:=0;wy:=0;zx:=0;zy:=0;
	else
	  scale:=sqrt((bqp*x)^2+(bsr*y)^2);
	  wx:=x/scale;
	  wy:=y/scale;
	  zx:=-y*srqp/scale;
	  zy:=x*qpsr/scale;
	fi;

    ap:=Ax[j,k]+2*x/3+wx/3; bp:=Ay[j,k]+2*y/3+wy/3;
    Ax[j,k+1]:=ap;
    Ay[j,k+1]:=bp;
    Ax[j,k+2]:=ap-wx/6+zx/6;
    Ay[j,k+2]:=bp-wy/6+zy/6;
    Ax[j,k+3]:=ap;
    Ay[j,k+3]:=bp;
    Ax[j,k+4]:=ap-wx/6-zx/6;
    Ay[j,k+4]:=bp-wy/6-zy/6;
    Ax[j,k+5]:=ap;
    Ay[j,k+5]:=bp;
    od;od
end:

#main difference in computevects2 from computevectshf2 is trapping
#for errors when evaluating input procedures.

computevects2:=proc(Ax,Ay,f1,f2,p,q,r,s,m,n)
local i,j,k,maxlen,len,a,b,x,y,ap,bp,
 scale,wx,wy,zx,zy,blocks,dx,dy;
options `Copyright 1991 by Daniel Schwalbe`;
dx:=(q-p)/m: dy:=(s-r)/n:
if m > n then blocks:= m else blocks:= n fi;
maxlen:=0;
a:=p;
for j from 0 to m do
  b:=r;
  for i from 0 to n do
    x:= traperror(f1(0,a,b));
	if x = lasterror then x:=0; y:=0 fi;
    y:= traperror(f2(0,a,b));
	if y = lasterror then x:=0; y:=0 fi;
	if not type([x,y],[numeric,numeric]) then
	  x:=0; y:=0 fi;
    len:=((blocks/(q-p)*x)^2+(blocks/(s-r)*y)^2)^(1/2);
    if len > maxlen then maxlen:=len fi;
      b:=b+dy od; a:=a+dx od;
a:=p;
if blocks < 20 then blocks:=20 fi;
for j from 0 to m do
  b:=r;
  for i from 0 to n do
    k:=6*i;
    x:= traperror(f1(0,a,b));
	if x = lasterror then x:=0; y:=0 else
    y:= traperror(f2(0,a,b));
	if y = lasterror then x:=0; y:=0 fi;
    x:=x/maxlen;y:=y/maxlen;
	fi;
	if not type([x,y],[numeric,numeric]) then
	  x:=0; y:=0 fi;
    if x=0 and y=0 then
      scale:= 0 else
      scale:=(((blocks*x/(q-p))^2+(blocks*y/(s-r))^2))^(-1/2) fi;
    wx:= scale*x;   wy:=scale*y;
 #zx and zy are components of perpendicular vector
    zx:= -scale*y/(s-r)*(q-p); zy:= scale*x/(q-p)*(s-r) ;
    Ax[j,k]:=a;
    Ay[j,k]:=b;
    ap:=a+2*x/3+wx/3; bp:=b+2*y/3+wy/3;
    Ax[j,k+1]:=ap;
    Ay[j,k+1]:=bp;
    Ax[j,k+2]:=ap-wx/6+zx/6;
    Ay[j,k+2]:=bp-wy/6+zy/6;
    Ax[j,k+3]:=ap;
    Ay[j,k+3]:=bp;
    Ax[j,k+4]:=ap-wx/6-zx/6;
    Ay[j,k+4]:=bp-wy/6-zy/6;
    Ax[j,k+5]:=ap;
    Ay[j,k+5]:=bp;
    b:=b+dy od; a:=a+dx od
end:

phaserk2:=proc(A,f1,f2,CC,h0,nsteps,npts)
local pt,px,py,p,q,r,s,j,h,
   k11,k21,k31,k41,k12,k22,k32,k42,yt,yx,yy;
options `Copyright 1991 by Daniel Schwalbe`;
h:=h0/npts:
pt:=CC[1]: px:=A[1,0]: py:=A[2,0]:
p:=CC[2]:q:=CC[3]:r:=CC[4]:s:=CC[5]:
j:=0:
while j<nsteps and px<=q and px>=p and
  py<=s and py>=r do
    j:=j+1:
	
  to npts  do
    k11:=f1(pt,px,py):
    k12:=f2(pt,px,py):
	yt:=pt+h/2;
	yx:=px+h*k11/2;
	yy:=py+h*k12/2;
    k21:=f1(yt,yx,yy):
    k22:=f2(yt,yx,yy):
	yx:=px+h*k21/2;
	yy:=py+h*k22/2;
    k31:=f1(yt,yx,yy):
    k32:=f2(yt,yx,yy):
	yt:=pt+h;
	yx:=px+h*k31;
	yy:=py+h*k32;
    k41:=f1(yt,yx,yy):
    k42:=f2(yt,yx,yy):
    pt:=yt;
    px:= px+h*(k11+2*k21+2*k31+k41)/6; 
    py:= py+h*(k12+2*k22+2*k32+k42)/6; 
  od:
	
    A[1,j]:=px: A[2,j]:=py: od:
CC[6]:=j:

j:=0:h:=-h;
pt:=CC[1]: px:=A[1,0]: py:=A[2,0]:
while j>-nsteps and px<=q and px>=p and
  py<=s and py>=r do
    j:=j-1:
	
  to npts  do
    k11:=f1(pt,px,py):
    k12:=f2(pt,px,py):
	yt:=pt+h/2;
	yx:=px+h*k11/2;
	yy:=py+h*k12/2;
    k21:=f1(yt,yx,yy):
    k22:=f2(yt,yx,yy):
	yx:=px+h*k21/2;
	yy:=py+h*k22/2;
    k31:=f1(yt,yx,yy):
    k32:=f2(yt,yx,yy):
	yt:=pt+h;
	yx:=px+h*k31;
	yy:=py+h*k32;
    k41:=f1(yt,yx,yy):
    k42:=f2(yt,yx,yy):
    pt:=yt;
    px:= px+h*(k11+2*k21+2*k31+k41)/6; 
    py:= py+h*(k12+2*k22+2*k32+k42)/6; 
  od:
	
    A[1,j]:=px: A[2,j]:=py: od:
CC[7]:=j:

end:



macro( INITCOND = {[constant,constant],[constant,constant,constant]} );
macro(_body=`plots/ODE/body`);

phaseplot:=proc(F,H,V)
local i,j,p,q,r,s,phasevects,initlist,precision,
  pv,init,A,flow,f1,f2,f,npts,nsteps,m,n,step,
  h,cargs,carg,Ax,Ay,CC,flows,parm:
options `Copyright 1991 by Daniel Schwalbe`;
step:=proc(x) if x<0 then 0 else 1 fi end;


if type(F,list) then 
    parm:=op(1,eval(F[1]));
    f:=traperror([subs(Heaviside=step,F[1](parm)),
	subs(Heaviside=step,F[2](parm))]);
	if f<>lasterror then 
      f1:=subs('_arg'=parm,'_body'=f[1],proc(_arg) _body end);
      f2:=subs('_arg'=parm,'_body'=f[2],proc(_arg) _body end);
	else
	  f1:=subs('_arg'=parm,'_body'=g(parm),
	         g=eval(F[1]),Heaviside=step,
			 proc(_arg) _body end);
	  f2:=subs('_arg'=parm,'_body'=g(parm),
	         g=eval(F[2]),Heaviside=step,
			 proc(_arg) _body end); fi;
else
    parm:=op(1,eval(F)); 
    f:=[subs(Heaviside=step,F(parm)[1]),
	    subs(Heaviside=step,F(parm)[2])];
    f1:=subs('_arg'=parm,'_body'=f[1],proc(_arg) _body end);
    f2:=subs('_arg'=parm,'_body'=f[2],proc(_arg) _body end);
fi;

if not type([H,V],[constant..constant,constant..constant]) then
   ERROR(`2nd and 3rd arguments (domains) must be real ranges`) fi;
step:=proc(z) if z<0 then 0 else 1 fi end:
nsteps:=50:npts:=1;h:=0.1:m:=16:n:=14:
precision := max(Digits,trunc(evalhf(Digits)));

p:=evalf(lhs(H),precision); q:=evalf(rhs(H),precision);
r:=evalf(lhs(V),precision); s:=evalf(rhs(V),precision);
cargs:={args[4..nargs]};
initlist:={};
for carg in cargs do
  if type(carg,'set') then
    if type(carg,'set(INITCOND)') then initlist := carg
    else ERROR(`initial conditions must be a set of points of type`,INITCOND)
    fi;
  elif type(carg,'list') then 
    if type(carg,INITCOND) then initlist:={carg}
	else ERROR(`initial condition must be of type`,INITCOND)
    fi;
  elif not type(carg,equation) then
	ERROR(`optional arguments must be equations`)
  elif op(1,carg)='stepsize' then
     h:=evalf(op(2,carg),precision);
     if not type(h,numeric) then
	ERROR(`stepsize must be a numerical value`)
     fi;
  elif op(1,carg)='numsteps' then
     nsteps:=op(2,carg);
     if not type(nsteps,posint) then
	ERROR(`numsteps must be a positive integer`)
     fi;
  elif op(1,carg)='iterations' then
     npts:=op(2,carg);
     if not type(npts,posint) then
	ERROR(`number of iterations must be a positive integer`)
     fi;
  elif op(1,carg)='grid' then
     if not type(rhs(carg),[integer,integer]) then
	ERROR(`grid must be a list of two positive integers`) fi;
     m:=op(2,carg)[1];
     n:=op(2,carg)[2]
  else ERROR(`unknown option`,lhs(carg))
  fi
od;

if m>0 and n>0 then
  Ax:=array(0..m,0..6*n+5);
  Ay:=array(0..m,0..6*n+5);
  pv:=traperror(evalhf(computevectshf2(var(Ax),var(Ay),
     f1,f2,p,q,r,s,m,n)));
  if lasterror=pv then #print(`computevectshf2 failed`,pv);
    computevects2(Ax,Ay,f1,f2,p,q,r,s,m,n); fi;
  phasevects:=seq(seq(CURVE(BLACK,LINE,
      [Ax[j,6*i],Ay[j,6*i],Ax[j,6*i+1],Ay[j,6*i+1],
      Ax[j,6*i+2],Ay[j,6*i+2],Ax[j,6*i+3],Ay[j,6*i+3],
      Ax[j,6*i+4],Ay[j,6*i+4],Ax[j,6*i+5],Ay[j,6*i+5]]),
        j=0..m),i=0..n):
  else phasevects:= NULL fi;

i:=1;
for init in initlist do
  A:=array(1..2,-nsteps..nsteps);
  CC:=array(1..7):
  if nops(init)=3 then
    CC[1]:=evalf(init[1]);
    A[1,0]:=evalf(init[2]);
    A[2,0]:=evalf(init[3]);
  else
    CC[1]:=0;
    A[1,0]:=evalf(init[1]);
    A[2,0]:=evalf(init[2]) fi;
  CC[2]:=p: CC[3]:=q: CC[4]:=r: CC[5]:=s:
  pv:=traperror(evalhf(
     phaserk2(var(A),f1,f2,var(CC),h,nsteps,npts)));
  if pv = lasterror then 
    pv:=traperror(phaserk2(A,f1,f2,CC,h,nsteps,npts)) fi;
  if pv = lasterror then 
    print(`solution not found for `,init,pv) 
  else	
    flow[i]:= CURVE(RED,LINE,
	   [seq(op([A[1,j],A[2,j]]),j=round(CC[7])..round(CC[6]))]);
    i:=i+1; 
  fi ;
  
od;

flows:=seq(flow[j],j=1..i-1):
m:=max(m,16);n:=max(n,12);

PLOT(AXIS(HORIZONTAL,SCALE(4,LINEAR,NUMBER)),
 AXIS(VERTICAL,SCALE(4,LINEAR,NUMBER)),
 RANGE(HORIZONTAL,p-(q-p)/m..q+(q-p)/m),
 RANGE(VERTICAL,r-(s-r)/n..s+(s-r)/n),
 phasevects,flows):
end:

macro( INITCOND = INITCOND );
macro(_body=_body);


##################################################

picard:=proc(f,init,n)
 local parray,ptemp,t0,y0:
options `Copyright 1991 by Daniel Schwalbe`;
 t0:=init[1]; y0:=init[2];
 ptemp:=y0+int(f('s',y0),'s'=t0..'t');
 parray:=ptemp;
 from 2 to n do
 ptemp:=y0+int(f('s',subs('t'='s',ptemp)),'s'=t0..'t');
 parray:=parray,sort(ptemp) od;
 parray;
end:

##################################################



#save `ODE.m`;
#quit
