[2696] | 1 | function [T, Y, RCNTRL, ICNTRL, RSTAT, ISTAT] = ... |
---|
| 2 | Rosenbrock(Function, Tspan, Y0, Options, RCNTRL, ICNTRL) |
---|
| 3 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 4 | % Rosenbrock - Implementation of several Rosenbrock methods: |
---|
| 5 | % * Ros2 |
---|
| 6 | % * Ros3 |
---|
| 7 | % * Ros4 |
---|
| 8 | % * Rodas3 |
---|
| 9 | % * Rodas4 |
---|
| 10 | % |
---|
| 11 | % Solves the system y'=F(t,y) using a Rosenbrock method defined by: |
---|
| 12 | % |
---|
| 13 | % G = 1/(H*gamma(1)) - Jac(t0,Y0) |
---|
| 14 | % T_i = t0 + Alpha(i)*H |
---|
| 15 | % Y_i = Y0 + \sum_{j=1}^{i-1} A(i,j)*K_j |
---|
| 16 | % G * K_i = Fun( T_i, Y_i ) + \sum_{j=1}^S C(i,j)/H * K_j + |
---|
| 17 | % gamma(i)*dF/dT(t0, Y0) |
---|
| 18 | % Y1 = Y0 + \sum_{j=1}^S M(j)*K_j |
---|
| 19 | % |
---|
| 20 | % For details on Rosenbrock methods and their implementation consult: |
---|
| 21 | % E. Hairer and G. Wanner |
---|
| 22 | % "Solving ODEs II. Stiff and differential-algebraic problems". |
---|
| 23 | % Springer series in computational mathematics, Springer-Verlag, 1996. |
---|
| 24 | % The codes contained in the book inspired this implementation. |
---|
| 25 | % |
---|
| 26 | % MATLAB implementation (C) John C. Linford (jlinford@vt.edu). |
---|
| 27 | % Virginia Polytechnic Institute and State University |
---|
| 28 | % November, 2009 |
---|
| 29 | % |
---|
| 30 | % Based on the Fortran90 implementation (C) Adrian Sandu, August 2004 |
---|
| 31 | % and revised by Philipp Miehe and Adrian Sandu, May 2006. |
---|
| 32 | % Virginia Polytechnic Institute and State University |
---|
| 33 | % Contact: sandu@cs.vt.edu |
---|
| 34 | % |
---|
| 35 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 36 | % Input Arguments : |
---|
| 37 | % The first four arguments are similar to the input arguments of |
---|
| 38 | % MATLAB's ODE solvers |
---|
| 39 | % Function - A function handle for the ODE function |
---|
| 40 | % Tspan - The time space to integrate |
---|
| 41 | % Y0 - Initial value |
---|
| 42 | % Options - ODE solver options created by odeset(): |
---|
| 43 | % AbsTol, InitialStep, Jacobian, MaxStep, and RelTol |
---|
| 44 | % are considered. Other options are ignored. |
---|
| 45 | % 'Jacobian' must be set. |
---|
| 46 | % RCNTRL - real value input parameters (explained below) |
---|
| 47 | % ICNTRL - integer input parameters (explained below) |
---|
| 48 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 49 | % Output Arguments: |
---|
| 50 | % The first two arguments are similar to the output arguments of |
---|
| 51 | % MATLAB's ODE solvers |
---|
| 52 | % T - A vector of final integration times. |
---|
| 53 | % Y - A matrix of function values. Y(T(i),:) is the value of |
---|
| 54 | % the function at the ith output time. |
---|
| 55 | % RCNTRL - real value input parameters (explained below) |
---|
| 56 | % ICNTRL - integer input parameters (explained below) |
---|
| 57 | % RSTAT - real output parameters (explained below) |
---|
| 58 | % ISTAT - integer output parameters (explained below) |
---|
| 59 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 60 | % |
---|
| 61 | % RCNTRL and ICNTRL on input and output: |
---|
| 62 | % |
---|
| 63 | % Note: For input parameters equal to zero the default values of the |
---|
| 64 | % corresponding variables are used. |
---|
| 65 | % |
---|
| 66 | % ICNTRL(1) = 1: F = F(y) Independent of T (AUTONOMOUS) |
---|
| 67 | % = 0: F = F(t,y) Depends on T (NON-AUTONOMOUS) |
---|
| 68 | % |
---|
| 69 | % ICNTRL(2) = 0: AbsTol, RelTol are N-dimensional vectors |
---|
| 70 | % = 1: AbsTol, RelTol are scalars |
---|
| 71 | % |
---|
| 72 | % ICNTRL(3) -> selection of a particular Rosenbrock method |
---|
| 73 | % = 0 : Rodas3 (default) |
---|
| 74 | % = 1 : Ros2 |
---|
| 75 | % = 2 : Ros3 |
---|
| 76 | % = 3 : Ros4 |
---|
| 77 | % = 4 : Rodas3 |
---|
| 78 | % = 5 : Rodas4 |
---|
| 79 | % |
---|
| 80 | % ICNTRL(4) -> maximum number of integration steps |
---|
| 81 | % For ICNTRL(4)=0) the default value of 100000 is used |
---|
| 82 | % |
---|
| 83 | % RCNTRL(1) -> Hmin, lower bound for the integration step size |
---|
| 84 | % It is strongly recommended to keep Hmin = ZERO |
---|
| 85 | % RCNTRL(2) -> Hmax, upper bound for the integration step size |
---|
| 86 | % RCNTRL(3) -> Hstart, starting value for the integration step size |
---|
| 87 | % |
---|
| 88 | % RCNTRL(4) -> FacMin, lower bound on step decrease factor (default=0.2) |
---|
| 89 | % RCNTRL(5) -> FacMax, upper bound on step increase factor (default=6) |
---|
| 90 | % RCNTRL(6) -> FacRej, step decrease factor after multiple rejections |
---|
| 91 | % (default=0.1) |
---|
| 92 | % RCNTRL(7) -> FacSafe, by which the new step is slightly smaller |
---|
| 93 | % than the predicted value (default=0.9) |
---|
| 94 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 95 | % |
---|
| 96 | % RSTAT and ISTAT on output: |
---|
| 97 | % |
---|
| 98 | % ISTAT(1) -> No. of function calls |
---|
| 99 | % ISTAT(2) -> No. of jacobian calls |
---|
| 100 | % ISTAT(3) -> No. of steps |
---|
| 101 | % ISTAT(4) -> No. of accepted steps |
---|
| 102 | % ISTAT(5) -> No. of rejected steps (except at very beginning) |
---|
| 103 | % ISTAT(6) -> No. of LU decompositions |
---|
| 104 | % ISTAT(7) -> No. of forward/backward substitutions |
---|
| 105 | % ISTAT(8) -> No. of singular matrix decompositions |
---|
| 106 | % |
---|
| 107 | % RSTAT(1) -> Texit, the time corresponding to the |
---|
| 108 | % computed Y upon return |
---|
| 109 | % RSTAT(2) -> Hexit, last accepted step before exit |
---|
| 110 | % RSTAT(3) -> Hnew, last predicted step (not yet taken) |
---|
| 111 | % For multiple restarts, use Hnew as Hstart |
---|
| 112 | % in the subsequent run |
---|
| 113 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 114 | |
---|
| 115 | %~~~> Statistics on the work performed by the Rosenbrock method |
---|
| 116 | global Nfun Njac Nstp Nacc Nrej Ndec Nsol Nsng Ntexit Nhexit Nhnew |
---|
| 117 | global ISTATUS RSTATUS |
---|
| 118 | |
---|
| 119 | % Parse ODE options |
---|
| 120 | AbsTol = odeget(Options, 'AbsTol'); |
---|
| 121 | if isempty(AbsTol) |
---|
| 122 | AbsTol = 1.0e-3; |
---|
| 123 | end |
---|
| 124 | Hstart = odeget(Options, 'InitialStep'); |
---|
| 125 | if isempty(Hstart) |
---|
| 126 | Hstart = 0; |
---|
| 127 | end |
---|
| 128 | Jacobian = odeget(Options, 'Jacobian'); |
---|
| 129 | if isempty(Jacobian) |
---|
| 130 | error('A Jacobian function is required.'); |
---|
| 131 | end |
---|
| 132 | Hmax = odeget(Options, 'MaxStep'); |
---|
| 133 | if isempty(Hmax) |
---|
| 134 | Hmax = 0; |
---|
| 135 | end |
---|
| 136 | RelTol = odeget(Options, 'RelTol'); |
---|
| 137 | if isempty(RelTol) |
---|
| 138 | RelTol = 1.0e-4; |
---|
| 139 | end |
---|
| 140 | |
---|
| 141 | % Initialize statistics |
---|
| 142 | Nfun=1; Njac=2; Nstp=3; Nacc=4; Nrej=5; Ndec=6; Nsol=7; Nsng=8; |
---|
| 143 | Ntexit=1; Nhexit=2; Nhnew=3; |
---|
| 144 | RSTATUS = zeros(20,1); |
---|
| 145 | ISTATUS = zeros(20,1); |
---|
| 146 | |
---|
| 147 | % Get problem size |
---|
| 148 | steps = length(Tspan); |
---|
| 149 | N = length(Y0); |
---|
| 150 | |
---|
| 151 | % Initialize tolerances |
---|
| 152 | ATOL = ones(N,1)*AbsTol; |
---|
| 153 | RTOL = ones(N,1)*RelTol; |
---|
| 154 | |
---|
| 155 | % Initialize step |
---|
| 156 | RCNTRL(2) = max(0, Hmax); |
---|
| 157 | RCNTRL(3) = max(0, Hstart); |
---|
| 158 | |
---|
| 159 | % Integrate |
---|
| 160 | Y = zeros(steps,N); |
---|
| 161 | T = zeros(steps,1); |
---|
| 162 | Y(1,:) = Y0; |
---|
| 163 | T(1) = Tspan(1); |
---|
| 164 | for i=2:steps |
---|
| 165 | [T(i), Y(i,:), IERR] = ... |
---|
| 166 | RosenbrockIntegrator(N, Y(i-1,:), T(i-1), Tspan(i), ... |
---|
| 167 | Function, Jacobian, ATOL, RTOL, RCNTRL, ICNTRL); |
---|
| 168 | |
---|
| 169 | if IERR < 0 |
---|
| 170 | error(['Rosenbrock exited with IERR=',num2str(IERR)]); |
---|
| 171 | end |
---|
| 172 | |
---|
| 173 | end |
---|
| 174 | |
---|
| 175 | ISTAT = ISTATUS; |
---|
| 176 | RSTAT = RSTATUS; |
---|
| 177 | |
---|
| 178 | return |
---|
| 179 | |
---|
| 180 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 181 | function [ T, Y, IERR ] = RosenbrockIntegrator(N, Y, Tstart, Tend, ... |
---|
| 182 | OdeFunction, OdeJacobian, AbsTol, RelTol, RCNTRL, ICNTRL) |
---|
| 183 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 184 | % Advances the ODE system defined by OdeFunction, OdeJacobian and Y |
---|
| 185 | % from time=Tstart to time=Tend |
---|
| 186 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 187 | |
---|
| 188 | %~~~> Statistics on the work performed by the Rosenbrock method |
---|
| 189 | global Nfun Njac Nstp Nacc Nrej Ndec Nsol Nsng Ntexit Nhexit Nhnew Ntotal |
---|
| 190 | global ISTATUS RSTATUS |
---|
| 191 | |
---|
| 192 | %~~~> Parameters |
---|
| 193 | global ZERO ONE DeltaMin |
---|
| 194 | ZERO = 0.0; ONE = 1.0; DeltaMin = 1.0E-5; |
---|
| 195 | |
---|
| 196 | %~~~> Initialize statistics |
---|
| 197 | ISTATUS(1:8) = 0; |
---|
| 198 | RSTATUS(1:3) = ZERO; |
---|
| 199 | |
---|
| 200 | %~~~> Autonomous or time dependent ODE. Default is time dependent. |
---|
| 201 | Autonomous = ~(ICNTRL(1) == 0); |
---|
| 202 | |
---|
| 203 | %~~~> For Scalar tolerances (ICNTRL(2).NE.0) the code uses AbsTol(1) and RelTol(1) |
---|
| 204 | % For Vector tolerances (ICNTRL(2) == 0) the code uses AbsTol(1:N) and RelTol(1:N) |
---|
| 205 | if ICNTRL(2) == 0 |
---|
| 206 | VectorTol = true; |
---|
| 207 | UplimTol = N; |
---|
| 208 | else |
---|
| 209 | VectorTol = false; |
---|
| 210 | UplimTol = 1; |
---|
| 211 | end |
---|
| 212 | |
---|
| 213 | %~~~> Initialize the particular Rosenbrock method selected |
---|
| 214 | switch (ICNTRL(3)) |
---|
| 215 | case (1) |
---|
| 216 | ros_Param = Ros2; |
---|
| 217 | case (2) |
---|
| 218 | ros_Param = Ros3; |
---|
| 219 | case (3) |
---|
| 220 | ros_Param = Ros4; |
---|
| 221 | case {0,4} |
---|
| 222 | ros_Param = Rodas3; |
---|
| 223 | case (5) |
---|
| 224 | ros_Param = Rodas4; |
---|
| 225 | otherwise |
---|
| 226 | disp(['Unknown Rosenbrock method: ICNTRL(3)=',num2str(ICNTRL(3))]); |
---|
| 227 | IERR = ros_ErrorMsg(-2,Tstart,ZERO); |
---|
| 228 | return |
---|
| 229 | end |
---|
| 230 | ros_S = ros_Param{1}; |
---|
| 231 | rosMethod = ros_Param{2}; |
---|
| 232 | ros_A = ros_Param{3}; |
---|
| 233 | ros_C = ros_Param{4}; |
---|
| 234 | ros_M = ros_Param{5}; |
---|
| 235 | ros_E = ros_Param{6}; |
---|
| 236 | ros_Alpha = ros_Param{7}; |
---|
| 237 | ros_Gamma = ros_Param{8}; |
---|
| 238 | ros_ELO = ros_Param{9}; |
---|
| 239 | ros_NewF = ros_Param{10}; |
---|
| 240 | ros_Name = ros_Param{11}; |
---|
| 241 | |
---|
| 242 | %~~~> The maximum number of steps admitted |
---|
| 243 | if ICNTRL(4) == 0 |
---|
| 244 | Max_no_steps = 200000; |
---|
| 245 | elseif ICNTRL(4) > 0 |
---|
| 246 | Max_no_steps=ICNTRL(4); |
---|
| 247 | else |
---|
| 248 | disp(['User-selected max no. of steps: ICNTRL(4)=',num2str(ICNTRL(4))]); |
---|
| 249 | IERR = ros_ErrorMsg(-1,Tstart,ZERO); |
---|
| 250 | return |
---|
| 251 | end |
---|
| 252 | |
---|
| 253 | %~~~> Unit roundoff (1+Roundoff>1) |
---|
| 254 | Roundoff = eps; |
---|
| 255 | |
---|
| 256 | %~~~> Lower bound on the step size: (positive value) |
---|
| 257 | if RCNTRL(1) == ZERO |
---|
| 258 | Hmin = ZERO; |
---|
| 259 | elseif RCNTRL(1) > ZERO |
---|
| 260 | Hmin = RCNTRL(1); |
---|
| 261 | else |
---|
| 262 | disp(['User-selected Hmin: RCNTRL(1)=', num2str(RCNTRL(1))]); |
---|
| 263 | IERR = ros_ErrorMsg(-3,Tstart,ZERO); |
---|
| 264 | return |
---|
| 265 | end |
---|
| 266 | |
---|
| 267 | %~~~> Upper bound on the step size: (positive value) |
---|
| 268 | if RCNTRL(2) == ZERO |
---|
| 269 | Hmax = abs(Tend-Tstart); |
---|
| 270 | elseif RCNTRL(2) > ZERO |
---|
| 271 | Hmax = min(abs(RCNTRL(2)),abs(Tend-Tstart)); |
---|
| 272 | else |
---|
| 273 | disp(['User-selected Hmax: RCNTRL(2)=', num2str(RCNTRL(2))]); |
---|
| 274 | IERR = ros_ErrorMsg(-3,Tstart,ZERO); |
---|
| 275 | return |
---|
| 276 | end |
---|
| 277 | |
---|
| 278 | %~~~> Starting step size: (positive value) |
---|
| 279 | if RCNTRL(3) == ZERO |
---|
| 280 | Hstart = max(Hmin,DeltaMin); |
---|
| 281 | elseif RCNTRL(3) > ZERO |
---|
| 282 | Hstart = min(abs(RCNTRL(3)),abs(Tend-Tstart)); |
---|
| 283 | else |
---|
| 284 | disp(['User-selected Hstart: RCNTRL(3)=', num2str(RCNTRL(3))]); |
---|
| 285 | IERR = ros_ErrorMsg(-3,Tstart,ZERO); |
---|
| 286 | return |
---|
| 287 | end |
---|
| 288 | |
---|
| 289 | %~~~> Step size can be changed s.t. FacMin < Hnew/Hold < FacMax |
---|
| 290 | if RCNTRL(4) == ZERO |
---|
| 291 | FacMin = 0.2; |
---|
| 292 | elseif RCNTRL(4) > ZERO |
---|
| 293 | FacMin = RCNTRL(4); |
---|
| 294 | else |
---|
| 295 | disp(['User-selected FacMin: RCNTRL(4)=', num2str(RCNTRL(4))]); |
---|
| 296 | IERR = ros_ErrorMsg(-4,Tstart,ZERO); |
---|
| 297 | return |
---|
| 298 | end |
---|
| 299 | if RCNTRL(5) == ZERO |
---|
| 300 | FacMax = 6.0; |
---|
| 301 | elseif RCNTRL(5) > ZERO |
---|
| 302 | FacMax = RCNTRL(5); |
---|
| 303 | else |
---|
| 304 | disp(['User-selected FacMax: RCNTRL(5)=', num2str(RCNTRL(5))]); |
---|
| 305 | IERR = ros_ErrorMsg(-4,Tstart,ZERO); |
---|
| 306 | return |
---|
| 307 | end |
---|
| 308 | |
---|
| 309 | %~~~> FacRej: Factor to decrease step after 2 succesive rejections |
---|
| 310 | if RCNTRL(6) == ZERO |
---|
| 311 | FacRej = 0.1; |
---|
| 312 | elseif RCNTRL(6) > ZERO |
---|
| 313 | FacRej = RCNTRL(6); |
---|
| 314 | else |
---|
| 315 | disp(['User-selected FacRej: RCNTRL(6)=', num2str(RCNTRL(6))]); |
---|
| 316 | IERR = ros_ErrorMsg(-4,Tstart,ZERO); |
---|
| 317 | return |
---|
| 318 | end |
---|
| 319 | |
---|
| 320 | %~~~> FacSafe: Safety Factor in the computation of new step size |
---|
| 321 | if RCNTRL(7) == ZERO |
---|
| 322 | FacSafe = 0.9; |
---|
| 323 | elseif RCNTRL(7) > ZERO |
---|
| 324 | FacSafe = RCNTRL(7); |
---|
| 325 | else |
---|
| 326 | disp(['User-selected FacSafe: RCNTRL(7)=', num2str(RCNTRL(7))]); |
---|
| 327 | IERR = ros_ErrorMsg(-4,Tstart,ZERO); |
---|
| 328 | return |
---|
| 329 | end |
---|
| 330 | |
---|
| 331 | %~~~> Check if tolerances are reasonable |
---|
| 332 | for i=1:UplimTol |
---|
| 333 | if (AbsTol(i) <= ZERO) || (RelTol(i) <= 10.0*Roundoff) || (RelTol(i) >= 1.0) |
---|
| 334 | disp([' AbsTol(',i,') = ',num2str(AbsTol(i))]); |
---|
| 335 | disp([' RelTol(',i,') = ',num2str(RelTol(i))]); |
---|
| 336 | IERR = ros_ErrorMsg(-5,Tstart,ZERO); |
---|
| 337 | return |
---|
| 338 | end |
---|
| 339 | end |
---|
| 340 | |
---|
| 341 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 342 | % Template for the implementation of a generic Rosenbrock method |
---|
| 343 | % defined by ros_S (no of stages) |
---|
| 344 | % and its coefficients ros_{A,C,M,E,Alpha,Gamma} |
---|
| 345 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 346 | |
---|
| 347 | % ~~~~ Local variables |
---|
| 348 | Ynew = zeros(N,1); |
---|
| 349 | Fcn0 = zeros(N,1); |
---|
| 350 | Fcn = zeros(N,1); |
---|
| 351 | K = zeros(N,ros_S); |
---|
| 352 | dFdT = zeros(N,1); |
---|
| 353 | Jac0 = sparse(N,N); |
---|
| 354 | Ghimj = sparse(N,N); |
---|
| 355 | Yerr = zeros(N,1); |
---|
| 356 | Pivot = zeros(N,1); |
---|
| 357 | |
---|
| 358 | %~~~> Initial preparations |
---|
| 359 | T = Tstart; |
---|
| 360 | RSTATUS(Nhexit) = ZERO; |
---|
| 361 | H = min(max(abs(Hmin), abs(Hstart)), abs(Hmax)); |
---|
| 362 | if abs(H) <= 10.0*Roundoff |
---|
| 363 | H = DeltaMin; |
---|
| 364 | end |
---|
| 365 | |
---|
| 366 | if Tend >= Tstart |
---|
| 367 | Direction = +1; |
---|
| 368 | else |
---|
| 369 | Direction = -1; |
---|
| 370 | end |
---|
| 371 | H = Direction*H; |
---|
| 372 | |
---|
| 373 | RejectLastH = false; |
---|
| 374 | RejectMoreH = false; |
---|
| 375 | |
---|
| 376 | %~~~> Time loop begins below |
---|
| 377 | |
---|
| 378 | while ( (Direction > 0) && ((T-Tend)+Roundoff <= ZERO) || ... |
---|
| 379 | (Direction < 0) && ((Tend-T)+Roundoff <= ZERO) ) |
---|
| 380 | |
---|
| 381 | if ISTATUS(Nstp) > Max_no_steps % Too many steps |
---|
| 382 | IERR = ros_ErrorMsg(-6,T,H); |
---|
| 383 | return |
---|
| 384 | end |
---|
| 385 | if ((T+0.1*H) == T) || (H <= Roundoff) % Step size too small |
---|
| 386 | IERR = ros_ErrorMsg(-7,T,H); |
---|
| 387 | return |
---|
| 388 | end |
---|
| 389 | |
---|
| 390 | %~~~> Limit H if necessary to avoid going beyond Tend |
---|
| 391 | H = min(H, abs(Tend-T)); |
---|
| 392 | |
---|
| 393 | %~~~> Compute the function at current time |
---|
| 394 | Fcn0 = OdeFunction(T, Y); |
---|
| 395 | ISTATUS(Nfun) = ISTATUS(Nfun) + 1; |
---|
| 396 | |
---|
| 397 | %~~~> Compute the function derivative with respect to T |
---|
| 398 | if ~Autonomous |
---|
| 399 | dFdT = ros_FunTimeDerivative (T, OdeFunction, Roundoff, Y, Fcn0); |
---|
| 400 | end |
---|
| 401 | |
---|
| 402 | %~~~> Compute the Jacobian at current time |
---|
| 403 | Jac0 = OdeJacobian(T,Y); |
---|
| 404 | ISTATUS(Njac) = ISTATUS(Njac) + 1; |
---|
| 405 | |
---|
| 406 | %~~~> Repeat step calculation until current step accepted |
---|
| 407 | while(true) |
---|
| 408 | |
---|
| 409 | [H, Ghimj, Pivot, Singular] = ... |
---|
| 410 | ros_PrepareMatrix(N, H, Direction, ros_Gamma(1), Jac0); |
---|
| 411 | |
---|
| 412 | % Not calculating LU decomposition in ros_PrepareMatrix anymore |
---|
| 413 | % so don't need to check if the matrix is singular |
---|
| 414 | % |
---|
| 415 | %if Singular % More than 5 consecutive failed decompositions |
---|
| 416 | % IERR = ros_ErrorMsg(-8,T,H); |
---|
| 417 | % return |
---|
| 418 | %end |
---|
| 419 | |
---|
| 420 | % For the 1st istage the function has been computed previously |
---|
| 421 | Fcn = Fcn0; |
---|
| 422 | K(:,1) = Fcn; |
---|
| 423 | if (~Autonomous) && (ros_Gamma(1) ~= ZERO) |
---|
| 424 | HG = Direction*H*ros_Gamma(1); |
---|
| 425 | K(:,1) = K(:,1) + HG * dFdT; |
---|
| 426 | end |
---|
| 427 | K(:,1) = ros_Solve(Ghimj, Pivot, K(:,1)); |
---|
| 428 | |
---|
| 429 | %~~~> Compute the remaining stages |
---|
| 430 | for istage=2:ros_S |
---|
| 431 | |
---|
| 432 | % istage>1 and a new function evaluation is needed |
---|
| 433 | if ros_NewF(istage) |
---|
| 434 | Ynew = Y; |
---|
| 435 | for j=1:istage-1 |
---|
| 436 | Ynew = Ynew + ros_A((istage-1)*(istage-2)/2+j)*K(:,j)'; |
---|
| 437 | end |
---|
| 438 | Tau = T + ros_Alpha(istage)*Direction*H; |
---|
| 439 | Fcn = OdeFunction(Tau,Ynew); |
---|
| 440 | ISTATUS(Nfun) = ISTATUS(Nfun) + 1; |
---|
| 441 | end |
---|
| 442 | |
---|
| 443 | K(:,istage) = Fcn; |
---|
| 444 | for j=1:istage-1 |
---|
| 445 | HC = ros_C((istage-1)*(istage-2)/2+j)/(Direction*H); |
---|
| 446 | K(:,istage) = K(:,istage) + HC * K(:,j); |
---|
| 447 | end |
---|
| 448 | if (~Autonomous) && (ros_Gamma(istage) ~= ZERO) |
---|
| 449 | HG = Direction*H*ros_Gamma(istage); |
---|
| 450 | K(:,istage) = K(:,istage) + HG * dFdT; |
---|
| 451 | end |
---|
| 452 | |
---|
| 453 | % Linear system is solved here with MATLAB '\' operator |
---|
| 454 | % instead of LU decomposition. |
---|
| 455 | K(:,istage) = ros_Solve(Ghimj, Pivot, K(:,istage)); |
---|
| 456 | |
---|
| 457 | end |
---|
| 458 | |
---|
| 459 | %~~~> Compute the new solution |
---|
| 460 | Ynew = Y; |
---|
| 461 | for j=1:ros_S |
---|
| 462 | Ynew = Ynew + ros_M(j) * K(:,j)'; |
---|
| 463 | end |
---|
| 464 | |
---|
| 465 | %~~~> Compute the error estimation |
---|
| 466 | Yerr = zeros(N,1); |
---|
| 467 | for j=1:ros_S |
---|
| 468 | Yerr = Yerr + ros_E(j) * K(:,j); |
---|
| 469 | end |
---|
| 470 | Err = ros_ErrorNorm(N, Y, Ynew, Yerr, AbsTol, RelTol, VectorTol); |
---|
| 471 | |
---|
| 472 | %~~~> New step size is bounded by FacMin <= Hnew/H <= FacMax |
---|
| 473 | Fac = min(FacMax, max(FacMin, FacSafe/(Err^(ONE/ros_ELO)))); |
---|
| 474 | Hnew = H*Fac; |
---|
| 475 | |
---|
| 476 | %~~~> Check the error magnitude and adjust step size |
---|
| 477 | ISTATUS(Nstp) = ISTATUS(Nstp) + 1; |
---|
| 478 | if (Err <= ONE) || (H <= Hmin) %~~~> Accept step |
---|
| 479 | ISTATUS(Nacc) = ISTATUS(Nacc) + 1; |
---|
| 480 | Y = Ynew; |
---|
| 481 | T = T + Direction*H; |
---|
| 482 | Hnew = max(Hmin, min(Hnew, Hmax)); |
---|
| 483 | if RejectLastH % No step size increase after a rejected step |
---|
| 484 | Hnew = min(Hnew, H); |
---|
| 485 | end |
---|
| 486 | RSTATUS(Nhexit) = H; |
---|
| 487 | RSTATUS(Nhnew) = Hnew; |
---|
| 488 | RSTATUS(Ntexit) = T; |
---|
| 489 | RejectLastH = false; |
---|
| 490 | RejectMoreH = false; |
---|
| 491 | H = Hnew; |
---|
| 492 | break % EXIT THE LOOP: WHILE STEP NOT ACCEPTED |
---|
| 493 | else %~~~> Reject step |
---|
| 494 | if RejectMoreH |
---|
| 495 | Hnew = H*FacRej; |
---|
| 496 | end |
---|
| 497 | RejectMoreH = RejectLastH; |
---|
| 498 | RejectLastH = true; |
---|
| 499 | H = Hnew; |
---|
| 500 | if ISTATUS(Nacc) >= 1 |
---|
| 501 | ISTATUS(Nrej) = ISTATUS(Nrej) + 1; |
---|
| 502 | end |
---|
| 503 | end % Err <= 1 |
---|
| 504 | |
---|
| 505 | end % UntilAccepted |
---|
| 506 | |
---|
| 507 | end % TimeLoop |
---|
| 508 | |
---|
| 509 | %~~~> Succesful exit |
---|
| 510 | IERR = 1; %~~~> The integration was successful |
---|
| 511 | |
---|
| 512 | return |
---|
| 513 | |
---|
| 514 | |
---|
| 515 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 516 | function [ ErrNorm ] = ros_ErrorNorm(N, Y, Ynew, Yerr, AbsTol, RelTol, VectorTol) |
---|
| 517 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 518 | %~~~> Computes the "scaled norm" of the error vector Yerr |
---|
| 519 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 520 | |
---|
| 521 | Err = 0; |
---|
| 522 | for i=1:N |
---|
| 523 | Ymax = max(abs(Y(i)), abs(Ynew(i))); |
---|
| 524 | if VectorTol |
---|
| 525 | Scale = AbsTol(i) + RelTol(i)*Ymax; |
---|
| 526 | else |
---|
| 527 | Scale = AbsTol(1) + RelTol(1)*Ymax; |
---|
| 528 | end |
---|
| 529 | Err = Err + (Yerr(i)/Scale)^2; |
---|
| 530 | end |
---|
| 531 | Err = sqrt(Err/N); |
---|
| 532 | ErrNorm = max(Err,1.0e-10); |
---|
| 533 | |
---|
| 534 | return |
---|
| 535 | |
---|
| 536 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 537 | function [ dFdT ] = ros_FunTimeDerivative (T, OdeFunction, Roundoff, Y, Fcn0) |
---|
| 538 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 539 | %~~~> The time partial derivative of the function by finite differences |
---|
| 540 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 541 | |
---|
| 542 | %~~~> Statistics on the work performed by the Rosenbrock method |
---|
| 543 | global Nfun ISTATUS |
---|
| 544 | |
---|
| 545 | %~~~> Parameters |
---|
| 546 | global ZERO ONE DeltaMin |
---|
| 547 | |
---|
| 548 | Delta = sqrt(Roundoff) * max(DeltaMin, abs(T)); |
---|
| 549 | dFdT = OdeFunction(T+Delta,Y); |
---|
| 550 | ISTATUS(Nfun) = ISTATUS(Nfun) + 1; |
---|
| 551 | dFdT = (dFdT - Fcn0) / Delta; |
---|
| 552 | |
---|
| 553 | return |
---|
| 554 | |
---|
| 555 | |
---|
| 556 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 557 | function [H, Ghimj, Pivot, Singular] = ... |
---|
| 558 | ros_PrepareMatrix(N, H, Direction, gam, Jac0) |
---|
| 559 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 560 | % Prepares the LHS matrix for stage calculations |
---|
| 561 | % 1. Construct Ghimj = 1/(H*ham) - Jac0 |
---|
| 562 | % "(Gamma H) Inverse Minus Jacobian" |
---|
| 563 | % 2. LU decomposition not performed here because |
---|
| 564 | % MATLAB solves the linear system with '\'. |
---|
| 565 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 566 | |
---|
| 567 | Ghimj = -Jac0; |
---|
| 568 | ghinv = 1.0/(Direction*H*gam); |
---|
| 569 | for i=1:N |
---|
| 570 | Ghimj(i,i) = Ghimj(i,i)+ghinv; |
---|
| 571 | end |
---|
| 572 | |
---|
| 573 | Pivot(1) = 1; |
---|
| 574 | Singular = false; |
---|
| 575 | |
---|
| 576 | return |
---|
| 577 | |
---|
| 578 | |
---|
| 579 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 580 | function [ IERR ] = ros_ErrorMsg(Code, T, H) |
---|
| 581 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 582 | % Handles all error messages |
---|
| 583 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 584 | |
---|
| 585 | IERR = Code; |
---|
| 586 | |
---|
| 587 | disp('Forced exit from Rosenbrock due to the following error:'); |
---|
| 588 | |
---|
| 589 | switch (Code) |
---|
| 590 | case (-1) |
---|
| 591 | disp('--> Improper value for maximal no of steps'); |
---|
| 592 | case (-2) |
---|
| 593 | disp('--> Selected Rosenbrock method not implemented'); |
---|
| 594 | case (-3) |
---|
| 595 | disp('--> Hmin/Hmax/Hstart must be positive'); |
---|
| 596 | case (-4) |
---|
| 597 | disp('--> FacMin/FacMax/FacRej must be positive'); |
---|
| 598 | case (-5) |
---|
| 599 | disp('--> Improper tolerance values'); |
---|
| 600 | case (-6) |
---|
| 601 | disp('--> No of steps exceeds maximum bound'); |
---|
| 602 | case (-7) |
---|
| 603 | disp('--> Step size too small: T + 10*H = T or H < Roundoff'); |
---|
| 604 | case (-8) |
---|
| 605 | disp('--> Matrix is repeatedly singular'); |
---|
| 606 | otherwise |
---|
| 607 | disp(['Unknown Error code: ', num2str(Code)]); |
---|
| 608 | end |
---|
| 609 | |
---|
| 610 | disp(['T=',num2str(T),' and H=',num2str(H)]); |
---|
| 611 | |
---|
| 612 | return |
---|
| 613 | |
---|
| 614 | |
---|
| 615 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 616 | function [X] = ros_Solve(JVS, Pivot, X) |
---|
| 617 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 618 | % Template for the forward/backward substitution |
---|
| 619 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 620 | |
---|
| 621 | global Nsol ISTATUS |
---|
| 622 | |
---|
| 623 | X = JVS\X; |
---|
| 624 | Pivot(1) = 1; |
---|
| 625 | ISTATUS(Nsol) = ISTATUS(Nsol) + 1; |
---|
| 626 | |
---|
| 627 | return |
---|
| 628 | |
---|
| 629 | |
---|
| 630 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 631 | function [ params ] = Ros2() |
---|
| 632 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 633 | % --- AN L-STABLE METHOD, 2 stages, order 2 |
---|
| 634 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 635 | |
---|
| 636 | g = 1.0 + 1.0/sqrt(2.0); |
---|
| 637 | rosMethod = 1; |
---|
| 638 | %~~~> Name of the method |
---|
| 639 | ros_Name = 'ROS-2'; |
---|
| 640 | %~~~> Number of stages |
---|
| 641 | ros_S = 2; |
---|
| 642 | |
---|
| 643 | %~~~> The coefficient matrices A and C are strictly lower triangular. |
---|
| 644 | % The lower triangular (subdiagonal) elements are stored in row-wise order: |
---|
| 645 | % A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc. |
---|
| 646 | % The general mapping formula is: |
---|
| 647 | % A(i,j) = ros_A( (i-1)*(i-2)/2 + j ) |
---|
| 648 | % C(i,j) = ros_C( (i-1)*(i-2)/2 + j ) |
---|
| 649 | |
---|
| 650 | ros_A(1) = (1.0)/g; |
---|
| 651 | ros_C(1) = (-2.0)/g; |
---|
| 652 | %~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE) |
---|
| 653 | % or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE) |
---|
| 654 | ros_NewF(1) = true; |
---|
| 655 | ros_NewF(2) = true; |
---|
| 656 | %~~~> M_i = Coefficients for new step solution |
---|
| 657 | ros_M(1)= (3.0)/(2.0*g); |
---|
| 658 | ros_M(2)= (1.0)/(2.0*g); |
---|
| 659 | % E_i = Coefficients for error estimator |
---|
| 660 | ros_E(1) = 1.0/(2.0*g); |
---|
| 661 | ros_E(2) = 1.0/(2.0*g); |
---|
| 662 | %~~~> ros_ELO = estimator of local order - the minimum between the |
---|
| 663 | % main and the embedded scheme orders plus one |
---|
| 664 | ros_ELO = 2.0; |
---|
| 665 | %~~~> Y_stage_i ~ Y( T + H*Alpha_i ) |
---|
| 666 | ros_Alpha(1) = 0.0; |
---|
| 667 | ros_Alpha(2) = 1.0; |
---|
| 668 | %~~~> Gamma_i = \sum_j gamma_{i,j} |
---|
| 669 | ros_Gamma(1) = g; |
---|
| 670 | ros_Gamma(2) =-g; |
---|
| 671 | |
---|
| 672 | params = { ros_S, rosMethod, ros_A, ros_C, ros_M, ros_E, ... |
---|
| 673 | ros_Alpha, ros_Gamma, ros_ELO, ros_NewF, ros_Name }; |
---|
| 674 | |
---|
| 675 | return |
---|
| 676 | |
---|
| 677 | |
---|
| 678 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 679 | function [ params ] = Ros3() |
---|
| 680 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 681 | % --- AN L-STABLE METHOD, 3 stages, order 3, 2 function evaluations |
---|
| 682 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 683 | |
---|
| 684 | rosMethod = 2; |
---|
| 685 | %~~~> Name of the method |
---|
| 686 | ros_Name = 'ROS-3'; |
---|
| 687 | %~~~> Number of stages |
---|
| 688 | ros_S = 3; |
---|
| 689 | |
---|
| 690 | %~~~> The coefficient matrices A and C are strictly lower triangular. |
---|
| 691 | % The lower triangular (subdiagonal) elements are stored in row-wise order: |
---|
| 692 | % A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc. |
---|
| 693 | % The general mapping formula is: |
---|
| 694 | % A(i,j) = ros_A( (i-1)*(i-2)/2 + j ) |
---|
| 695 | % C(i,j) = ros_C( (i-1)*(i-2)/2 + j ) |
---|
| 696 | |
---|
| 697 | ros_A(1)= 1.0; |
---|
| 698 | ros_A(2)= 1.0; |
---|
| 699 | ros_A(3)= 0.0; |
---|
| 700 | |
---|
| 701 | ros_C(1) = -0.10156171083877702091975600115545E+01; |
---|
| 702 | ros_C(2) = 0.40759956452537699824805835358067E+01; |
---|
| 703 | ros_C(3) = 0.92076794298330791242156818474003E+01; |
---|
| 704 | %~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE) |
---|
| 705 | % or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE) |
---|
| 706 | ros_NewF(1) = true; |
---|
| 707 | ros_NewF(2) = true; |
---|
| 708 | ros_NewF(3) = false; |
---|
| 709 | %~~~> M_i = Coefficients for new step solution |
---|
| 710 | ros_M(1) = 0.1E+01; |
---|
| 711 | ros_M(2) = 0.61697947043828245592553615689730E+01; |
---|
| 712 | ros_M(3) = -0.42772256543218573326238373806514; |
---|
| 713 | % E_i = Coefficients for error estimator |
---|
| 714 | ros_E(1) = 0.5; |
---|
| 715 | ros_E(2) = -0.29079558716805469821718236208017E+01; |
---|
| 716 | ros_E(3) = 0.22354069897811569627360909276199; |
---|
| 717 | %~~~> ros_ELO = estimator of local order - the minimum between the |
---|
| 718 | % main and the embedded scheme orders plus 1 |
---|
| 719 | ros_ELO = 3.0; |
---|
| 720 | %~~~> Y_stage_i ~ Y( T + H*Alpha_i ) |
---|
| 721 | ros_Alpha(1)= 0.0; |
---|
| 722 | ros_Alpha(2)= 0.43586652150845899941601945119356; |
---|
| 723 | ros_Alpha(3)= 0.43586652150845899941601945119356; |
---|
| 724 | %~~~> Gamma_i = \sum_j gamma_{i,j} |
---|
| 725 | ros_Gamma(1)= 0.43586652150845899941601945119356; |
---|
| 726 | ros_Gamma(2)= 0.24291996454816804366592249683314; |
---|
| 727 | ros_Gamma(3)= 0.21851380027664058511513169485832E+01; |
---|
| 728 | |
---|
| 729 | params = { ros_S, rosMethod, ros_A, ros_C, ros_M, ros_E, ... |
---|
| 730 | ros_Alpha, ros_Gamma, ros_ELO, ros_NewF, ros_Name }; |
---|
| 731 | |
---|
| 732 | return |
---|
| 733 | |
---|
| 734 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 735 | |
---|
| 736 | |
---|
| 737 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 738 | function [ params ] = Ros4() |
---|
| 739 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 740 | % L-STABLE ROSENBROCK METHOD OF ORDER 4, WITH 4 STAGES |
---|
| 741 | % L-STABLE EMBEDDED ROSENBROCK METHOD OF ORDER 3 |
---|
| 742 | % |
---|
| 743 | % E. HAIRER AND G. WANNER, SOLVING ORDINARY DIFFERENTIAL |
---|
| 744 | % EQUATIONS II. STIFF AND DIFFERENTIAL-ALGEBRAIC PROBLEMS. |
---|
| 745 | % SPRINGER SERIES IN COMPUTATIONAL MATHEMATICS, |
---|
| 746 | % SPRINGER-VERLAG (1990) |
---|
| 747 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 748 | |
---|
| 749 | rosMethod = 3; |
---|
| 750 | %~~~> Name of the method |
---|
| 751 | ros_Name = 'ROS-4'; |
---|
| 752 | %~~~> Number of stages |
---|
| 753 | ros_S = 4; |
---|
| 754 | |
---|
| 755 | %~~~> The coefficient matrices A and C are strictly lower triangular. |
---|
| 756 | % The lower triangular (subdiagonal) elements are stored in row-wise order: |
---|
| 757 | % A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc. |
---|
| 758 | % The general mapping formula is: |
---|
| 759 | % A(i,j) = ros_A( (i-1)*(i-2)/2 + j ) |
---|
| 760 | % C(i,j) = ros_C( (i-1)*(i-2)/2 + j ) |
---|
| 761 | |
---|
| 762 | ros_A(1) = 0.2000000000000000E+01; |
---|
| 763 | ros_A(2) = 0.1867943637803922E+01; |
---|
| 764 | ros_A(3) = 0.2344449711399156; |
---|
| 765 | ros_A(4) = ros_A(2); |
---|
| 766 | ros_A(5) = ros_A(3); |
---|
| 767 | ros_A(6) = 0.0; |
---|
| 768 | |
---|
| 769 | ros_C(1) =-0.7137615036412310E+01; |
---|
| 770 | ros_C(2) = 0.2580708087951457E+01; |
---|
| 771 | ros_C(3) = 0.6515950076447975; |
---|
| 772 | ros_C(4) =-0.2137148994382534E+01; |
---|
| 773 | ros_C(5) =-0.3214669691237626; |
---|
| 774 | ros_C(6) =-0.6949742501781779; |
---|
| 775 | %~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE) |
---|
| 776 | % or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE) |
---|
| 777 | ros_NewF(1) = true; |
---|
| 778 | ros_NewF(2) = true; |
---|
| 779 | ros_NewF(3) = true; |
---|
| 780 | ros_NewF(4) = false; |
---|
| 781 | %~~~> M_i = Coefficients for new step solution |
---|
| 782 | ros_M(1) = 0.2255570073418735E+01; |
---|
| 783 | ros_M(2) = 0.2870493262186792; |
---|
| 784 | ros_M(3) = 0.4353179431840180; |
---|
| 785 | ros_M(4) = 0.1093502252409163E+01; |
---|
| 786 | %~~~> E_i = Coefficients for error estimator |
---|
| 787 | ros_E(1) =-0.2815431932141155; |
---|
| 788 | ros_E(2) =-0.7276199124938920E-01; |
---|
| 789 | ros_E(3) =-0.1082196201495311; |
---|
| 790 | ros_E(4) =-0.1093502252409163E+01; |
---|
| 791 | %~~~> ros_ELO = estimator of local order - the minimum between the |
---|
| 792 | % main and the embedded scheme orders plus 1 |
---|
| 793 | ros_ELO = 4.0; |
---|
| 794 | %~~~> Y_stage_i ~ Y( T + H*Alpha_i ) |
---|
| 795 | ros_Alpha(1) = 0.0; |
---|
| 796 | ros_Alpha(2) = 0.1145640000000000E+01; |
---|
| 797 | ros_Alpha(3) = 0.6552168638155900; |
---|
| 798 | ros_Alpha(4) = ros_Alpha(3); |
---|
| 799 | %~~~> Gamma_i = \sum_j gamma_{i,j} |
---|
| 800 | ros_Gamma(1) = 0.5728200000000000; |
---|
| 801 | ros_Gamma(2) =-0.1769193891319233E+01; |
---|
| 802 | ros_Gamma(3) = 0.7592633437920482; |
---|
| 803 | ros_Gamma(4) =-0.1049021087100450; |
---|
| 804 | |
---|
| 805 | params = { ros_S, rosMethod, ros_A, ros_C, ros_M, ros_E, ... |
---|
| 806 | ros_Alpha, ros_Gamma, ros_ELO, ros_NewF, ros_Name }; |
---|
| 807 | |
---|
| 808 | return |
---|
| 809 | |
---|
| 810 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 811 | function [ params ] = Rodas3() |
---|
| 812 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 813 | % --- A STIFFLY-STABLE METHOD, 4 stages, order 3 |
---|
| 814 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 815 | |
---|
| 816 | rosMethod = 4; |
---|
| 817 | %~~~> Name of the method |
---|
| 818 | ros_Name = 'RODAS-3'; |
---|
| 819 | %~~~> Number of stages |
---|
| 820 | ros_S = 4; |
---|
| 821 | |
---|
| 822 | %~~~> The coefficient matrices A and C are strictly lower triangular. |
---|
| 823 | % The lower triangular (subdiagonal) elements are stored in row-wise order: |
---|
| 824 | % A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc. |
---|
| 825 | % The general mapping formula is: |
---|
| 826 | % A(i,j) = ros_A( (i-1)*(i-2)/2 + j ) |
---|
| 827 | % C(i,j) = ros_C( (i-1)*(i-2)/2 + j ) |
---|
| 828 | |
---|
| 829 | ros_A(1) = 0.0; |
---|
| 830 | ros_A(2) = 2.0; |
---|
| 831 | ros_A(3) = 0.0; |
---|
| 832 | ros_A(4) = 2.0; |
---|
| 833 | ros_A(5) = 0.0; |
---|
| 834 | ros_A(6) = 1.0; |
---|
| 835 | |
---|
| 836 | ros_C(1) = 4.0; |
---|
| 837 | ros_C(2) = 1.0; |
---|
| 838 | ros_C(3) =-1.0; |
---|
| 839 | ros_C(4) = 1.0; |
---|
| 840 | ros_C(5) =-1.0; |
---|
| 841 | ros_C(6) =-(8.0/3.0); |
---|
| 842 | |
---|
| 843 | %~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE) |
---|
| 844 | % or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE) |
---|
| 845 | ros_NewF(1) = true; |
---|
| 846 | ros_NewF(2) = false; |
---|
| 847 | ros_NewF(3) = true; |
---|
| 848 | ros_NewF(4) = true; |
---|
| 849 | %~~~> M_i = Coefficients for new step solution |
---|
| 850 | ros_M(1) = 2.0; |
---|
| 851 | ros_M(2) = 0.0; |
---|
| 852 | ros_M(3) = 1.0; |
---|
| 853 | ros_M(4) = 1.0; |
---|
| 854 | %~~~> E_i = Coefficients for error estimator |
---|
| 855 | ros_E(1) = 0.0; |
---|
| 856 | ros_E(2) = 0.0; |
---|
| 857 | ros_E(3) = 0.0; |
---|
| 858 | ros_E(4) = 1.0; |
---|
| 859 | %~~~> ros_ELO = estimator of local order - the minimum between the |
---|
| 860 | % main and the embedded scheme orders plus 1 |
---|
| 861 | ros_ELO = 3.0; |
---|
| 862 | %~~~> Y_stage_i ~ Y( T + H*Alpha_i ) |
---|
| 863 | ros_Alpha(1) = 0.0; |
---|
| 864 | ros_Alpha(2) = 0.0; |
---|
| 865 | ros_Alpha(3) = 1.0; |
---|
| 866 | ros_Alpha(4) = 1.0; |
---|
| 867 | %~~~> Gamma_i = \sum_j gamma_{i,j} |
---|
| 868 | ros_Gamma(1) = 0.5; |
---|
| 869 | ros_Gamma(2) = 1.5; |
---|
| 870 | ros_Gamma(3) = 0.0; |
---|
| 871 | ros_Gamma(4) = 0.0; |
---|
| 872 | |
---|
| 873 | params = { ros_S, rosMethod, ros_A, ros_C, ros_M, ros_E, ... |
---|
| 874 | ros_Alpha, ros_Gamma, ros_ELO, ros_NewF, ros_Name }; |
---|
| 875 | |
---|
| 876 | return |
---|
| 877 | |
---|
| 878 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 879 | function [ params ] = Rodas4() |
---|
| 880 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 881 | % STIFFLY-STABLE ROSENBROCK METHOD OF ORDER 4, WITH 6 STAGES |
---|
| 882 | % |
---|
| 883 | % E. HAIRER AND G. WANNER, SOLVING ORDINARY DIFFERENTIAL |
---|
| 884 | % EQUATIONS II. STIFF AND DIFFERENTIAL-ALGEBRAIC PROBLEMS. |
---|
| 885 | % SPRINGER SERIES IN COMPUTATIONAL MATHEMATICS, |
---|
| 886 | % SPRINGER-VERLAG (1996) |
---|
| 887 | %~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 888 | |
---|
| 889 | rosMethod = 5; |
---|
| 890 | %~~~> Name of the method |
---|
| 891 | ros_Name = 'RODAS-4'; |
---|
| 892 | %~~~> Number of stages |
---|
| 893 | ros_S = 6; |
---|
| 894 | |
---|
| 895 | %~~~> Y_stage_i ~ Y( T + H*Alpha_i ) |
---|
| 896 | ros_Alpha(1) = 0.000; |
---|
| 897 | ros_Alpha(2) = 0.386; |
---|
| 898 | ros_Alpha(3) = 0.210; |
---|
| 899 | ros_Alpha(4) = 0.630; |
---|
| 900 | ros_Alpha(5) = 1.000; |
---|
| 901 | ros_Alpha(6) = 1.000; |
---|
| 902 | |
---|
| 903 | %~~~> Gamma_i = \sum_j gamma_{i,j} |
---|
| 904 | ros_Gamma(1) = 0.2500000000000000; |
---|
| 905 | ros_Gamma(2) =-0.1043000000000000; |
---|
| 906 | ros_Gamma(3) = 0.1035000000000000; |
---|
| 907 | ros_Gamma(4) =-0.3620000000000023E-01; |
---|
| 908 | ros_Gamma(5) = 0.0; |
---|
| 909 | ros_Gamma(6) = 0.0; |
---|
| 910 | |
---|
| 911 | %~~~> The coefficient matrices A and C are strictly lower triangular. |
---|
| 912 | % The lower triangular (subdiagonal) elements are stored in row-wise order: |
---|
| 913 | % A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc. |
---|
| 914 | % The general mapping formula is: A(i,j) = ros_A( (i-1)*(i-2)/2 + j ) |
---|
| 915 | % C(i,j) = ros_C( (i-1)*(i-2)/2 + j ) |
---|
| 916 | |
---|
| 917 | ros_A(1) = 0.1544000000000000E+01; |
---|
| 918 | ros_A(2) = 0.9466785280815826; |
---|
| 919 | ros_A(3) = 0.2557011698983284; |
---|
| 920 | ros_A(4) = 0.3314825187068521E+01; |
---|
| 921 | ros_A(5) = 0.2896124015972201E+01; |
---|
| 922 | ros_A(6) = 0.9986419139977817; |
---|
| 923 | ros_A(7) = 0.1221224509226641E+01; |
---|
| 924 | ros_A(8) = 0.6019134481288629E+01; |
---|
| 925 | ros_A(9) = 0.1253708332932087E+02; |
---|
| 926 | ros_A(10) =-0.6878860361058950; |
---|
| 927 | ros_A(11) = ros_A(7); |
---|
| 928 | ros_A(12) = ros_A(8); |
---|
| 929 | ros_A(13) = ros_A(9); |
---|
| 930 | ros_A(14) = ros_A(10); |
---|
| 931 | ros_A(15) = 1.0; |
---|
| 932 | |
---|
| 933 | ros_C(1) =-0.5668800000000000E+01; |
---|
| 934 | ros_C(2) =-0.2430093356833875E+01; |
---|
| 935 | ros_C(3) =-0.2063599157091915; |
---|
| 936 | ros_C(4) =-0.1073529058151375; |
---|
| 937 | ros_C(5) =-0.9594562251023355E+01; |
---|
| 938 | ros_C(6) =-0.2047028614809616E+02; |
---|
| 939 | ros_C(7) = 0.7496443313967647E+01; |
---|
| 940 | ros_C(8) =-0.1024680431464352E+02; |
---|
| 941 | ros_C(9) =-0.3399990352819905E+02; |
---|
| 942 | ros_C(10) = 0.1170890893206160E+02; |
---|
| 943 | ros_C(11) = 0.8083246795921522E+01; |
---|
| 944 | ros_C(12) =-0.7981132988064893E+01; |
---|
| 945 | ros_C(13) =-0.3152159432874371E+02; |
---|
| 946 | ros_C(14) = 0.1631930543123136E+02; |
---|
| 947 | ros_C(15) =-0.6058818238834054E+01; |
---|
| 948 | |
---|
| 949 | %~~~> M_i = Coefficients for new step solution |
---|
| 950 | ros_M(1) = ros_A(7); |
---|
| 951 | ros_M(2) = ros_A(8); |
---|
| 952 | ros_M(3) = ros_A(9); |
---|
| 953 | ros_M(4) = ros_A(10); |
---|
| 954 | ros_M(5) = 1.0; |
---|
| 955 | ros_M(6) = 1.0; |
---|
| 956 | |
---|
| 957 | %~~~> E_i = Coefficients for error estimator |
---|
| 958 | ros_E(1) = 0.0; |
---|
| 959 | ros_E(2) = 0.0; |
---|
| 960 | ros_E(3) = 0.0; |
---|
| 961 | ros_E(4) = 0.0; |
---|
| 962 | ros_E(5) = 0.0; |
---|
| 963 | ros_E(6) = 1.0; |
---|
| 964 | |
---|
| 965 | %~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE) |
---|
| 966 | % or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE) |
---|
| 967 | ros_NewF(1) = true; |
---|
| 968 | ros_NewF(2) = true; |
---|
| 969 | ros_NewF(3) = true; |
---|
| 970 | ros_NewF(4) = true; |
---|
| 971 | ros_NewF(5) = true; |
---|
| 972 | ros_NewF(6) = true; |
---|
| 973 | |
---|
| 974 | %~~~> ros_ELO = estimator of local order - the minimum between the |
---|
| 975 | % main and the embedded scheme orders plus 1 |
---|
| 976 | ros_ELO = 4.0; |
---|
| 977 | |
---|
| 978 | params = { ros_S, rosMethod, ros_A, ros_C, ros_M, ros_E, ... |
---|
| 979 | ros_Alpha, ros_Gamma, ros_ELO, ros_NewF, ros_Name }; |
---|
| 980 | |
---|
| 981 | return |
---|
| 982 | |
---|
| 983 | % End of INTEGRATE function |
---|
| 984 | % ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 985 | |
---|
| 986 | |
---|