Mercurial > octave
comparison scripts/ode/private/runge_kutta_45_dorpri.m @ 30893:e1788b1a315f
maint: Use "fcn" as preferred abbreviation for "function" in m-files.
* accumarray.m, accumdim.m, quadl.m, quadv.m, randi.m, structfun.m,
__is_function__.m, uigetfile.m, uimenu.m, uiputfile.m, doc_cache_create.m,
colorspace_conversion_input_check.m, imageIO.m, argnames.m, vectorize.m,
vectorize.m, normest1.m, inputname.m, nthargout.m, display_info_file.m,
decic.m, ode15i.m, ode15s.m, ode23.m, ode23s.m, ode45.m, odeset.m,
check_default_input.m, integrate_adaptive.m, ode_event_handler.m,
runge_kutta_23.m, runge_kutta_23s.m, runge_kutta_45_dorpri.m,
runge_kutta_interpolate.m, starting_stepsize.m, __all_opts__.m, fminbnd.m,
fminsearch.m, fminunc.m, fsolve.m, fzero.m, sqp.m, fplot.m, plotyy.m,
__bar__.m, __ezplot__.m, flat_entry.html, profexport.m, movfun.m, bicg.m,
bicgstab.m, cgs.m, eigs.m, gmres.m, pcg.m, __alltohandles__.m, __sprand__.m,
qmr.m, tfqmr.m, dump_demos.m:
Replace "func", "fun", "fn" in documentation and variable names with "fcn".
author | Rik <rik@octave.org> |
---|---|
date | Mon, 04 Apr 2022 18:14:56 -0700 |
parents | 796f54d4ddbf |
children | 597f3ee61a48 |
comparison
equal
deleted
inserted
replaced
30892:1a3cc2811090 | 30893:e1788b1a315f |
---|---|
22 ## <https://www.gnu.org/licenses/>. | 22 ## <https://www.gnu.org/licenses/>. |
23 ## | 23 ## |
24 ######################################################################## | 24 ######################################################################## |
25 | 25 |
26 ## -*- texinfo -*- | 26 ## -*- texinfo -*- |
27 ## @deftypefn {} {[@var{t_next}, @var{x_next}] =} runge_kutta_45_dorpri (@var{@@fun}, @var{t}, @var{x}, @var{dt}) | 27 ## @deftypefn {} {[@var{t_next}, @var{x_next}] =} runge_kutta_45_dorpri (@var{@@fcn}, @var{t}, @var{x}, @var{dt}) |
28 ## @deftypefnx {} {[@var{t_next}, @var{x_next}] =} runge_kutta_45_dorpri (@var{@@fun}, @var{t}, @var{x}, @var{dt}, @var{options}) | 28 ## @deftypefnx {} {[@var{t_next}, @var{x_next}] =} runge_kutta_45_dorpri (@var{@@fcn}, @var{t}, @var{x}, @var{dt}, @var{options}) |
29 ## @deftypefnx {} {[@var{t_next}, @var{x_next}] =} runge_kutta_45_dorpri (@var{@@fun}, @var{t}, @var{x}, @var{dt}, @var{options}, @var{k_vals}) | 29 ## @deftypefnx {} {[@var{t_next}, @var{x_next}] =} runge_kutta_45_dorpri (@var{@@fcn}, @var{t}, @var{x}, @var{dt}, @var{options}, @var{k_vals}) |
30 ## @deftypefnx {} {[@var{t_next}, @var{x_next}] =} runge_kutta_45_dorpri (@var{@@fun}, @var{t}, @var{x}, @var{dt}, @var{options}, @var{k_vals}, @var{t_next}) | 30 ## @deftypefnx {} {[@var{t_next}, @var{x_next}] =} runge_kutta_45_dorpri (@var{@@fcn}, @var{t}, @var{x}, @var{dt}, @var{options}, @var{k_vals}, @var{t_next}) |
31 ## @deftypefnx {} {[@var{t_next}, @var{x_next}, @var{x_est}] =} runge_kutta_45_dorpri (@dots{}) | 31 ## @deftypefnx {} {[@var{t_next}, @var{x_next}, @var{x_est}] =} runge_kutta_45_dorpri (@dots{}) |
32 ## @deftypefnx {} {[@var{t_next}, @var{x_next}, @var{x_est}, @var{k_vals_out}] =} runge_kutta_45_dorpri (@dots{}) | 32 ## @deftypefnx {} {[@var{t_next}, @var{x_next}, @var{x_est}, @var{k_vals_out}] =} runge_kutta_45_dorpri (@dots{}) |
33 ## | 33 ## |
34 ## This function can be used to integrate a system of ODEs with a given initial | 34 ## This function can be used to integrate a system of ODEs with a given initial |
35 ## condition @var{x} from @var{t} to @var{t+dt} with the Dormand-Prince method. | 35 ## condition @var{x} from @var{t} to @var{t+dt} with the Dormand-Prince method. |
65 ## | 65 ## |
66 ## Fourth output parameter is matrix containing the Runge-Kutta evaluations | 66 ## Fourth output parameter is matrix containing the Runge-Kutta evaluations |
67 ## to use in an FSAL scheme or for dense output. | 67 ## to use in an FSAL scheme or for dense output. |
68 ## @end deftypefn | 68 ## @end deftypefn |
69 | 69 |
70 function [t_next, x_next, x_est, k] = runge_kutta_45_dorpri (fun, t, x, dt, | 70 function [t_next, x_next, x_est, k] = runge_kutta_45_dorpri (fcn, t, x, dt, |
71 options = [], | 71 options = [], |
72 k_vals = [], | 72 k_vals = [], |
73 t_next = t + dt) | 73 t_next = t + dt) |
74 | 74 |
75 ## Reference: Hairer, Ernst; Nørsett, Syvert Paul; Wanner, Gerhard (2008), | 75 ## Reference: Hairer, Ernst; Nørsett, Syvert Paul; Wanner, Gerhard (2008), |
98 endif | 98 endif |
99 | 99 |
100 if (! isempty (k_vals)) # k values from previous step are passed | 100 if (! isempty (k_vals)) # k values from previous step are passed |
101 k(:,1) = k_vals(:,end); # FSAL property | 101 k(:,1) = k_vals(:,end); # FSAL property |
102 else | 102 else |
103 k(:,1) = feval (fun, t, x, args{:}); | 103 k(:,1) = feval (fcn, t, x, args{:}); |
104 endif | 104 endif |
105 | 105 |
106 k(:,2) = feval (fun, s(2), x + k(:,1) * aa(2, 1).' , args{:}); | 106 k(:,2) = feval (fcn, s(2), x + k(:,1) * aa(2, 1).' , args{:}); |
107 k(:,3) = feval (fun, s(3), x + k(:,1:2) * aa(3, 1:2).', args{:}); | 107 k(:,3) = feval (fcn, s(3), x + k(:,1:2) * aa(3, 1:2).', args{:}); |
108 k(:,4) = feval (fun, s(4), x + k(:,1:3) * aa(4, 1:3).', args{:}); | 108 k(:,4) = feval (fcn, s(4), x + k(:,1:3) * aa(4, 1:3).', args{:}); |
109 k(:,5) = feval (fun, s(5), x + k(:,1:4) * aa(5, 1:4).', args{:}); | 109 k(:,5) = feval (fcn, s(5), x + k(:,1:4) * aa(5, 1:4).', args{:}); |
110 k(:,6) = feval (fun, s(6), x + k(:,1:5) * aa(6, 1:5).', args{:}); | 110 k(:,6) = feval (fcn, s(6), x + k(:,1:5) * aa(6, 1:5).', args{:}); |
111 | 111 |
112 ## compute new time and new values for the unknowns | 112 ## compute new time and new values for the unknowns |
113 ## t_next = t + dt; | 113 ## t_next = t + dt; |
114 x_next = x + k(:,1:6) * cc(:); # 5th order approximation | 114 x_next = x + k(:,1:6) * cc(:); # 5th order approximation |
115 | 115 |
116 ## if the estimation of the error is required | 116 ## if the estimation of the error is required |
117 if (nargout >= 3) | 117 if (nargout >= 3) |
118 ## new solution to be compared with the previous one | 118 ## new solution to be compared with the previous one |
119 k(:,7) = feval (fun, t_next, x_next, args{:}); | 119 k(:,7) = feval (fcn, t_next, x_next, args{:}); |
120 cc_prime = dt * c_prime; | 120 cc_prime = dt * c_prime; |
121 x_est = x + k * cc_prime(:); | 121 x_est = x + k * cc_prime(:); |
122 endif | 122 endif |
123 | 123 |
124 endfunction | 124 endfunction |