Original WRF subgrid support version from John Michalakes without fire
[wrffire.git] / wrfv2_fire / chem / KPP / kpp / kpp-2.1 / int / rosenbrock.f90
blob6946f9384ba9bd24a52c11591bb24634125df387
1 MODULE KPP_ROOT_Integrator
3 IMPLICIT NONE
4 PUBLIC
5 SAVE
6 !~~~> Statistics on the work performed by the Rosenbrock method
7 INTEGER :: Nfun,Njac,Nstp,Nacc,Nrej,Ndec,Nsol,Nsng
8 INTEGER, PARAMETER :: ifun=1, ijac=2, istp=3, iacc=4, &
9 irej=5, idec=6, isol=7, isng=8, itexit=1, ihexit=2
12 ! description of the error numbers IERR
13 CHARACTER(LEN=50), PARAMETER, DIMENSION(-8:1) :: IERR_NAMES = (/ &
14 'Matrix is repeatedly singular ', & ! -8
15 'Step size too small ', & ! -7
16 'No of steps exceeds maximum bound ', & ! -6
17 'Improper tolerance values ', & ! -5
18 'FacMin/FacMax/FacRej must be positive ', & ! -4
19 'Hmin/Hmax/Hstart must be positive ', & ! -3
20 'Selected Rosenbrock method not implemented ', & ! -2
21 'Improper value for maximal no of steps ', & ! -1
22 ' ', & ! 0 (not used)
23 'Success ' /) ! 1
25 CONTAINS
27 SUBROUTINE INTEGRATE( TIN, TOUT, &
28 ICNTRL_U, RCNTRL_U, ISTATUS_U, RSTATUS_U, IERR_U )
30 USE KPP_ROOT_Parameters
31 USE KPP_ROOT_Global
32 IMPLICIT NONE
34 KPP_REAL, INTENT(IN) :: TIN ! Start Time
35 KPP_REAL, INTENT(IN) :: TOUT ! End Time
36 ! Optional input parameters and statistics
37 INTEGER, INTENT(IN), OPTIONAL :: ICNTRL_U(20)
38 KPP_REAL, INTENT(IN), OPTIONAL :: RCNTRL_U(20)
39 INTEGER, INTENT(OUT), OPTIONAL :: ISTATUS_U(20)
40 KPP_REAL, INTENT(OUT), OPTIONAL :: RSTATUS_U(20)
41 INTEGER, INTENT(OUT), OPTIONAL :: IERR_U
43 INTEGER :: N_stp, N_acc, N_rej, N_sng
44 SAVE N_stp, N_acc, N_rej, N_sng
45 INTEGER :: i, IERR
46 KPP_REAL :: RCNTRL(20), RSTATUS(20)
47 INTEGER :: ICNTRL(20), ISTATUS(20)
49 ICNTRL(:) = 0
50 RCNTRL(:) = 0.0_dp
51 ISTATUS(:) = 0
52 RSTATUS(:) = 0.0_dp
54 ! If optional parameters are given, and if they are >0,
55 ! then they overwrite default settings.
56 IF (PRESENT(ICNTRL_U)) THEN
57 WHERE(ICNTRL_U(:) > 0) ICNTRL(:) = ICNTRL_U(:)
58 END IF
59 IF (PRESENT(RCNTRL_U)) THEN
60 WHERE(RCNTRL_U(:) > 0) RCNTRL(:) = RCNTRL_U(:)
61 END IF
63 CALL Rosenbrock(VAR,TIN,TOUT, &
64 ATOL,RTOL, &
65 RCNTRL,ICNTRL,RSTATUS,ISTATUS,IERR)
67 STEPMIN = RCNTRL(ihexit)
68 ! if optional parameters are given for output they to return information
69 IF (PRESENT(ISTATUS_U)) ISTATUS_U(:) = ISTATUS(:)
70 IF (PRESENT(RSTATUS_U)) RSTATUS_U(:) = RSTATUS(:)
71 IF (PRESENT(IERR_U)) IERR_U = IERR
73 END SUBROUTINE INTEGRATE
75 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
76 SUBROUTINE Rosenbrock(Y,Tstart,Tend, &
77 AbsTol,RelTol, &
78 RCNTRL,ICNTRL,RSTATUS,ISTATUS,IERR)
79 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
81 ! Solves the system y'=F(t,y) using a Rosenbrock method defined by:
83 ! G = 1/(H*gamma(1)) - Jac(t0,Y0)
84 ! T_i = t0 + Alpha(i)*H
85 ! Y_i = Y0 + \sum_{j=1}^{i-1} A(i,j)*K_j
86 ! G * K_i = Fun( T_i, Y_i ) + \sum_{j=1}^S C(i,j)/H * K_j +
87 ! gamma(i)*dF/dT(t0, Y0)
88 ! Y1 = Y0 + \sum_{j=1}^S M(j)*K_j
90 ! For details on Rosenbrock methods and their implementation consult:
91 ! E. Hairer and G. Wanner
92 ! "Solving ODEs II. Stiff and differential-algebraic problems".
93 ! Springer series in computational mathematics, Springer-Verlag, 1996.
94 ! The codes contained in the book inspired this implementation.
96 ! (C) Adrian Sandu, August 2004
97 ! Virginia Polytechnic Institute and State University
98 ! Contact: sandu@cs.vt.edu
99 ! This implementation is part of KPP - the Kinetic PreProcessor
100 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
102 !~~~> INPUT ARGUMENTS:
104 !- Y(NVAR) = vector of initial conditions (at T=Tstart)
105 !- [Tstart,Tend] = time range of integration
106 ! (if Tstart>Tend the integration is performed backwards in time)
107 !- RelTol, AbsTol = user precribed accuracy
108 !- SUBROUTINE Fun( T, Y, Ydot ) = ODE function,
109 ! returns Ydot = Y' = F(T,Y)
110 !- SUBROUTINE Jac( T, Y, Jcb ) = Jacobian of the ODE function,
111 ! returns Jcb = dFun/dY
112 !- ICNTRL(1:20) = integer inputs parameters
113 !- RCNTRL(1:20) = real inputs parameters
114 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
116 !~~~> OUTPUT ARGUMENTS:
118 !- Y(NVAR) -> vector of final states (at T->Tend)
119 !- ISTATUS(1:20) -> integer output parameters
120 !- RSTATUS(1:20) -> real output parameters
121 !- IERR -> job status upon return
122 ! success (positive value) or
123 ! failure (negative value)
124 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
126 !~~~> INPUT PARAMETERS:
128 ! Note: For input parameters equal to zero the default values of the
129 ! corresponding variables are used.
131 ! ICNTRL(1) = 1: F = F(y) Independent of T (AUTONOMOUS)
132 ! = 0: F = F(t,y) Depends on T (NON-AUTONOMOUS)
134 ! ICNTRL(2) = 0: AbsTol, RelTol are NVAR-dimensional vectors
135 ! = 1: AbsTol, RelTol are scalars
137 ! ICNTRL(3) -> selection of a particular Rosenbrock method
138 ! = 0 : default method is Rodas3
139 ! = 1 : method is Ros2
140 ! = 2 : method is Ros3
141 ! = 3 : method is Ros4
142 ! = 4 : method is Rodas3
143 ! = 5: method is Rodas4
145 ! ICNTRL(4) -> maximum number of integration steps
146 ! For ICNTRL(4)=0) the default value of 100000 is used
148 ! RCNTRL(1) -> Hmin, lower bound for the integration step size
149 ! It is strongly recommended to keep Hmin = ZERO
150 ! RCNTRL(2) -> Hmax, upper bound for the integration step size
151 ! RCNTRL(3) -> Hstart, starting value for the integration step size
153 ! RCNTRL(4) -> FacMin, lower bound on step decrease factor (default=0.2)
154 ! RCNTRL(5) -> FacMax, upper bound on step increase factor (default=6)
155 ! RCNTRL(6) -> FacRej, step decrease factor after multiple rejections
156 ! (default=0.1)
157 ! RCNTRL(7) -> FacSafe, by which the new step is slightly smaller
158 ! than the predicted value (default=0.9)
159 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
161 !~~~> OUTPUT PARAMETERS:
163 ! Note: each call to Rosenbrock adds the current no. of fcn calls
164 ! to previous value of ISTATUS(1), and similar for the other params.
165 ! Set ISTATUS(1:20) = 0 before call to avoid this accumulation.
167 ! ISTATUS(1) = No. of function calls
168 ! ISTATUS(2) = No. of jacobian calls
169 ! ISTATUS(3) = No. of steps
170 ! ISTATUS(4) = No. of accepted steps
171 ! ISTATUS(5) = No. of rejected steps (except at the beginning)
172 ! ISTATUS(6) = No. of LU decompositions
173 ! ISTATUS(7) = No. of forward/backward substitutions
174 ! ISTATUS(8) = No. of singular matrix decompositions
176 ! RSTATUS(1) -> Texit, the time corresponding to the
177 ! computed Y upon return
178 ! RSTATUS(2) -> Hexit, last accepted step before exit
179 ! For multiple restarts, use Hexit as Hstart in the following run
180 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
182 USE KPP_ROOT_Parameters
183 USE KPP_ROOT_LinearAlgebra
184 IMPLICIT NONE
186 !~~~> Arguments
187 KPP_REAL, INTENT(INOUT) :: Y(NVAR)
188 KPP_REAL, INTENT(IN) :: Tstart,Tend
189 KPP_REAL, INTENT(IN) :: AbsTol(NVAR),RelTol(NVAR)
190 INTEGER, INTENT(IN) :: ICNTRL(20)
191 KPP_REAL, INTENT(IN) :: RCNTRL(20)
192 INTEGER, INTENT(INOUT) :: ISTATUS(20)
193 KPP_REAL, INTENT(INOUT) :: RSTATUS(20)
194 INTEGER, INTENT(OUT) :: IERR
195 !~~~> The method parameters
196 INTEGER, PARAMETER :: Smax = 6
197 INTEGER :: Method, ros_S
198 KPP_REAL, DIMENSION(Smax) :: ros_M, ros_E, ros_Alpha, ros_Gamma
199 KPP_REAL, DIMENSION(Smax*(Smax-1)/2) :: ros_A, ros_C
200 KPP_REAL :: ros_ELO
201 LOGICAL, DIMENSION(Smax) :: ros_NewF
202 CHARACTER(LEN=12) :: ros_Name
203 !~~~> Local variables
204 KPP_REAL :: Roundoff, FacMin, FacMax, FacRej, FacSafe
205 KPP_REAL :: Hmin, Hmax, Hstart, Hexit
206 KPP_REAL :: Texit
207 INTEGER :: i, UplimTol, Max_no_steps
208 LOGICAL :: Autonomous, VectorTol
209 !~~~> Parameters
210 KPP_REAL, PARAMETER :: ZERO = 0.0_dp, ONE = 1.0_dp
211 KPP_REAL, PARAMETER :: DeltaMin = 1.0E-5_dp
213 !~~~> Initialize statistics
214 Nfun = ISTATUS(ifun)
215 Njac = ISTATUS(ijac)
216 Nstp = ISTATUS(istp)
217 Nacc = ISTATUS(iacc)
218 Nrej = ISTATUS(irej)
219 Ndec = ISTATUS(idec)
220 Nsol = ISTATUS(isol)
221 Nsng = ISTATUS(isng)
223 !~~~> Autonomous or time dependent ODE. Default is time dependent.
224 Autonomous = .NOT.(ICNTRL(1) == 0)
226 !~~~> For Scalar tolerances (ICNTRL(2).NE.0) the code uses AbsTol(1) and RelTol(1)
227 ! For Vector tolerances (ICNTRL(2) == 0) the code uses AbsTol(1:NVAR) and RelTol(1:NVAR)
228 IF (ICNTRL(2) == 0) THEN
229 VectorTol = .TRUE.
230 UplimTol = NVAR
231 ELSE
232 VectorTol = .FALSE.
233 UplimTol = 1
234 END IF
236 !~~~> The particular Rosenbrock method chosen
237 IF (ICNTRL(3) == 0) THEN
238 Method = 4
239 ELSEIF ( (ICNTRL(3) >= 1).AND.(ICNTRL(3) <= 5) ) THEN
240 Method = ICNTRL(3)
241 ELSE
242 PRINT * , 'User-selected Rosenbrock method: ICNTRL(3)=', Method
243 CALL ros_ErrorMsg(-2,Tstart,ZERO,IERR)
244 RETURN
245 END IF
247 !~~~> The maximum number of steps admitted
248 IF (ICNTRL(4) == 0) THEN
249 Max_no_steps = 100000
250 ELSEIF (ICNTRL(4) > 0) THEN
251 Max_no_steps=ICNTRL(4)
252 ELSE
253 PRINT * ,'User-selected max no. of steps: ICNTRL(4)=',ICNTRL(4)
254 CALL ros_ErrorMsg(-1,Tstart,ZERO,IERR)
255 RETURN
256 END IF
258 !~~~> Unit roundoff (1+Roundoff>1)
259 Roundoff = WLAMCH('E')
261 !~~~> Lower bound on the step size: (positive value)
262 IF (RCNTRL(1) == ZERO) THEN
263 Hmin = ZERO
264 ELSEIF (RCNTRL(1) > ZERO) THEN
265 Hmin = RCNTRL(1)
266 ELSE
267 PRINT * , 'User-selected Hmin: RCNTRL(1)=', RCNTRL(1)
268 CALL ros_ErrorMsg(-3,Tstart,ZERO,IERR)
269 RETURN
270 END IF
271 !~~~> Upper bound on the step size: (positive value)
272 IF (RCNTRL(2) == ZERO) THEN
273 Hmax = ABS(Tend-Tstart)
274 ELSEIF (RCNTRL(2) > ZERO) THEN
275 Hmax = MIN(ABS(RCNTRL(2)),ABS(Tend-Tstart))
276 ELSE
277 PRINT * , 'User-selected Hmax: RCNTRL(2)=', RCNTRL(2)
278 CALL ros_ErrorMsg(-3,Tstart,ZERO,IERR)
279 RETURN
280 END IF
281 !~~~> Starting step size: (positive value)
282 IF (RCNTRL(3) == ZERO) THEN
283 Hstart = MAX(Hmin,DeltaMin)
284 ELSEIF (RCNTRL(3) > ZERO) THEN
285 Hstart = MIN(ABS(RCNTRL(3)),ABS(Tend-Tstart))
286 ELSE
287 PRINT * , 'User-selected Hstart: RCNTRL(3)=', RCNTRL(3)
288 CALL ros_ErrorMsg(-3,Tstart,ZERO,IERR)
289 RETURN
290 END IF
291 !~~~> Step size can be changed s.t. FacMin < Hnew/Hexit < FacMax
292 IF (RCNTRL(4) == ZERO) THEN
293 FacMin = 0.2_dp
294 ELSEIF (RCNTRL(4) > ZERO) THEN
295 FacMin = RCNTRL(4)
296 ELSE
297 PRINT * , 'User-selected FacMin: RCNTRL(4)=', RCNTRL(4)
298 CALL ros_ErrorMsg(-4,Tstart,ZERO,IERR)
299 RETURN
300 END IF
301 IF (RCNTRL(5) == ZERO) THEN
302 FacMax = 6.0_dp
303 ELSEIF (RCNTRL(5) > ZERO) THEN
304 FacMax = RCNTRL(5)
305 ELSE
306 PRINT * , 'User-selected FacMax: RCNTRL(5)=', RCNTRL(5)
307 CALL ros_ErrorMsg(-4,Tstart,ZERO,IERR)
308 RETURN
309 END IF
310 !~~~> FacRej: Factor to decrease step after 2 succesive rejections
311 IF (RCNTRL(6) == ZERO) THEN
312 FacRej = 0.1_dp
313 ELSEIF (RCNTRL(6) > ZERO) THEN
314 FacRej = RCNTRL(6)
315 ELSE
316 PRINT * , 'User-selected FacRej: RCNTRL(6)=', RCNTRL(6)
317 CALL ros_ErrorMsg(-4,Tstart,ZERO,IERR)
318 RETURN
319 END IF
320 !~~~> FacSafe: Safety Factor in the computation of new step size
321 IF (RCNTRL(7) == ZERO) THEN
322 FacSafe = 0.9_dp
323 ELSEIF (RCNTRL(7) > ZERO) THEN
324 FacSafe = RCNTRL(7)
325 ELSE
326 PRINT * , 'User-selected FacSafe: RCNTRL(7)=', RCNTRL(7)
327 CALL ros_ErrorMsg(-4,Tstart,ZERO,IERR)
328 RETURN
329 END IF
330 !~~~> Check if tolerances are reasonable
331 DO i=1,UplimTol
332 IF ( (AbsTol(i) <= ZERO) .OR. (RelTol(i) <= 10.0_dp*Roundoff) &
333 .OR. (RelTol(i) >= 1.0_dp) ) THEN
334 PRINT * , ' AbsTol(',i,') = ',AbsTol(i)
335 PRINT * , ' RelTol(',i,') = ',RelTol(i)
336 CALL ros_ErrorMsg(-5,Tstart,ZERO,IERR)
337 RETURN
338 END IF
339 END DO
342 !~~~> Initialize the particular Rosenbrock method
343 SELECT CASE (Method)
344 CASE (1)
345 CALL Ros2(ros_S, ros_A, ros_C, ros_M, ros_E, &
346 ros_Alpha, ros_Gamma, ros_NewF, ros_ELO, ros_Name)
347 CASE (2)
348 CALL Ros3(ros_S, ros_A, ros_C, ros_M, ros_E, &
349 ros_Alpha, ros_Gamma, ros_NewF, ros_ELO, ros_Name)
350 CASE (3)
351 CALL Ros4(ros_S, ros_A, ros_C, ros_M, ros_E, &
352 ros_Alpha, ros_Gamma, ros_NewF, ros_ELO, ros_Name)
353 CASE (4)
354 CALL Rodas3(ros_S, ros_A, ros_C, ros_M, ros_E, &
355 ros_Alpha, ros_Gamma, ros_NewF, ros_ELO, ros_Name)
356 CASE (5)
357 CALL Rodas4(ros_S, ros_A, ros_C, ros_M, ros_E, &
358 ros_Alpha, ros_Gamma, ros_NewF, ros_ELO, ros_Name)
359 CASE DEFAULT
360 PRINT * , 'Unknown Rosenbrock method: ICNTRL(4)=', Method
361 CALL ros_ErrorMsg(-2,Tstart,ZERO,IERR)
362 RETURN
363 END SELECT
365 !~~~> CALL Rosenbrock method
366 CALL ros_Integrator(Y,Tstart,Tend,Texit, &
367 AbsTol, RelTol, &
368 ! Rosenbrock method coefficients
369 ros_S, ros_M, ros_E, ros_A, ros_C, &
370 ros_Alpha, ros_Gamma, ros_ELO, ros_NewF, &
371 ! Integration parameters
372 Autonomous, VectorTol, Max_no_steps, &
373 Roundoff, Hmin, Hmax, Hstart, Hexit, &
374 FacMin, FacMax, FacRej, FacSafe, &
375 ! Error indicator
376 IERR)
379 !~~~> Collect run statistics
380 ISTATUS(ifun) = Nfun
381 ISTATUS(ijac) = Njac
382 ISTATUS(istp) = Nstp
383 ISTATUS(iacc) = Nacc
384 ISTATUS(irej) = Nrej
385 ISTATUS(idec) = Ndec
386 ISTATUS(isol) = Nsol
387 ISTATUS(isng) = Nsng
388 !~~~> Last T and H
389 RSTATUS(itexit) = Texit
390 RSTATUS(ihexit) = Hexit
392 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
393 CONTAINS ! SUBROUTINES internal to Rosenbrock
394 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
396 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
397 SUBROUTINE ros_ErrorMsg(Code,T,H,IERR)
398 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
399 ! Handles all error messages
400 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
402 KPP_REAL, INTENT(IN) :: T, H
403 INTEGER, INTENT(IN) :: Code
404 INTEGER, INTENT(OUT) :: IERR
406 IERR = Code
407 PRINT * , &
408 'Forced exit from Rosenbrock due to the following error:'
409 IF ((Code>=-8).AND.(Code<=-1)) THEN
410 PRINT *, IERR_NAMES(Code)
411 ELSE
412 PRINT *, 'Unknown Error code: ', Code
413 ENDIF
415 PRINT *, "T=", T, "and H=", H
417 END SUBROUTINE ros_ErrorMsg
419 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
420 SUBROUTINE ros_Integrator (Y, Tstart, Tend, T, &
421 AbsTol, RelTol, &
422 !~~~> Rosenbrock method coefficients
423 ros_S, ros_M, ros_E, ros_A, ros_C, &
424 ros_Alpha, ros_Gamma, ros_ELO, ros_NewF, &
425 !~~~> Integration parameters
426 Autonomous, VectorTol, Max_no_steps, &
427 Roundoff, Hmin, Hmax, Hstart, Hexit, &
428 FacMin, FacMax, FacRej, FacSafe, &
429 !~~~> Error indicator
430 IERR )
431 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
432 ! Template for the implementation of a generic Rosenbrock method
433 ! defined by ros_S (no of stages)
434 ! and its coefficients ros_{A,C,M,E,Alpha,Gamma}
435 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
437 IMPLICIT NONE
439 !~~~> Input: the initial condition at Tstart; Output: the solution at T
440 KPP_REAL, INTENT(INOUT) :: Y(NVAR)
441 !~~~> Input: integration interval
442 KPP_REAL, INTENT(IN) :: Tstart,Tend
443 !~~~> Output: time at which the solution is returned (T=Tend if success)
444 KPP_REAL, INTENT(OUT) :: T
445 !~~~> Input: tolerances
446 KPP_REAL, INTENT(IN) :: AbsTol(NVAR), RelTol(NVAR)
447 !~~~> Input: The Rosenbrock method parameters
448 INTEGER, INTENT(IN) :: ros_S
449 KPP_REAL, INTENT(IN) :: ros_M(ros_S), ros_E(ros_S), &
450 ros_Alpha(ros_S), ros_A(ros_S*(ros_S-1)/2), &
451 ros_Gamma(ros_S), ros_C(ros_S*(ros_S-1)/2), ros_ELO
452 LOGICAL, INTENT(IN) :: ros_NewF(ros_S)
453 !~~~> Input: integration parameters
454 LOGICAL, INTENT(IN) :: Autonomous, VectorTol
455 KPP_REAL, INTENT(IN) :: Hstart, Hmin, Hmax
456 INTEGER, INTENT(IN) :: Max_no_steps
457 KPP_REAL, INTENT(IN) :: Roundoff, FacMin, FacMax, FacRej, FacSafe
458 !~~~> Output: last accepted step
459 KPP_REAL, INTENT(OUT) :: Hexit
460 !~~~> Output: Error indicator
461 INTEGER, INTENT(OUT) :: IERR
462 ! ~~~~ Local variables
463 KPP_REAL :: Ynew(NVAR), Fcn0(NVAR), Fcn(NVAR)
464 KPP_REAL :: K(NVAR*ros_S), dFdT(NVAR)
465 #ifdef FULL_ALGEBRA
466 KPP_REAL :: Jac0(NVAR,NVAR), Ghimj(NVAR,NVAR)
467 #else
468 KPP_REAL :: Jac0(LU_NONZERO), Ghimj(LU_NONZERO)
469 #endif
470 KPP_REAL :: H, Hnew, HC, HG, Fac, Tau
471 KPP_REAL :: Err, Yerr(NVAR)
472 INTEGER :: Pivot(NVAR), Direction, ioffset, j, istage
473 LOGICAL :: RejectLastH, RejectMoreH, Singular
474 !~~~> Local parameters
475 KPP_REAL, PARAMETER :: ZERO = 0.0_dp, ONE = 1.0_dp
476 KPP_REAL, PARAMETER :: DeltaMin = 1.0E-5_dp
477 !~~~> Locally called functions
478 ! KPP_REAL WLAMCH
479 ! EXTERNAL WLAMCH
480 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
483 !~~~> Initial preparations
484 T = Tstart
485 Hexit = 0.0_dp
486 H = MIN(Hstart,Hmax)
487 IF (ABS(H) <= 10.0_dp*Roundoff) H = DeltaMin
489 IF (Tend >= Tstart) THEN
490 Direction = +1
491 ELSE
492 Direction = -1
493 END IF
495 RejectLastH=.FALSE.
496 RejectMoreH=.FALSE.
498 !~~~> Time loop begins below
500 TimeLoop: DO WHILE ( (Direction > 0).AND.((T-Tend)+Roundoff <= ZERO) &
501 .OR. (Direction < 0).AND.((Tend-T)+Roundoff <= ZERO) )
503 IF ( Nstp > Max_no_steps ) THEN ! Too many steps
504 CALL ros_ErrorMsg(-6,T,H,IERR)
505 RETURN
506 END IF
507 IF ( ((T+0.1_dp*H) == T).OR.(H <= Roundoff) ) THEN ! Step size too small
508 CALL ros_ErrorMsg(-7,T,H,IERR)
509 RETURN
510 END IF
512 !~~~> Limit H if necessary to avoid going beyond Tend
513 Hexit = H
514 H = MIN(H,ABS(Tend-T))
516 !~~~> Compute the function at current time
517 CALL FunTemplate(T,Y,Fcn0)
519 !~~~> Compute the function derivative with respect to T
520 IF (.NOT.Autonomous) THEN
521 CALL ros_FunTimeDerivative ( T, Roundoff, Y, &
522 Fcn0, dFdT )
523 END IF
525 !~~~> Compute the Jacobian at current time
526 CALL JacTemplate(T,Y,Jac0)
528 !~~~> Repeat step calculation until current step accepted
529 UntilAccepted: DO
531 CALL ros_PrepareMatrix(H,Direction,ros_Gamma(1), &
532 Jac0,Ghimj,Pivot,Singular)
533 IF (Singular) THEN ! More than 5 consecutive failed decompositions
534 CALL ros_ErrorMsg(-8,T,H,IERR)
535 RETURN
536 END IF
538 !~~~> Compute the stages
539 Stage: DO istage = 1, ros_S
541 ! Current istage offset. Current istage vector is K(ioffset+1:ioffset+NVAR)
542 ioffset = NVAR*(istage-1)
544 ! For the 1st istage the function has been computed previously
545 IF ( istage == 1 ) THEN
546 CALL WCOPY(NVAR,Fcn0,1,Fcn,1)
547 ! istage>1 and a new function evaluation is needed at the current istage
548 ELSEIF ( ros_NewF(istage) ) THEN
549 CALL WCOPY(NVAR,Y,1,Ynew,1)
550 DO j = 1, istage-1
551 CALL WAXPY(NVAR,ros_A((istage-1)*(istage-2)/2+j), &
552 K(NVAR*(j-1)+1),1,Ynew,1)
553 END DO
554 Tau = T + ros_Alpha(istage)*Direction*H
555 CALL FunTemplate(Tau,Ynew,Fcn)
556 END IF ! if istage == 1 elseif ros_NewF(istage)
557 CALL WCOPY(NVAR,Fcn,1,K(ioffset+1),1)
558 DO j = 1, istage-1
559 HC = ros_C((istage-1)*(istage-2)/2+j)/(Direction*H)
560 CALL WAXPY(NVAR,HC,K(NVAR*(j-1)+1),1,K(ioffset+1),1)
561 END DO
562 IF ((.NOT. Autonomous).AND.(ros_Gamma(istage).NE.ZERO)) THEN
563 HG = Direction*H*ros_Gamma(istage)
564 CALL WAXPY(NVAR,HG,dFdT,1,K(ioffset+1),1)
565 END IF
566 CALL ros_Solve(Ghimj, Pivot, K(ioffset+1))
568 END DO Stage
571 !~~~> Compute the new solution
572 CALL WCOPY(NVAR,Y,1,Ynew,1)
573 DO j=1,ros_S
574 CALL WAXPY(NVAR,ros_M(j),K(NVAR*(j-1)+1),1,Ynew,1)
575 END DO
577 !~~~> Compute the error estimation
578 CALL WSCAL(NVAR,ZERO,Yerr,1)
579 DO j=1,ros_S
580 CALL WAXPY(NVAR,ros_E(j),K(NVAR*(j-1)+1),1,Yerr,1)
581 END DO
582 Err = ros_ErrorNorm ( Y, Ynew, Yerr, AbsTol, RelTol, VectorTol )
584 !~~~> New step size is bounded by FacMin <= Hnew/H <= FacMax
585 Fac = MIN(FacMax,MAX(FacMin,FacSafe/Err**(ONE/ros_ELO)))
586 Hnew = H*Fac
588 !~~~> Check the error magnitude and adjust step size
589 Nstp = Nstp+1
590 IF ( (Err <= ONE).OR.(H <= Hmin) ) THEN !~~~> Accept step
591 Nacc = Nacc+1
592 CALL WCOPY(NVAR,Ynew,1,Y,1)
593 T = T + Direction*H
594 Hnew = MAX(Hmin,MIN(Hnew,Hmax))
595 IF (RejectLastH) THEN ! No step size increase after a rejected step
596 Hnew = MIN(Hnew,H)
597 END IF
598 RejectLastH = .FALSE.
599 RejectMoreH = .FALSE.
600 H = Hnew
601 EXIT UntilAccepted ! EXIT THE LOOP: WHILE STEP NOT ACCEPTED
602 ELSE !~~~> Reject step
603 IF (RejectMoreH) THEN
604 Hnew = H*FacRej
605 END IF
606 RejectMoreH = RejectLastH
607 RejectLastH = .TRUE.
608 H = Hnew
609 IF (Nacc >= 1) THEN
610 Nrej = Nrej+1
611 END IF
612 END IF ! Err <= 1
614 END DO UntilAccepted
616 END DO TimeLoop
618 !~~~> Succesful exit
619 IERR = 1 !~~~> The integration was successful
621 END SUBROUTINE ros_Integrator
624 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
625 KPP_REAL FUNCTION ros_ErrorNorm ( Y, Ynew, Yerr, &
626 AbsTol, RelTol, VectorTol )
627 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
628 !~~~> Computes the "scaled norm" of the error vector Yerr
629 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
630 IMPLICIT NONE
632 ! Input arguments
633 KPP_REAL, INTENT(IN) :: Y(NVAR), Ynew(NVAR), &
634 Yerr(NVAR), AbsTol(NVAR), RelTol(NVAR)
635 LOGICAL, INTENT(IN) :: VectorTol
636 ! Local variables
637 KPP_REAL :: Err, Scale, Ymax
638 INTEGER :: i
639 KPP_REAL, PARAMETER :: ZERO = 0.0_dp
641 Err = ZERO
642 DO i=1,NVAR
643 Ymax = MAX(ABS(Y(i)),ABS(Ynew(i)))
644 IF (VectorTol) THEN
645 Scale = AbsTol(i)+RelTol(i)*Ymax
646 ELSE
647 Scale = AbsTol(1)+RelTol(1)*Ymax
648 END IF
649 Err = Err+(Yerr(i)/Scale)**2
650 END DO
651 Err = SQRT(Err/NVAR)
653 ros_ErrorNorm = Err
655 END FUNCTION ros_ErrorNorm
658 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
659 SUBROUTINE ros_FunTimeDerivative ( T, Roundoff, Y, &
660 Fcn0, dFdT )
661 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
662 !~~~> The time partial derivative of the function by finite differences
663 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
664 IMPLICIT NONE
666 !~~~> Input arguments
667 KPP_REAL, INTENT(IN) :: T, Roundoff, Y(NVAR), Fcn0(NVAR)
668 !~~~> Output arguments
669 KPP_REAL, INTENT(OUT) :: dFdT(NVAR)
670 !~~~> Local variables
671 KPP_REAL :: Delta
672 KPP_REAL, PARAMETER :: ONE = 1.0_dp, DeltaMin = 1.0E-6_dp
674 Delta = SQRT(Roundoff)*MAX(DeltaMin,ABS(T))
675 CALL FunTemplate(T+Delta,Y,dFdT)
676 CALL WAXPY(NVAR,(-ONE),Fcn0,1,dFdT,1)
677 CALL WSCAL(NVAR,(ONE/Delta),dFdT,1)
679 END SUBROUTINE ros_FunTimeDerivative
682 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
683 SUBROUTINE ros_PrepareMatrix ( H, Direction, gam, &
684 Jac0, Ghimj, Pivot, Singular )
685 ! --- --- --- --- --- --- --- --- --- --- --- --- ---
686 ! Prepares the LHS matrix for stage calculations
687 ! 1. Construct Ghimj = 1/(H*ham) - Jac0
688 ! "(Gamma H) Inverse Minus Jacobian"
689 ! 2. Repeat LU decomposition of Ghimj until successful.
690 ! -half the step size if LU decomposition fails and retry
691 ! -exit after 5 consecutive fails
692 ! --- --- --- --- --- --- --- --- --- --- --- --- ---
693 IMPLICIT NONE
695 !~~~> Input arguments
696 #ifdef FULL_ALGEBRA
697 KPP_REAL, INTENT(IN) :: Jac0(NVAR,NVAR)
698 #else
699 KPP_REAL, INTENT(IN) :: Jac0(LU_NONZERO)
700 #endif
701 KPP_REAL, INTENT(IN) :: gam
702 INTEGER, INTENT(IN) :: Direction
703 !~~~> Output arguments
704 #ifdef FULL_ALGEBRA
705 KPP_REAL, INTENT(OUT) :: Ghimj(NVAR,NVAR)
706 #else
707 KPP_REAL, INTENT(OUT) :: Ghimj(LU_NONZERO)
708 #endif
709 LOGICAL, INTENT(OUT) :: Singular
710 INTEGER, INTENT(OUT) :: Pivot(NVAR)
711 !~~~> Inout arguments
712 KPP_REAL, INTENT(INOUT) :: H ! step size is decreased when LU fails
713 !~~~> Local variables
714 INTEGER :: i, ising, Nconsecutive
715 KPP_REAL :: ghinv
716 KPP_REAL, PARAMETER :: ONE = 1.0_dp, HALF = 0.5_dp
718 Nconsecutive = 0
719 Singular = .TRUE.
721 DO WHILE (Singular)
723 !~~~> Construct Ghimj = 1/(H*gam) - Jac0
724 #ifdef FULL_ALGEBRA
725 CALL WCOPY(NVAR*NVAR,Jac0,1,Ghimj,1)
726 CALL WSCAL(NVAR*NVAR,(-ONE),Ghimj,1)
727 ghinv = ONE/(Direction*H*gam)
728 DO i=1,NVAR
729 Ghimj(i,i) = Ghimj(i,i)+ghinv
730 END DO
731 #else
732 CALL WCOPY(LU_NONZERO,Jac0,1,Ghimj,1)
733 CALL WSCAL(LU_NONZERO,(-ONE),Ghimj,1)
734 ghinv = ONE/(Direction*H*gam)
735 DO i=1,NVAR
736 Ghimj(LU_DIAG(i)) = Ghimj(LU_DIAG(i))+ghinv
737 END DO
738 #endif
739 !~~~> Compute LU decomposition
740 CALL ros_Decomp( Ghimj, Pivot, ising )
741 IF (ising == 0) THEN
742 !~~~> If successful done
743 Singular = .FALSE.
744 ELSE ! ising .ne. 0
745 !~~~> If unsuccessful half the step size; if 5 consecutive fails then return
746 Nsng = Nsng+1
747 Nconsecutive = Nconsecutive+1
748 Singular = .TRUE.
749 PRINT*,'Warning: LU Decomposition returned ising = ',ising
750 IF (Nconsecutive <= 5) THEN ! Less than 5 consecutive failed decompositions
751 H = H*HALF
752 ELSE ! More than 5 consecutive failed decompositions
753 RETURN
754 END IF ! Nconsecutive
755 END IF ! ising
757 END DO ! WHILE Singular
759 END SUBROUTINE ros_PrepareMatrix
762 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
763 SUBROUTINE ros_Decomp( A, Pivot, ising )
764 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
765 ! Template for the LU decomposition
766 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
767 IMPLICIT NONE
768 !~~~> Inout variables
769 KPP_REAL, INTENT(INOUT) :: A(LU_NONZERO)
770 !~~~> Output variables
771 INTEGER, INTENT(OUT) :: Pivot(NVAR), ising
773 #ifdef FULL_ALGEBRA
774 CALL DGETRF( NVAR, NVAR, A, NVAR, Pivot, ising )
775 #else
776 CALL KppDecomp ( A, ising )
777 Pivot(1) = 1
778 #endif
779 Ndec = Ndec + 1
781 END SUBROUTINE ros_Decomp
784 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
785 SUBROUTINE ros_Solve( A, Pivot, b )
786 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
787 ! Template for the forward/backward substitution (using pre-computed LU decomposition)
788 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
789 IMPLICIT NONE
790 !~~~> Input variables
791 KPP_REAL, INTENT(IN) :: A(LU_NONZERO)
792 INTEGER, INTENT(IN) :: Pivot(NVAR)
793 !~~~> InOut variables
794 KPP_REAL, INTENT(INOUT) :: b(NVAR)
796 #ifdef FULL_ALGEBRA
797 CALL DGETRS( 'N', NVAR , 1, A, NVAR, Pivot, b, NVAR, 0 )
798 #else
799 CALL KppSolve( A, b )
800 #endif
802 Nsol = Nsol+1
804 END SUBROUTINE ros_Solve
808 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
809 SUBROUTINE Ros2 (ros_S,ros_A,ros_C,ros_M,ros_E,ros_Alpha,&
810 ros_Gamma,ros_NewF,ros_ELO,ros_Name)
811 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
812 ! --- AN L-STABLE METHOD, 2 stages, order 2
813 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
815 IMPLICIT NONE
817 INTEGER, PARAMETER :: S=2
818 INTEGER, INTENT(OUT) :: ros_S
819 KPP_REAL, DIMENSION(S), INTENT(OUT) :: ros_M,ros_E,ros_Alpha,ros_Gamma
820 KPP_REAL, DIMENSION(S*(S-1)/2), INTENT(OUT) :: ros_A, ros_C
821 KPP_REAL, INTENT(OUT) :: ros_ELO
822 LOGICAL, DIMENSION(S), INTENT(OUT) :: ros_NewF
823 CHARACTER(LEN=12), INTENT(OUT) :: ros_Name
824 DOUBLE PRECISION g
826 g = 1.0_dp + 1.0_dp/SQRT(2.0_dp)
828 !~~~> Name of the method
829 ros_Name = 'ROS-2'
830 !~~~> Number of stages
831 ros_S = S
833 !~~~> The coefficient matrices A and C are strictly lower triangular.
834 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
835 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
836 ! The general mapping formula is:
837 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
838 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
840 ros_A(1) = (1.0_dp)/g
841 ros_C(1) = (-2.0_dp)/g
842 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
843 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
844 ros_NewF(1) = .TRUE.
845 ros_NewF(2) = .TRUE.
846 !~~~> M_i = Coefficients for new step solution
847 ros_M(1)= (3.0_dp)/(2.0_dp*g)
848 ros_M(2)= (1.0_dp)/(2.0_dp*g)
849 ! E_i = Coefficients for error estimator
850 ros_E(1) = 1.0_dp/(2.0_dp*g)
851 ros_E(2) = 1.0_dp/(2.0_dp*g)
852 !~~~> ros_ELO = estimator of local order - the minimum between the
853 ! main and the embedded scheme orders plus one
854 ros_ELO = 2.0_dp
855 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
856 ros_Alpha(1) = 0.0_dp
857 ros_Alpha(2) = 1.0_dp
858 !~~~> Gamma_i = \sum_j gamma_{i,j}
859 ros_Gamma(1) = g
860 ros_Gamma(2) =-g
862 END SUBROUTINE Ros2
865 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
866 SUBROUTINE Ros3 (ros_S,ros_A,ros_C,ros_M,ros_E,ros_Alpha,&
867 ros_Gamma,ros_NewF,ros_ELO,ros_Name)
868 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
869 ! --- AN L-STABLE METHOD, 3 stages, order 3, 2 function evaluations
870 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
872 IMPLICIT NONE
874 INTEGER, PARAMETER :: S=3
875 INTEGER, INTENT(OUT) :: ros_S
876 KPP_REAL, DIMENSION(S), INTENT(OUT) :: ros_M,ros_E,ros_Alpha,ros_Gamma
877 KPP_REAL, DIMENSION(S*(S-1)/2), INTENT(OUT) :: ros_A, ros_C
878 KPP_REAL, INTENT(OUT) :: ros_ELO
879 LOGICAL, DIMENSION(S), INTENT(OUT) :: ros_NewF
880 CHARACTER(LEN=12), INTENT(OUT) :: ros_Name
882 !~~~> Name of the method
883 ros_Name = 'ROS-3'
884 !~~~> Number of stages
885 ros_S = S
887 !~~~> The coefficient matrices A and C are strictly lower triangular.
888 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
889 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
890 ! The general mapping formula is:
891 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
892 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
894 ros_A(1)= 1.0_dp
895 ros_A(2)= 1.0_dp
896 ros_A(3)= 0.0_dp
898 ros_C(1) = -0.10156171083877702091975600115545E+01_dp
899 ros_C(2) = 0.40759956452537699824805835358067E+01_dp
900 ros_C(3) = 0.92076794298330791242156818474003E+01_dp
901 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
902 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
903 ros_NewF(1) = .TRUE.
904 ros_NewF(2) = .TRUE.
905 ros_NewF(3) = .FALSE.
906 !~~~> M_i = Coefficients for new step solution
907 ros_M(1) = 0.1E+01_dp
908 ros_M(2) = 0.61697947043828245592553615689730E+01_dp
909 ros_M(3) = -0.42772256543218573326238373806514E+00_dp
910 ! E_i = Coefficients for error estimator
911 ros_E(1) = 0.5E+00_dp
912 ros_E(2) = -0.29079558716805469821718236208017E+01_dp
913 ros_E(3) = 0.22354069897811569627360909276199E+00_dp
914 !~~~> ros_ELO = estimator of local order - the minimum between the
915 ! main and the embedded scheme orders plus 1
916 ros_ELO = 3.0_dp
917 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
918 ros_Alpha(1)= 0.0E+00_dp
919 ros_Alpha(2)= 0.43586652150845899941601945119356E+00_dp
920 ros_Alpha(3)= 0.43586652150845899941601945119356E+00_dp
921 !~~~> Gamma_i = \sum_j gamma_{i,j}
922 ros_Gamma(1)= 0.43586652150845899941601945119356E+00_dp
923 ros_Gamma(2)= 0.24291996454816804366592249683314E+00_dp
924 ros_Gamma(3)= 0.21851380027664058511513169485832E+01_dp
926 END SUBROUTINE Ros3
928 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
931 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
932 SUBROUTINE Ros4 (ros_S,ros_A,ros_C,ros_M,ros_E,ros_Alpha,&
933 ros_Gamma,ros_NewF,ros_ELO,ros_Name)
934 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
935 ! L-STABLE ROSENBROCK METHOD OF ORDER 4, WITH 4 STAGES
936 ! L-STABLE EMBEDDED ROSENBROCK METHOD OF ORDER 3
938 ! E. HAIRER AND G. WANNER, SOLVING ORDINARY DIFFERENTIAL
939 ! EQUATIONS II. STIFF AND DIFFERENTIAL-ALGEBRAIC PROBLEMS.
940 ! SPRINGER SERIES IN COMPUTATIONAL MATHEMATICS,
941 ! SPRINGER-VERLAG (1990)
942 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
944 IMPLICIT NONE
946 INTEGER, PARAMETER :: S=4
947 INTEGER, INTENT(OUT) :: ros_S
948 KPP_REAL, DIMENSION(4), INTENT(OUT) :: ros_M,ros_E,ros_Alpha,ros_Gamma
949 KPP_REAL, DIMENSION(6), INTENT(OUT) :: ros_A, ros_C
950 KPP_REAL, INTENT(OUT) :: ros_ELO
951 LOGICAL, DIMENSION(4), INTENT(OUT) :: ros_NewF
952 CHARACTER(LEN=12), INTENT(OUT) :: ros_Name
953 DOUBLE PRECISION g
955 !~~~> Name of the method
956 ros_Name = 'ROS-4'
957 !~~~> Number of stages
958 ros_S = S
960 !~~~> The coefficient matrices A and C are strictly lower triangular.
961 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
962 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
963 ! The general mapping formula is:
964 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
965 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
967 ros_A(1) = 0.2000000000000000E+01_dp
968 ros_A(2) = 0.1867943637803922E+01_dp
969 ros_A(3) = 0.2344449711399156E+00_dp
970 ros_A(4) = ros_A(2)
971 ros_A(5) = ros_A(3)
972 ros_A(6) = 0.0_dp
974 ros_C(1) =-0.7137615036412310E+01_dp
975 ros_C(2) = 0.2580708087951457E+01_dp
976 ros_C(3) = 0.6515950076447975E+00_dp
977 ros_C(4) =-0.2137148994382534E+01_dp
978 ros_C(5) =-0.3214669691237626E+00_dp
979 ros_C(6) =-0.6949742501781779E+00_dp
980 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
981 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
982 ros_NewF(1) = .TRUE.
983 ros_NewF(2) = .TRUE.
984 ros_NewF(3) = .TRUE.
985 ros_NewF(4) = .FALSE.
986 !~~~> M_i = Coefficients for new step solution
987 ros_M(1) = 0.2255570073418735E+01_dp
988 ros_M(2) = 0.2870493262186792E+00_dp
989 ros_M(3) = 0.4353179431840180E+00_dp
990 ros_M(4) = 0.1093502252409163E+01_dp
991 !~~~> E_i = Coefficients for error estimator
992 ros_E(1) =-0.2815431932141155E+00_dp
993 ros_E(2) =-0.7276199124938920E-01_dp
994 ros_E(3) =-0.1082196201495311E+00_dp
995 ros_E(4) =-0.1093502252409163E+01_dp
996 !~~~> ros_ELO = estimator of local order - the minimum between the
997 ! main and the embedded scheme orders plus 1
998 ros_ELO = 4.0_dp
999 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
1000 ros_Alpha(1) = 0.0_dp
1001 ros_Alpha(2) = 0.1145640000000000E+01_dp
1002 ros_Alpha(3) = 0.6552168638155900E+00_dp
1003 ros_Alpha(4) = ros_Alpha(3)
1004 !~~~> Gamma_i = \sum_j gamma_{i,j}
1005 ros_Gamma(1) = 0.5728200000000000E+00_dp
1006 ros_Gamma(2) =-0.1769193891319233E+01_dp
1007 ros_Gamma(3) = 0.7592633437920482E+00_dp
1008 ros_Gamma(4) =-0.1049021087100450E+00_dp
1010 END SUBROUTINE Ros4
1012 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1013 SUBROUTINE Rodas3 (ros_S,ros_A,ros_C,ros_M,ros_E,ros_Alpha,&
1014 ros_Gamma,ros_NewF,ros_ELO,ros_Name)
1015 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1016 ! --- A STIFFLY-STABLE METHOD, 4 stages, order 3
1017 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1019 IMPLICIT NONE
1021 INTEGER, PARAMETER :: S=4
1022 INTEGER, INTENT(OUT) :: ros_S
1023 KPP_REAL, DIMENSION(S), INTENT(OUT) :: ros_M,ros_E,ros_Alpha,ros_Gamma
1024 KPP_REAL, DIMENSION(S*(S-1)/2), INTENT(OUT) :: ros_A, ros_C
1025 KPP_REAL, INTENT(OUT) :: ros_ELO
1026 LOGICAL, DIMENSION(S), INTENT(OUT) :: ros_NewF
1027 CHARACTER(LEN=12), INTENT(OUT) :: ros_Name
1028 DOUBLE PRECISION g
1030 !~~~> Name of the method
1031 ros_Name = 'RODAS-3'
1032 !~~~> Number of stages
1033 ros_S = S
1035 !~~~> The coefficient matrices A and C are strictly lower triangular.
1036 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
1037 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
1038 ! The general mapping formula is:
1039 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
1040 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
1042 ros_A(1) = 0.0E+00_dp
1043 ros_A(2) = 2.0E+00_dp
1044 ros_A(3) = 0.0E+00_dp
1045 ros_A(4) = 2.0E+00_dp
1046 ros_A(5) = 0.0E+00_dp
1047 ros_A(6) = 1.0E+00_dp
1049 ros_C(1) = 4.0E+00_dp
1050 ros_C(2) = 1.0E+00_dp
1051 ros_C(3) =-1.0E+00_dp
1052 ros_C(4) = 1.0E+00_dp
1053 ros_C(5) =-1.0E+00_dp
1054 ros_C(6) =-(8.0E+00_dp/3.0E+00_dp)
1056 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
1057 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
1058 ros_NewF(1) = .TRUE.
1059 ros_NewF(2) = .FALSE.
1060 ros_NewF(3) = .TRUE.
1061 ros_NewF(4) = .TRUE.
1062 !~~~> M_i = Coefficients for new step solution
1063 ros_M(1) = 2.0E+00_dp
1064 ros_M(2) = 0.0E+00_dp
1065 ros_M(3) = 1.0E+00_dp
1066 ros_M(4) = 1.0E+00_dp
1067 !~~~> E_i = Coefficients for error estimator
1068 ros_E(1) = 0.0E+00_dp
1069 ros_E(2) = 0.0E+00_dp
1070 ros_E(3) = 0.0E+00_dp
1071 ros_E(4) = 1.0E+00_dp
1072 !~~~> ros_ELO = estimator of local order - the minimum between the
1073 ! main and the embedded scheme orders plus 1
1074 ros_ELO = 3.0E+00_dp
1075 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
1076 ros_Alpha(1) = 0.0E+00_dp
1077 ros_Alpha(2) = 0.0E+00_dp
1078 ros_Alpha(3) = 1.0E+00_dp
1079 ros_Alpha(4) = 1.0E+00_dp
1080 !~~~> Gamma_i = \sum_j gamma_{i,j}
1081 ros_Gamma(1) = 0.5E+00_dp
1082 ros_Gamma(2) = 1.5E+00_dp
1083 ros_Gamma(3) = 0.0E+00_dp
1084 ros_Gamma(4) = 0.0E+00_dp
1086 END SUBROUTINE Rodas3
1088 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1089 SUBROUTINE Rodas4 (ros_S,ros_A,ros_C,ros_M,ros_E,ros_Alpha,&
1090 ros_Gamma,ros_NewF,ros_ELO,ros_Name)
1091 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1092 ! STIFFLY-STABLE ROSENBROCK METHOD OF ORDER 4, WITH 6 STAGES
1094 ! E. HAIRER AND G. WANNER, SOLVING ORDINARY DIFFERENTIAL
1095 ! EQUATIONS II. STIFF AND DIFFERENTIAL-ALGEBRAIC PROBLEMS.
1096 ! SPRINGER SERIES IN COMPUTATIONAL MATHEMATICS,
1097 ! SPRINGER-VERLAG (1996)
1098 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1100 IMPLICIT NONE
1102 INTEGER, PARAMETER :: S=6
1103 INTEGER, INTENT(OUT) :: ros_S
1104 KPP_REAL, DIMENSION(S), INTENT(OUT) :: ros_M,ros_E,ros_Alpha,ros_Gamma
1105 KPP_REAL, DIMENSION(S*(S-1)/2), INTENT(OUT) :: ros_A, ros_C
1106 KPP_REAL, INTENT(OUT) :: ros_ELO
1107 LOGICAL, DIMENSION(S), INTENT(OUT) :: ros_NewF
1108 CHARACTER(LEN=12), INTENT(OUT) :: ros_Name
1109 DOUBLE PRECISION g
1111 !~~~> Name of the method
1112 ros_Name = 'RODAS-4'
1113 !~~~> Number of stages
1114 ros_S = 6
1116 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
1117 ros_Alpha(1) = 0.000_dp
1118 ros_Alpha(2) = 0.386_dp
1119 ros_Alpha(3) = 0.210_dp
1120 ros_Alpha(4) = 0.630_dp
1121 ros_Alpha(5) = 1.000_dp
1122 ros_Alpha(6) = 1.000_dp
1124 !~~~> Gamma_i = \sum_j gamma_{i,j}
1125 ros_Gamma(1) = 0.2500000000000000E+00_dp
1126 ros_Gamma(2) =-0.1043000000000000E+00_dp
1127 ros_Gamma(3) = 0.1035000000000000E+00_dp
1128 ros_Gamma(4) =-0.3620000000000023E-01_dp
1129 ros_Gamma(5) = 0.0_dp
1130 ros_Gamma(6) = 0.0_dp
1132 !~~~> The coefficient matrices A and C are strictly lower triangular.
1133 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
1134 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
1135 ! The general mapping formula is: A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
1136 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
1138 ros_A(1) = 0.1544000000000000E+01_dp
1139 ros_A(2) = 0.9466785280815826E+00_dp
1140 ros_A(3) = 0.2557011698983284E+00_dp
1141 ros_A(4) = 0.3314825187068521E+01_dp
1142 ros_A(5) = 0.2896124015972201E+01_dp
1143 ros_A(6) = 0.9986419139977817E+00_dp
1144 ros_A(7) = 0.1221224509226641E+01_dp
1145 ros_A(8) = 0.6019134481288629E+01_dp
1146 ros_A(9) = 0.1253708332932087E+02_dp
1147 ros_A(10) =-0.6878860361058950E+00_dp
1148 ros_A(11) = ros_A(7)
1149 ros_A(12) = ros_A(8)
1150 ros_A(13) = ros_A(9)
1151 ros_A(14) = ros_A(10)
1152 ros_A(15) = 1.0E+00_dp
1154 ros_C(1) =-0.5668800000000000E+01_dp
1155 ros_C(2) =-0.2430093356833875E+01_dp
1156 ros_C(3) =-0.2063599157091915E+00_dp
1157 ros_C(4) =-0.1073529058151375E+00_dp
1158 ros_C(5) =-0.9594562251023355E+01_dp
1159 ros_C(6) =-0.2047028614809616E+02_dp
1160 ros_C(7) = 0.7496443313967647E+01_dp
1161 ros_C(8) =-0.1024680431464352E+02_dp
1162 ros_C(9) =-0.3399990352819905E+02_dp
1163 ros_C(10) = 0.1170890893206160E+02_dp
1164 ros_C(11) = 0.8083246795921522E+01_dp
1165 ros_C(12) =-0.7981132988064893E+01_dp
1166 ros_C(13) =-0.3152159432874371E+02_dp
1167 ros_C(14) = 0.1631930543123136E+02_dp
1168 ros_C(15) =-0.6058818238834054E+01_dp
1170 !~~~> M_i = Coefficients for new step solution
1171 ros_M(1) = ros_A(7)
1172 ros_M(2) = ros_A(8)
1173 ros_M(3) = ros_A(9)
1174 ros_M(4) = ros_A(10)
1175 ros_M(5) = 1.0E+00_dp
1176 ros_M(6) = 1.0E+00_dp
1178 !~~~> E_i = Coefficients for error estimator
1179 ros_E(1) = 0.0E+00_dp
1180 ros_E(2) = 0.0E+00_dp
1181 ros_E(3) = 0.0E+00_dp
1182 ros_E(4) = 0.0E+00_dp
1183 ros_E(5) = 0.0E+00_dp
1184 ros_E(6) = 1.0E+00_dp
1186 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
1187 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
1188 ros_NewF(1) = .TRUE.
1189 ros_NewF(2) = .TRUE.
1190 ros_NewF(3) = .TRUE.
1191 ros_NewF(4) = .TRUE.
1192 ros_NewF(5) = .TRUE.
1193 ros_NewF(6) = .TRUE.
1195 !~~~> ros_ELO = estimator of local order - the minimum between the
1196 ! main and the embedded scheme orders plus 1
1197 ros_ELO = 4.0_dp
1199 END SUBROUTINE Rodas4
1201 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1202 ! End of the set of internal Rosenbrock subroutines
1203 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1204 END SUBROUTINE Rosenbrock
1205 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1208 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1209 SUBROUTINE FunTemplate( T, Y, Ydot )
1210 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1211 ! Template for the ODE function call.
1212 ! Updates the rate coefficients (and possibly the fixed species) at each call
1213 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1214 USE KPP_ROOT_Parameters
1215 USE KPP_ROOT_Global
1216 USE KPP_ROOT_Function
1217 USE KPP_ROOT_Rates
1218 !~~~> Input variables
1219 KPP_REAL :: T, Y(NVAR)
1220 !~~~> Output variables
1221 KPP_REAL :: Ydot(NVAR)
1222 !~~~> Local variables
1223 KPP_REAL :: Told
1225 Told = TIME
1226 TIME = T
1227 CALL Update_SUN()
1228 CALL Update_RCONST()
1229 CALL Fun( Y, FIX, RCONST, Ydot )
1230 TIME = Told
1232 Nfun = Nfun+1
1234 END SUBROUTINE FunTemplate
1237 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1238 SUBROUTINE JacTemplate( T, Y, Jcb )
1239 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1240 ! Template for the ODE Jacobian call.
1241 ! Updates the rate coefficients (and possibly the fixed species) at each call
1242 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1243 USE KPP_ROOT_Parameters
1244 USE KPP_ROOT_Global
1245 USE KPP_ROOT_Jacobian
1246 USE KPP_ROOT_LinearAlgebra
1247 USE KPP_ROOT_Rates
1248 !~~~> Input variables
1249 KPP_REAL :: T, Y(NVAR)
1250 !~~~> Output variables
1251 #ifdef FULL_ALGEBRA
1252 KPP_REAL :: JV(LU_NONZERO), Jcb(NVAR,NVAR)
1253 #else
1254 KPP_REAL :: Jcb(LU_NONZERO)
1255 #endif
1256 !~~~> Local variables
1257 KPP_REAL :: Told
1258 #ifdef FULL_ALGEBRA
1259 INTEGER :: i, j
1260 #endif
1262 Told = TIME
1263 TIME = T
1264 CALL Update_SUN()
1265 CALL Update_RCONST()
1266 #ifdef FULL_ALGEBRA
1267 CALL Jac_SP(Y, FIX, RCONST, JV)
1268 DO j=1,NVAR
1269 DO i=1,NVAR
1270 Jcb(i,j) = 0.0d0
1271 END DO
1272 END DO
1273 DO i=1,LU_NONZERO
1274 Jcb(LU_IROW(i),LU_ICOL(i)) = JV(i)
1275 END DO
1276 #else
1277 CALL Jac_SP( Y, FIX, RCONST, Jcb )
1278 #endif
1279 TIME = Told
1281 Njac = Njac+1
1283 END SUBROUTINE JacTemplate
1285 END MODULE KPP_ROOT_Integrator