1 MODULE KPP_ROOT_Integrator
4 USE KPP_ROOT_Parameters
6 USE KPP_ROOT_LinearAlgebra
16 !~~~> Statistics on the work performed by the Rosenbrock method
17 INTEGER :: Nfun
,Njac
,Nstp
,Nacc
,Nrej
,Ndec
,Nsol
,Nsng
18 INTEGER, PARAMETER :: ifun
=1, ijac
=2, istp
=3, iacc
=4, &
19 irej
=5, idec
=6, isol
=7, isng
=8, &
22 KPP_REAL
, PARAMETER :: ZERO
= 0.0d0, ONE
= 1.0d0
25 CONTAINS ! Functions in the module KPP_ROOT_Integrator
27 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
28 SUBROUTINE INTEGRATE_TLM( NTLM
, Y
, Y_tlm
, TIN
, TOUT
, &
29 ICNTRL_U
, RCNTRL_U
, ISTATUS_U
, RSTATUS_U
)
30 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
33 !~~~> Y - Concentrations
35 !~~~> NTLM - No. of sensitivity coefficients
37 !~~~> Y_tlm - Sensitivities of concentrations
38 ! Note: Y_tlm (1:NVAR,j) contains sensitivities of
39 ! Y(1:NVAR) w.r.t. the j-th parameter, j=1...NTLM
40 KPP_REAL
:: Y_tlm(NVAR
,NTLM
)
41 KPP_REAL
, INTENT(IN
) :: TIN
! TIN - Start Time
42 KPP_REAL
, INTENT(IN
) :: TOUT
! TOUT - End Time
43 !~~~> Optional input parameters and statistics
44 INTEGER, INTENT(IN
), OPTIONAL
:: ICNTRL_U(20)
45 KPP_REAL
, INTENT(IN
), OPTIONAL
:: RCNTRL_U(20)
46 INTEGER, INTENT(OUT
), OPTIONAL
:: ISTATUS_U(20)
47 KPP_REAL
, INTENT(OUT
), OPTIONAL
:: RSTATUS_U(20)
49 INTEGER, SAVE :: N_stp
, N_acc
, N_rej
, N_sng
, IERR
50 KPP_REAL
:: RCNTRL(20), RSTATUS(20)
51 INTEGER :: ICNTRL(20), ISTATUS(20)
56 RSTATUS(1:10) = 0.0_dp
58 ICNTRL(1) = 0 ! non-autonomous
59 ICNTRL(2) = 1 ! vector tolerances
60 RCNTRL(3) = STEPMIN
! starting step
61 ICNTRL(4) = 5 ! choice of the method
63 ! if optional parameters are given, and if they are >=0, then they overwrite default settings
64 IF (PRESENT(ICNTRL_U
)) THEN
65 WHERE(ICNTRL_U(:) >= 0) ICNTRL(1:20) = ICNTRL_U(:)
67 IF (PRESENT(RCNTRL_U
)) THEN
68 WHERE(RCNTRL_U(:) >= 0) RCNTRL(1:20) = RCNTRL_U(:)
71 CALL RosenbrockTLM(VAR
, NTLM
, Y_tlm
, &
74 FunTemplate
,JacTemplate
,HessTemplate
, &
75 RCNTRL
,ICNTRL
,RSTATUS
,ISTATUS
,IERR
)
77 ! N_stp = N_stp + ICNTRL(istp)
78 ! N_acc = N_acc + ICNTRL(iacc)
79 ! N_rej = N_rej + ICNTRL(irej)
80 ! N_sng = N_sng + ICNTRL(isng)
81 ! PRINT*,'Step=',N_stp,' Acc=',N_acc,' Rej=',N_rej, &
85 print *,'Rosenbrock: Unsucessful step at T=', &
86 TIN
,' (IERR=',IERR
,')'
89 STEPMIN
= RCNTRL(ihexit
)
90 ! if optional parameters are given for output they return information
91 IF (PRESENT(ISTATUS_U
)) ISTATUS_U(:) = ISTATUS(1:20)
92 IF (PRESENT(RSTATUS_U
)) RSTATUS_U(:) = RSTATUS(1:20)
94 END SUBROUTINE INTEGRATE_TLM
97 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
98 SUBROUTINE RosenbrockTLM(Y
,NTLM
,Y_tlm
,&
101 RCNTRL
,ICNTRL
,RSTATUS
,ISTATUS
,IERR
)
102 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
104 ! TLM = Tangent Linear Model of a Rosenbrock Method
106 ! Solves the system y'=F(t,y) using a Rosenbrock method defined by:
108 ! G = 1/(H*gamma(1)) - Jac(t0,Y0)
109 ! T_i = t0 + Alpha(i)*H
110 ! Y_i = Y0 + \sum_{j=1}^{i-1} A(i,j)*K_j
111 ! G * K_i = Fun( T_i, Y_i ) + \sum_{j=1}^S C(i,j)/H * K_j +
112 ! gamma(i)*dF/dT(t0, Y0)
113 ! Y1 = Y0 + \sum_{j=1}^S M(j)*K_j
115 ! For details on Rosenbrock methods and their implementation consult:
116 ! E. Hairer and G. Wanner
117 ! "Solving ODEs II. Stiff and differential-algebraic problems".
118 ! Springer series in computational mathematics, Springer-Verlag, 1996.
119 ! The codes contained in the book inspired this implementation.
121 ! (C) Adrian Sandu, August 2004
122 ! Virginia Polytechnic Institute and State University
123 ! Contact: sandu@cs.vt.edu
124 ! This implementation is part of KPP - the Kinetic PreProcessor
125 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
127 !~~~> INPUT ARGUMENTS:
129 !- Y(NVAR) -> vector of initial conditions (at T=Tstart)
130 ! NTLM -> dimension of linearized system,
131 ! i.e. the number of sensitivity coefficients
132 !- Y_tlm(NVAR*NTLM) -> vector of initial sensitivity conditions (at T=Tstart)
133 !- [Tstart,Tend] -> time range of integration
134 ! (if Tstart>Tend the integration is performed backwards in time)
135 !- RelTol, AbsTol -> user precribed accuracy
136 !- SUBROUTINE Fun( T, Y, Ydot ) -> ODE function,
137 ! returns Ydot = Y' = F(T,Y)
138 !- SUBROUTINE Jac( T, Y, Jcb ) -> Jacobian of the ODE function,
139 ! returns Jcb = dF/dY
140 !- ICNTRL(1:10) -> integer inputs parameters
141 !- RCNTRL(1:10) -> real inputs parameters
142 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
144 !~~~> OUTPUT ARGUMENTS:
146 !- Y(NVAR) -> vector of final states (at T->Tend)
147 !- Y_tlm(NVAR*NTLM) -> vector of final sensitivities (at T=Tend)
148 !- ICNTRL(11:20) -> integer output parameters
149 !- RCNTRL(11:20) -> real output parameters
150 !- IERR -> job status upon return
151 ! - succes (positive value) or failure (negative value) -
153 ! = -1 : Improper value for maximal no of steps
154 ! = -2 : Selected Rosenbrock method not implemented
155 ! = -3 : Hmin/Hmax/Hstart must be positive
156 ! = -4 : FacMin/FacMax/FacRej must be positive
157 ! = -5 : Improper tolerance values
158 ! = -6 : No of steps exceeds maximum bound
159 ! = -7 : Step size too small
160 ! = -8 : Matrix is repeatedly singular
161 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
163 !~~~> INPUT PARAMETERS:
165 ! Note: For input parameters equal to zero the default values of the
166 ! corresponding variables are used.
168 ! ICNTRL(1) = 1: F = F(y) Independent of T (AUTONOMOUS)
169 ! = 0: F = F(t,y) Depends on T (NON-AUTONOMOUS)
170 ! ICNTRL(2) = 0: AbsTol, RelTol are NVAR-dimensional vectors
171 ! = 1: AbsTol, RelTol are scalars
172 ! ICNTRL(3) -> maximum number of integration steps
173 ! For ICNTRL(3)=0) the default value of 100000 is used
175 ! ICNTRL(4) -> selection of a particular Rosenbrock method
176 ! = 0 : default method is Rodas3
177 ! = 1 : method is Ros2
178 ! = 2 : method is Ros3
179 ! = 3 : method is Ros4
180 ! = 4 : method is Rodas3
181 ! = 5: method is Rodas4
183 ! RCNTRL(1) -> Hmin, lower bound for the integration step size
184 ! It is strongly recommended to keep Hmin = ZERO
185 ! RCNTRL(2) -> Hmax, upper bound for the integration step size
186 ! RCNTRL(3) -> Hstart, starting value for the integration step size
188 ! RCNTRL(4) -> FacMin, lower bound on step decrease factor (default=0.2)
189 ! RCNTRL(5) -> FacMin,upper bound on step increase factor (default=6)
190 ! RCNTRL(6) -> FacRej, step decrease factor after multiple rejections
192 ! RCNTRL(7) -> FacSafe, by which the new step is slightly smaller
193 ! than the predicted value (default=0.9)
194 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
196 !~~~> OUTPUT PARAMETERS:
198 ! Note: each call to Rosenbrock adds the corrent no. of fcn calls
199 ! to previous value of ISTATUS(1), and similar for the other params.
200 ! Set ISTATUS(1:10) = 0 before call to avoid this accumulation.
202 ! ISTATUS(1) = No. of function calls
203 ! ISTATUS(2) = No. of jacobian calls
204 ! ISTATUS(3) = No. of steps
205 ! ISTATUS(4) = No. of accepted steps
206 ! ISTATUS(5) = No. of rejected steps (except at the beginning)
207 ! ISTATUS(6) = No. of LU decompositions
208 ! ISTATUS(7) = No. of forward/backward substitutions
209 ! ISTATUS(8) = No. of singular matrix decompositions
211 ! RSTATUS(1) -> Texit, the time corresponding to the
212 ! computed Y upon return
213 ! RSTATUS(2) -> Hexit, last accepted step before exit
214 ! For multiple restarts, use Hexit as Hstart in the following run
215 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
220 KPP_REAL
, INTENT(INOUT
) :: Y(NVAR
)
221 INTEGER, INTENT(IN
) :: NTLM
222 KPP_REAL
, INTENT(INOUT
) :: Y_tlm(NVAR
,NTLM
)
223 KPP_REAL
, INTENT(IN
) :: Tstart
, Tend
224 KPP_REAL
, INTENT(IN
) :: AbsTol(NVAR
),RelTol(NVAR
)
225 INTEGER, INTENT(IN
) :: ICNTRL(10)
226 KPP_REAL
, INTENT(IN
) :: RCNTRL(10)
227 INTEGER, INTENT(INOUT
) :: ISTATUS(10)
228 KPP_REAL
, INTENT(INOUT
) :: RSTATUS(10)
229 INTEGER, INTENT(OUT
) :: IERR
230 !~~~> The method parameters
231 INTEGER, PARAMETER :: Smax
= 6
232 INTEGER :: Method
, ros_S
233 KPP_REAL
, DIMENSION(Smax
) :: ros_M
, ros_E
, ros_Alpha
, ros_Gamma
234 KPP_REAL
, DIMENSION(Smax
*(Smax
-1)/2) :: ros_A
, ros_C
236 LOGICAL, DIMENSION(Smax
) :: ros_NewF
237 CHARACTER(LEN
=12) :: ros_Name
238 !~~~> Local variables
239 KPP_REAL
:: Roundoff
, FacMin
, FacMax
, FacRej
, FacSafe
240 KPP_REAL
:: Hmin
, Hmax
, Hstart
, Hexit
242 INTEGER :: i
, UplimTol
, Max_no_steps
243 LOGICAL :: Autonomous
, VectorTol
245 KPP_REAL
, PARAMETER :: ZERO
= 0.0d0, ONE
= 1.0d0
246 KPP_REAL
, PARAMETER :: DeltaMin
= 1.0d-5
248 !~~~> Initialize statistics
258 !~~~> Autonomous or time dependent ODE. Default is time dependent.
259 Autonomous
= .NOT
.(ICNTRL(1) == 0)
261 !~~~> For Scalar tolerances (ICNTRL(2).NE.0) the code uses AbsTol(1) and RelTol(1)
262 ! For Vector tolerances (ICNTRL(2) == 0) the code uses AbsTol(1:NVAR) and RelTol(1:NVAR)
263 IF (ICNTRL(2) == 0) THEN
271 !~~~> The maximum number of steps admitted
272 IF (ICNTRL(3) == 0) THEN
273 Max_no_steps
= 100000
274 ELSEIF (Max_no_steps
> 0) THEN
275 Max_no_steps
=ICNTRL(3)
277 PRINT * ,'User-selected max no. of steps: ICNTRL(3)=',ICNTRL(3)
278 CALL ros_ErrorMsg(-1,Tstart
,ZERO
,IERR
)
282 !~~~> The particular Rosenbrock method chosen
283 IF (ICNTRL(4) == 0) THEN
285 ELSEIF ( (ICNTRL(4) >= 1).AND
.(ICNTRL(4) <= 5) ) THEN
288 PRINT * , 'User-selected Rosenbrock method: ICNTRL(4)=', Method
289 CALL ros_ErrorMsg(-2,Tstart
,ZERO
,IERR
)
293 !~~~> Unit roundoff (1+Roundoff>1)
294 Roundoff
= WLAMCH('E')
296 !~~~> Lower bound on the step size: (positive value)
297 IF (RCNTRL(1) == ZERO
) THEN
299 ELSEIF (RCNTRL(1) > ZERO
) THEN
302 PRINT * , 'User-selected Hmin: RCNTRL(1)=', RCNTRL(1)
303 CALL ros_ErrorMsg(-3,Tstart
,ZERO
,IERR
)
306 !~~~> Upper bound on the step size: (positive value)
307 IF (RCNTRL(2) == ZERO
) THEN
308 Hmax
= ABS(Tend
-Tstart
)
309 ELSEIF (RCNTRL(2) > ZERO
) THEN
310 Hmax
= MIN(ABS(RCNTRL(2)),ABS(Tend
-Tstart
))
312 PRINT * , 'User-selected Hmax: RCNTRL(2)=', RCNTRL(2)
313 CALL ros_ErrorMsg(-3,Tstart
,ZERO
,IERR
)
316 !~~~> Starting step size: (positive value)
317 IF (RCNTRL(3) == ZERO
) THEN
318 Hstart
= MAX(Hmin
,DeltaMin
)
319 ELSEIF (RCNTRL(3) > ZERO
) THEN
320 Hstart
= MIN(ABS(RCNTRL(3)),ABS(Tend
-Tstart
))
322 PRINT * , 'User-selected Hstart: RCNTRL(3)=', RCNTRL(3)
323 CALL ros_ErrorMsg(-3,Tstart
,ZERO
,IERR
)
326 !~~~> Step size can be changed s.t. FacMin < Hnew/Hexit < FacMax
327 IF (RCNTRL(4) == ZERO
) THEN
329 ELSEIF (RCNTRL(4) > ZERO
) THEN
332 PRINT * , 'User-selected FacMin: RCNTRL(4)=', RCNTRL(4)
333 CALL ros_ErrorMsg(-4,Tstart
,ZERO
,IERR
)
336 IF (RCNTRL(5) == ZERO
) THEN
338 ELSEIF (RCNTRL(5) > ZERO
) THEN
341 PRINT * , 'User-selected FacMax: RCNTRL(5)=', RCNTRL(5)
342 CALL ros_ErrorMsg(-4,Tstart
,ZERO
,IERR
)
345 !~~~> FacRej: Factor to decrease step after 2 succesive rejections
346 IF (RCNTRL(6) == ZERO
) THEN
348 ELSEIF (RCNTRL(6) > ZERO
) THEN
351 PRINT * , 'User-selected FacRej: RCNTRL(6)=', RCNTRL(6)
352 CALL ros_ErrorMsg(-4,Tstart
,ZERO
,IERR
)
355 !~~~> FacSafe: Safety Factor in the computation of new step size
356 IF (RCNTRL(7) == ZERO
) THEN
358 ELSEIF (RCNTRL(7) > ZERO
) THEN
361 PRINT * , 'User-selected FacSafe: RCNTRL(7)=', RCNTRL(7)
362 CALL ros_ErrorMsg(-4,Tstart
,ZERO
,IERR
)
365 !~~~> Check if tolerances are reasonable
367 IF ( (AbsTol(i
) <= ZERO
) .OR
. (RelTol(i
) <= 10.d0
*Roundoff
) &
368 .OR
. (RelTol(i
) >= 1.0d0) ) THEN
369 PRINT * , ' AbsTol(',i
,') = ',AbsTol(i
)
370 PRINT * , ' RelTol(',i
,') = ',RelTol(i
)
371 CALL ros_ErrorMsg(-5,Tstart
,ZERO
,IERR
)
377 !~~~> Initialize the particular Rosenbrock method
380 CALL Ros2(ros_S
, ros_A
, ros_C
, ros_M
, ros_E
, &
381 ros_Alpha
, ros_Gamma
, ros_NewF
, ros_ELO
, ros_Name
)
383 CALL Ros3(ros_S
, ros_A
, ros_C
, ros_M
, ros_E
, &
384 ros_Alpha
, ros_Gamma
, ros_NewF
, ros_ELO
, ros_Name
)
386 CALL Ros4(ros_S
, ros_A
, ros_C
, ros_M
, ros_E
, &
387 ros_Alpha
, ros_Gamma
, ros_NewF
, ros_ELO
, ros_Name
)
389 CALL Rodas3(ros_S
, ros_A
, ros_C
, ros_M
, ros_E
, &
390 ros_Alpha
, ros_Gamma
, ros_NewF
, ros_ELO
, ros_Name
)
392 CALL Rodas4(ros_S
, ros_A
, ros_C
, ros_M
, ros_E
, &
393 ros_Alpha
, ros_Gamma
, ros_NewF
, ros_ELO
, ros_Name
)
395 PRINT * , 'Unknown Rosenbrock method: ICNTRL(4)=', Method
396 CALL ros_ErrorMsg(-2,Tstart
,ZERO
,IERR
)
400 !~~~> CALL Rosenbrock method
401 CALL ros_TLM_Int(Y
, NTLM
, Y_tlm
, &
402 Tstart
, Tend
, Texit
, &
404 ! Rosenbrock method coefficients
405 ros_S
, ros_M
, ros_E
, ros_A
, ros_C
, &
406 ros_Alpha
, ros_Gamma
, ros_ELO
, ros_NewF
, &
407 ! Integration parameters
408 Autonomous
, VectorTol
, Max_no_steps
, &
409 Roundoff
, Hmin
, Hmax
, Hstart
, Hexit
, &
410 FacMin
, FacMax
, FacRej
, FacSafe
, &
415 !~~~> Collect run statistics
425 RSTATUS(itexit
) = Texit
426 RSTATUS(ihexit
) = Hexit
429 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
430 END SUBROUTINE RosenbrockTLM
431 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
434 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
435 SUBROUTINE ros_ErrorMsg(Code
,T
,H
,IERR
)
436 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
437 ! Handles all error messages
438 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
440 KPP_REAL
, INTENT(IN
) :: T
, H
441 INTEGER, INTENT(IN
) :: Code
442 INTEGER, INTENT(OUT
) :: IERR
446 'Forced exit from Rosenbrock due to the following error:'
450 PRINT * , '--> Improper value for maximal no of steps'
452 PRINT * , '--> Selected Rosenbrock method not implemented'
454 PRINT * , '--> Hmin/Hmax/Hstart must be positive'
456 PRINT * , '--> FacMin/FacMax/FacRej must be positive'
458 PRINT * , '--> Improper tolerance values'
460 PRINT * , '--> No of steps exceeds maximum bound'
462 PRINT * , '--> Step size too small: T + 10*H = T', &
465 PRINT * , '--> Matrix is repeatedly singular'
467 PRINT *, 'Unknown Error code: ', Code
470 PRINT *, "T=", T
, "and H=", H
472 END SUBROUTINE ros_ErrorMsg
476 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
477 SUBROUTINE ros_TLM_Int (Y
, NTLM
, Y_tlm
, &
480 !~~~> Rosenbrock method coefficients
481 ros_S
, ros_M
, ros_E
, ros_A
, ros_C
, &
482 ros_Alpha
, ros_Gamma
, ros_ELO
, ros_NewF
, &
483 !~~~> Integration parameters
484 Autonomous
, VectorTol
, Max_no_steps
, &
485 Roundoff
, Hmin
, Hmax
, Hstart
, Hexit
, &
486 FacMin
, FacMax
, FacRej
, FacSafe
, &
487 !~~~> Error indicator
489 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
490 ! Template for the implementation of a generic Rosenbrock method
491 ! defined by ros_S (no of stages)
492 ! and its coefficients ros_{A,C,M,E,Alpha,Gamma}
493 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
497 !~~~> Input: the initial condition at Tstart; Output: the solution at T
498 KPP_REAL
, INTENT(INOUT
) :: Y(NVAR
)
499 !~~~> Input: Number of sensitivity coefficients
500 INTEGER, INTENT(IN
) :: NTLM
501 !~~~> Input: the initial sensitivites at Tstart; Output: the sensitivities at T
502 KPP_REAL
, INTENT(INOUT
) :: Y_tlm(NVAR
,NTLM
)
503 !~~~> Input: integration interval
504 KPP_REAL
, INTENT(IN
) :: Tstart
,Tend
505 !~~~> Output: time at which the solution is returned (T=Tend if success)
506 KPP_REAL
, INTENT(OUT
) :: T
507 !~~~> Input: tolerances
508 KPP_REAL
, INTENT(IN
) :: AbsTol(NVAR
), RelTol(NVAR
)
509 !~~~> Input: The Rosenbrock method parameters
510 INTEGER, INTENT(IN
) :: ros_S
511 KPP_REAL
, INTENT(IN
) :: ros_M(ros_S
), ros_E(ros_S
), &
512 ros_Alpha(ros_S
), ros_A(ros_S
*(ros_S
-1)/2), &
513 ros_Gamma(ros_S
), ros_C(ros_S
*(ros_S
-1)/2), ros_ELO
514 LOGICAL, INTENT(IN
) :: ros_NewF(ros_S
)
515 !~~~> Input: integration parameters
516 LOGICAL, INTENT(IN
) :: Autonomous
, VectorTol
517 KPP_REAL
, INTENT(IN
) :: Hstart
, Hmin
, Hmax
518 INTEGER, INTENT(IN
) :: Max_no_steps
519 KPP_REAL
, INTENT(IN
) :: Roundoff
, FacMin
, FacMax
, FacRej
, FacSafe
520 !~~~> Output: last accepted step
521 KPP_REAL
, INTENT(OUT
) :: Hexit
522 !~~~> Output: Error indicator
523 INTEGER, INTENT(OUT
) :: IERR
524 ! ~~~~ Local variables
525 KPP_REAL
:: Ynew(NVAR
), Fcn0(NVAR
), Fcn(NVAR
)
526 KPP_REAL
:: K(NVAR
*ros_S
)
527 KPP_REAL
:: Ynew_tlm(NVAR
,NTLM
), Fcn0_tlm(NVAR
,NTLM
), Fcn_tlm(NVAR
,NTLM
)
528 KPP_REAL
:: K_tlm(NVAR
*ros_S
,NTLM
)
529 KPP_REAL
:: Hes0(NHESS
)
530 KPP_REAL
:: dFdT(NVAR
), dJdT(LU_NONZERO
)
531 KPP_REAL
:: Jac0(LU_NONZERO
), Jac(LU_NONZERO
), Ghimj(LU_NONZERO
)
532 KPP_REAL
:: H
, Hnew
, HC
, HG
, Fac
, Tau
533 KPP_REAL
:: Err
, Yerr(NVAR
)
534 INTEGER :: Pivot(NVAR
), Direction
, ioffset
, i
, j
, istage
, itlm
535 LOGICAL :: RejectLastH
, RejectMoreH
, Singular
536 !~~~> Local parameters
537 KPP_REAL
, PARAMETER :: ZERO
= 0.0d0, ONE
= 1.0d0
538 KPP_REAL
, PARAMETER :: DeltaMin
= 1.0d-5
539 !~~~> Locally called functions
542 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
545 !~~~> Initial preparations
549 IF (ABS(H
) <= 10.D0
*Roundoff
) H
= DeltaMin
551 IF (Tend
>= Tstart
) THEN
560 !~~~> Time loop begins below
562 TimeLoop
: DO WHILE ( (Direction
> 0).AND
.((T
-Tend
)+Roundoff
<= ZERO
) &
563 .OR
. (Direction
< 0).AND
.((Tend
-T
)+Roundoff
<= ZERO
) )
565 IF ( Nstp
> Max_no_steps
) THEN ! Too many steps
566 CALL ros_ErrorMsg(-6,T
,H
,IERR
)
569 IF ( ((T
+0.1d0*H
) == T
).OR
.(H
<= Roundoff
) ) THEN ! Step size too small
570 CALL ros_ErrorMsg(-7,T
,H
,IERR
)
574 !~~~> Limit H if necessary to avoid going beyond Tend
576 H
= MIN(H
,ABS(Tend
-T
))
578 !~~~> Compute the function at current time
579 CALL FunTemplate(T
,Y
,Fcn0
)
581 !~~~> Compute the Jacobian at current time
582 CALL JacTemplate(T
,Y
,Jac0
)
584 !~~~> Compute the Hessian at current time
585 CALL HessTemplate(T
,Y
,Hes0
)
587 !~~~> Compute the TLM function at current time
589 CALL Jac_SP_Vec ( Jac0
, Y_tlm(1,itlm
), Fcn0_tlm(1,itlm
) )
592 !~~~> Compute the function and Jacobian derivatives with respect to T
593 IF (.NOT
.Autonomous
) THEN
594 CALL ros_FunTimeDerivative ( T
, Roundoff
, Y
, Fcn0
, dFdT
)
595 CALL ros_JacTimeDerivative ( T
, Roundoff
, Y
, Jac0
, dJdT
)
598 !~~~> Repeat step calculation until current step accepted
601 CALL ros_PrepareMatrix(H
,Direction
,ros_Gamma(1),&
602 Jac0
,Ghimj
,Pivot
,Singular
)
603 IF (Singular
) THEN ! More than 5 consecutive failed decompositions
604 CALL ros_ErrorMsg(-8,T
,H
,IERR
)
608 !~~~> Compute the stages
609 Stage
: DO istage
= 1, ros_S
611 ! Current istage offset. Current istage vector is K(ioffset+1:ioffset+NVAR)
612 ioffset
= NVAR
*(istage
-1)
614 ! For the 1st istage the function has been computed previously
615 IF ( istage
== 1 ) THEN
616 CALL WCOPY(NVAR
,Fcn0
,1,Fcn
,1)
617 CALL WCOPY(NVAR
*NTLM
,Fcn0_tlm
,1,Fcn_tlm
,1)
618 ! istage>1 and a new function evaluation is needed at the current istage
619 ELSEIF ( ros_NewF(istage
) ) THEN
620 CALL WCOPY(NVAR
,Y
,1,Ynew
,1)
621 CALL WCOPY(NVAR
*NTLM
,Y_tlm
,1,Ynew_tlm
,1)
623 CALL WAXPY(NVAR
,ros_A((istage
-1)*(istage
-2)/2+j
), &
624 K(NVAR
*(j
-1)+1),1,Ynew
,1)
626 CALL WAXPY(NVAR
,ros_A((istage
-1)*(istage
-2)/2+j
), &
627 K_tlm(NVAR
*(j
-1)+1,itlm
),1,Ynew_tlm(1,itlm
),1)
630 Tau
= T
+ ros_Alpha(istage
)*Direction
*H
631 CALL FunTemplate(Tau
,Ynew
,Fcn
)
632 CALL JacTemplate(Tau
,Ynew
,Jac
)
634 CALL Jac_SP_Vec ( Jac
, Ynew_tlm(1,itlm
), Fcn_tlm(1,itlm
) )
636 END IF ! if istage == 1 elseif ros_NewF(istage)
637 CALL WCOPY(NVAR
,Fcn
,1,K(ioffset
+1),1)
639 CALL WCOPY(NVAR
,Fcn_tlm(1,itlm
),1,K_tlm(ioffset
+1,itlm
),1)
642 HC
= ros_C((istage
-1)*(istage
-2)/2+j
)/(Direction
*H
)
643 CALL WAXPY(NVAR
,HC
,K(NVAR
*(j
-1)+1),1,K(ioffset
+1),1)
645 CALL WAXPY(NVAR
,HC
,K_tlm(NVAR
*(j
-1)+1,itlm
),1,K_tlm(ioffset
+1,itlm
),1)
648 IF ((.NOT
. Autonomous
).AND
.(ros_Gamma(istage
).NE
.ZERO
)) THEN
649 HG
= Direction
*H
*ros_Gamma(istage
)
650 CALL WAXPY(NVAR
,HG
,dFdT
,1,K(ioffset
+1),1)
652 CALL Jac_SP_Vec ( dJdT
, Ynew_tlm(1,itlm
), Fcn_tlm(1,itlm
) )
653 CALL WAXPY(NVAR
,HG
,Fcn_tlm(1,itlm
),1,K_tlm(ioffset
+1,itlm
),1)
656 CALL ros_Solve(Ghimj
, Pivot
, K(ioffset
+1))
658 CALL Hess_Vec ( Hes0
, K(ioffset
+1), Y_tlm(1,itlm
), Fcn_tlm(1,itlm
) )
659 CALL WAXPY(NVAR
,ONE
,Fcn_tlm(1,itlm
),1,K_tlm(ioffset
+1,itlm
),1)
660 CALL ros_Solve(Ghimj
, Pivot
, K_tlm(ioffset
+1,itlm
))
666 !~~~> Compute the new solution
667 CALL WCOPY(NVAR
,Y
,1,Ynew
,1)
669 CALL WAXPY(NVAR
,ros_M(j
),K(NVAR
*(j
-1)+1),1,Ynew
,1)
672 CALL WCOPY(NVAR
,Y_tlm(1,itlm
),1,Ynew_tlm(1,itlm
),1)
674 CALL WAXPY(NVAR
,ros_M(j
),K_tlm(NVAR
*(j
-1)+1,itlm
),1,Ynew_tlm(1,itlm
),1)
678 !~~~> Compute the error estimation
679 CALL WSCAL(NVAR
,ZERO
,Yerr
,1)
681 CALL WAXPY(NVAR
,ros_E(j
),K(NVAR
*(j
-1)+1),1,Yerr
,1)
683 Err
= ros_ErrorNorm ( Y
, Ynew
, Yerr
, AbsTol
, RelTol
, VectorTol
)
685 !~~~> New step size is bounded by FacMin <= Hnew/H <= FacMax
686 Fac
= MIN(FacMax
,MAX(FacMin
,FacSafe
/Err
**(ONE
/ros_ELO
)))
689 !~~~> Check the error magnitude and adjust step size
691 IF ( (Err
<= ONE
).OR
.(H
<= Hmin
) ) THEN !~~~> Accept step
693 CALL WCOPY(NVAR
,Ynew
,1,Y
,1)
694 CALL WCOPY(NVAR
*NTLM
,Ynew_tlm
,1,Y_tlm
,1)
696 Hnew
= MAX(Hmin
,MIN(Hnew
,Hmax
))
697 IF (RejectLastH
) THEN ! No step size increase after a rejected step
700 RejectLastH
= .FALSE
.
701 RejectMoreH
= .FALSE
.
703 EXIT UntilAccepted
! EXIT THE LOOP: WHILE STEP NOT ACCEPTED
704 ELSE !~~~> Reject step
705 IF (RejectMoreH
) THEN
708 RejectMoreH
= RejectLastH
721 IERR
= 1 !~~~> The integration was successful
723 END SUBROUTINE ros_TLM_Int
726 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
727 KPP_REAL
FUNCTION ros_ErrorNorm ( Y
, Ynew
, Yerr
, &
728 AbsTol
, RelTol
, VectorTol
)
729 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
730 !~~~> Computes the "scaled norm" of the error vector Yerr
731 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
735 KPP_REAL
, INTENT(IN
) :: Y(NVAR
), Ynew(NVAR
), &
736 Yerr(NVAR
), AbsTol(NVAR
), RelTol(NVAR
)
737 LOGICAL, INTENT(IN
) :: VectorTol
739 KPP_REAL
:: Err
, Scale
, Ymax
741 KPP_REAL
, PARAMETER :: ZERO
= 0.0d0
745 Ymax
= MAX(ABS(Y(i
)),ABS(Ynew(i
)))
747 Scale
= AbsTol(i
)+RelTol(i
)*Ymax
749 Scale
= AbsTol(1)+RelTol(1)*Ymax
751 Err
= Err
+(Yerr(i
)/Scale
)**2
757 END FUNCTION ros_ErrorNorm
760 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
761 SUBROUTINE ros_FunTimeDerivative ( T
, Roundoff
, Y
, &
763 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
764 !~~~> The time partial derivative of the function by finite differences
765 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
768 !~~~> Input arguments
769 KPP_REAL
, INTENT(IN
) :: T
, Roundoff
, Y(NVAR
), Fcn0(NVAR
)
770 !~~~> Output arguments
771 KPP_REAL
, INTENT(OUT
) :: dFdT(NVAR
)
772 !~~~> Local variables
774 KPP_REAL
, PARAMETER :: ONE
= 1.0d0, DeltaMin
= 1.0d-6
776 Delta
= SQRT(Roundoff
)*MAX(DeltaMin
,ABS(T
))
777 CALL FunTemplate(T
+Delta
,Y
,dFdT
)
778 CALL WAXPY(NVAR
,(-ONE
),Fcn0
,1,dFdT
,1)
779 CALL WSCAL(NVAR
,(ONE
/Delta
),dFdT
,1)
781 END SUBROUTINE ros_FunTimeDerivative
784 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
785 SUBROUTINE ros_JacTimeDerivative ( T
, Roundoff
, Y
, &
787 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
788 !~~~> The time partial derivative of the Jacobian by finite differences
789 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
792 !~~~> Input arguments
793 KPP_REAL
, INTENT(IN
) :: T
, Roundoff
, Y(NVAR
), Jac0(LU_NONZERO
)
794 !~~~> Output arguments
795 KPP_REAL
, INTENT(OUT
) :: dJdT(LU_NONZERO
)
796 !~~~> Local variables
798 KPP_REAL
, PARAMETER :: ONE
= 1.0d0, DeltaMin
= 1.0d-6
800 Delta
= SQRT(Roundoff
)*MAX(DeltaMin
,ABS(T
))
801 CALL JacTemplate(T
+Delta
,Y
,dJdT
)
802 CALL WAXPY(LU_NONZERO
,(-ONE
),Jac0
,1,dJdT
,1)
803 CALL WSCAL(LU_NONZERO
,(ONE
/Delta
),dJdT
,1)
805 END SUBROUTINE ros_JacTimeDerivative
808 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
809 SUBROUTINE ros_PrepareMatrix ( H
, Direction
, gam
, &
810 Jac0
, Ghimj
, Pivot
, Singular
)
811 ! --- --- --- --- --- --- --- --- --- --- --- --- ---
812 ! Prepares the LHS matrix for stage calculations
813 ! 1. Construct Ghimj = 1/(H*ham) - Jac0
814 ! "(Gamma H) Inverse Minus Jacobian"
815 ! 2. Repeat LU decomposition of Ghimj until successful.
816 ! -half the step size if LU decomposition fails and retry
817 ! -exit after 5 consecutive fails
818 ! --- --- --- --- --- --- --- --- --- --- --- --- ---
821 !~~~> Input arguments
822 KPP_REAL
, INTENT(IN
) :: gam
, Jac0(LU_NONZERO
)
823 INTEGER, INTENT(IN
) :: Direction
824 !~~~> Output arguments
825 KPP_REAL
, INTENT(OUT
) :: Ghimj(LU_NONZERO
)
826 LOGICAL, INTENT(OUT
) :: Singular
827 INTEGER, INTENT(OUT
) :: Pivot(NVAR
)
828 !~~~> Inout arguments
829 KPP_REAL
, INTENT(INOUT
) :: H
! step size is decreased when LU fails
830 !~~~> Local variables
831 INTEGER :: i
, ising
, Nconsecutive
833 KPP_REAL
, PARAMETER :: ONE
= 1.0d0, HALF
= 0.5d0
840 !~~~> Construct Ghimj = 1/(H*ham) - Jac0
841 CALL WCOPY(LU_NONZERO
,Jac0
,1,Ghimj
,1)
842 CALL WSCAL(LU_NONZERO
,(-ONE
),Ghimj
,1)
843 ghinv
= ONE
/(Direction
*H
*gam
)
845 Ghimj(LU_DIAG(i
)) = Ghimj(LU_DIAG(i
))+ghinv
847 !~~~> Compute LU decomposition
848 CALL ros_Decomp( Ghimj
, Pivot
, ising
)
850 !~~~> If successful done
853 !~~~> If unsuccessful half the step size; if 5 consecutive fails then return
855 Nconsecutive
= Nconsecutive
+1
857 PRINT*,'Warning: LU Decomposition returned ising = ',ising
858 IF (Nconsecutive
<= 5) THEN ! Less than 5 consecutive failed decompositions
860 ELSE ! More than 5 consecutive failed decompositions
862 END IF ! Nconsecutive
865 END DO ! WHILE Singular
867 END SUBROUTINE ros_PrepareMatrix
870 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
871 SUBROUTINE ros_Decomp( A
, Pivot
, ising
)
872 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
873 ! Template for the LU decomposition
874 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
876 !~~~> Inout variables
877 KPP_REAL
, INTENT(INOUT
) :: A(LU_NONZERO
)
878 !~~~> Output variables
879 INTEGER, INTENT(OUT
) :: Pivot(NVAR
), ising
881 CALL KppDecomp ( A
, ising
)
882 !~~~> Note: for a full matrix use Lapack:
883 ! CALL DGETRF( NVAR, NVAR, A, NVAR, Pivot, ising )
888 END SUBROUTINE ros_Decomp
891 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
892 SUBROUTINE ros_Solve( A
, Pivot
, b
)
893 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
894 ! Template for the forward/backward substitution (using pre-computed LU decomposition)
895 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
897 !~~~> Input variables
898 KPP_REAL
, INTENT(IN
) :: A(LU_NONZERO
)
899 INTEGER, INTENT(IN
) :: Pivot(NVAR
)
900 !~~~> InOut variables
901 KPP_REAL
, INTENT(INOUT
) :: b(NVAR
)
903 CALL KppSolve( A
, b
)
904 !~~~> Note: for a full matrix use Lapack:
906 ! CALL DGETRS( 'N', NVAR , NRHS, A, NVAR, Pivot, b, NVAR, INFO )
910 END SUBROUTINE ros_Solve
914 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
915 SUBROUTINE Ros2 (ros_S
,ros_A
,ros_C
,ros_M
,ros_E
,ros_Alpha
,&
916 ros_Gamma
,ros_NewF
,ros_ELO
,ros_Name
)
917 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
918 ! --- AN L-STABLE METHOD, 2 stages, order 2
919 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
923 INTEGER, PARAMETER :: S
= 2
924 INTEGER, INTENT(OUT
) :: ros_S
925 KPP_REAL
, DIMENSION(S
), INTENT(OUT
) :: ros_M
,ros_E
,ros_Alpha
,ros_Gamma
926 KPP_REAL
, DIMENSION(S
*(S
-1)/2), INTENT(OUT
) :: ros_A
, ros_C
927 KPP_REAL
, INTENT(OUT
) :: ros_ELO
928 LOGICAL, DIMENSION(S
), INTENT(OUT
) :: ros_NewF
929 CHARACTER(LEN
=12), INTENT(OUT
) :: ros_Name
932 g
= 1.0d0 + 1.0d0/SQRT(2.0d0)
934 !~~~> Name of the method
936 !~~~> Number of stages
939 !~~~> The coefficient matrices A and C are strictly lower triangular.
940 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
941 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
942 ! The general mapping formula is:
943 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
944 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
948 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
949 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
952 !~~~> M_i = Coefficients for new step solution
953 ros_M(1)= (3.d0
)/(2.d0
*g
)
954 ros_M(2)= (1.d0
)/(2.d0
*g
)
955 ! E_i = Coefficients for error estimator
956 ros_E(1) = 1.d0
/(2.d0
*g
)
957 ros_E(2) = 1.d0
/(2.d0
*g
)
958 !~~~> ros_ELO = estimator of local order - the minimum between the
959 ! main and the embedded scheme orders plus one
961 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
964 !~~~> Gamma_i = \sum_j gamma_{i,j}
971 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
972 SUBROUTINE Ros3 (ros_S
,ros_A
,ros_C
,ros_M
,ros_E
,ros_Alpha
,&
973 ros_Gamma
,ros_NewF
,ros_ELO
,ros_Name
)
974 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
975 ! --- AN L-STABLE METHOD, 3 stages, order 3, 2 function evaluations
976 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
980 INTEGER, PARAMETER :: S
= 3
981 INTEGER, INTENT(OUT
) :: ros_S
982 KPP_REAL
, DIMENSION(S
), INTENT(OUT
) :: ros_M
,ros_E
,ros_Alpha
,ros_Gamma
983 KPP_REAL
, DIMENSION(S
*(S
-1)/2), INTENT(OUT
) :: ros_A
, ros_C
984 KPP_REAL
, INTENT(OUT
) :: ros_ELO
985 LOGICAL, DIMENSION(S
), INTENT(OUT
) :: ros_NewF
986 CHARACTER(LEN
=12), INTENT(OUT
) :: ros_Name
988 !~~~> Name of the method
990 !~~~> Number of stages
993 !~~~> The coefficient matrices A and C are strictly lower triangular.
994 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
995 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
996 ! The general mapping formula is:
997 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
998 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
1004 ros_C(1) = -0.10156171083877702091975600115545d+01
1005 ros_C(2) = 0.40759956452537699824805835358067d+01
1006 ros_C(3) = 0.92076794298330791242156818474003d+01
1007 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
1008 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
1009 ros_NewF(1) = .TRUE
.
1010 ros_NewF(2) = .TRUE
.
1011 ros_NewF(3) = .FALSE
.
1012 !~~~> M_i = Coefficients for new step solution
1014 ros_M(2) = 0.61697947043828245592553615689730d+01
1015 ros_M(3) = -0.42772256543218573326238373806514d+00
1016 ! E_i = Coefficients for error estimator
1018 ros_E(2) = -0.29079558716805469821718236208017d+01
1019 ros_E(3) = 0.22354069897811569627360909276199d+00
1020 !~~~> ros_ELO = estimator of local order - the minimum between the
1021 ! main and the embedded scheme orders plus 1
1023 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
1024 ros_Alpha(1)= 0.0d+00
1025 ros_Alpha(2)= 0.43586652150845899941601945119356d+00
1026 ros_Alpha(3)= 0.43586652150845899941601945119356d+00
1027 !~~~> Gamma_i = \sum_j gamma_{i,j}
1028 ros_Gamma(1)= 0.43586652150845899941601945119356d+00
1029 ros_Gamma(2)= 0.24291996454816804366592249683314d+00
1030 ros_Gamma(3)= 0.21851380027664058511513169485832d+01
1034 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1037 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1038 SUBROUTINE Ros4 (ros_S
,ros_A
,ros_C
,ros_M
,ros_E
,ros_Alpha
,&
1039 ros_Gamma
,ros_NewF
,ros_ELO
,ros_Name
)
1040 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1041 ! L-STABLE ROSENBROCK METHOD OF ORDER 4, WITH 4 STAGES
1042 ! L-STABLE EMBEDDED ROSENBROCK METHOD OF ORDER 3
1044 ! E. HAIRER AND G. WANNER, SOLVING ORDINARY DIFFERENTIAL
1045 ! EQUATIONS II. STIFF AND DIFFERENTIAL-ALGEBRAIC PROBLEMS.
1046 ! SPRINGER SERIES IN COMPUTATIONAL MATHEMATICS,
1047 ! SPRINGER-VERLAG (1990)
1048 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1052 INTEGER, PARAMETER :: S
= 4
1053 INTEGER, INTENT(OUT
) :: ros_S
1054 KPP_REAL
, DIMENSION(S
), INTENT(OUT
) :: ros_M
,ros_E
,ros_Alpha
,ros_Gamma
1055 KPP_REAL
, DIMENSION(S
*(S
-1)/2), INTENT(OUT
) :: ros_A
, ros_C
1056 KPP_REAL
, INTENT(OUT
) :: ros_ELO
1057 LOGICAL, DIMENSION(S
), INTENT(OUT
) :: ros_NewF
1058 CHARACTER(LEN
=12), INTENT(OUT
) :: ros_Name
1061 !~~~> Name of the method
1063 !~~~> Number of stages
1066 !~~~> The coefficient matrices A and C are strictly lower triangular.
1067 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
1068 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
1069 ! The general mapping formula is:
1070 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
1071 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
1073 ros_A(1) = 0.2000000000000000d+01
1074 ros_A(2) = 0.1867943637803922d+01
1075 ros_A(3) = 0.2344449711399156d+00
1080 ros_C(1) =-0.7137615036412310d+01
1081 ros_C(2) = 0.2580708087951457d+01
1082 ros_C(3) = 0.6515950076447975d+00
1083 ros_C(4) =-0.2137148994382534d+01
1084 ros_C(5) =-0.3214669691237626d+00
1085 ros_C(6) =-0.6949742501781779d+00
1086 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
1087 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
1088 ros_NewF(1) = .TRUE
.
1089 ros_NewF(2) = .TRUE
.
1090 ros_NewF(3) = .TRUE
.
1091 ros_NewF(4) = .FALSE
.
1092 !~~~> M_i = Coefficients for new step solution
1093 ros_M(1) = 0.2255570073418735d+01
1094 ros_M(2) = 0.2870493262186792d+00
1095 ros_M(3) = 0.4353179431840180d+00
1096 ros_M(4) = 0.1093502252409163d+01
1097 !~~~> E_i = Coefficients for error estimator
1098 ros_E(1) =-0.2815431932141155d+00
1099 ros_E(2) =-0.7276199124938920d-01
1100 ros_E(3) =-0.1082196201495311d+00
1101 ros_E(4) =-0.1093502252409163d+01
1102 !~~~> ros_ELO = estimator of local order - the minimum between the
1103 ! main and the embedded scheme orders plus 1
1105 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
1107 ros_Alpha(2) = 0.1145640000000000d+01
1108 ros_Alpha(3) = 0.6552168638155900d+00
1109 ros_Alpha(4) = ros_Alpha(3)
1110 !~~~> Gamma_i = \sum_j gamma_{i,j}
1111 ros_Gamma(1) = 0.5728200000000000d+00
1112 ros_Gamma(2) =-0.1769193891319233d+01
1113 ros_Gamma(3) = 0.7592633437920482d+00
1114 ros_Gamma(4) =-0.1049021087100450d+00
1118 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1119 SUBROUTINE Rodas3 (ros_S
,ros_A
,ros_C
,ros_M
,ros_E
,ros_Alpha
,&
1120 ros_Gamma
,ros_NewF
,ros_ELO
,ros_Name
)
1121 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1122 ! --- A STIFFLY-STABLE METHOD, 4 stages, order 3
1123 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1127 INTEGER, PARAMETER :: S
= 4
1128 INTEGER, INTENT(OUT
) :: ros_S
1129 KPP_REAL
, DIMENSION(S
), INTENT(OUT
) :: ros_M
,ros_E
,ros_Alpha
,ros_Gamma
1130 KPP_REAL
, DIMENSION(S
*(S
-1)/2), INTENT(OUT
) :: ros_A
, ros_C
1131 KPP_REAL
, INTENT(OUT
) :: ros_ELO
1132 LOGICAL, DIMENSION(S
), INTENT(OUT
) :: ros_NewF
1133 CHARACTER(LEN
=12), INTENT(OUT
) :: ros_Name
1136 !~~~> Name of the method
1137 ros_Name
= 'RODAS-3'
1138 !~~~> Number of stages
1141 !~~~> The coefficient matrices A and C are strictly lower triangular.
1142 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
1143 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
1144 ! The general mapping formula is:
1145 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
1146 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
1160 ros_C(6) =-(8.0d+00/3.0d+00)
1162 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
1163 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
1164 ros_NewF(1) = .TRUE
.
1165 ros_NewF(2) = .FALSE
.
1166 ros_NewF(3) = .TRUE
.
1167 ros_NewF(4) = .TRUE
.
1168 !~~~> M_i = Coefficients for new step solution
1173 !~~~> E_i = Coefficients for error estimator
1178 !~~~> ros_ELO = estimator of local order - the minimum between the
1179 ! main and the embedded scheme orders plus 1
1181 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
1182 ros_Alpha(1) = 0.0d+00
1183 ros_Alpha(2) = 0.0d+00
1184 ros_Alpha(3) = 1.0d+00
1185 ros_Alpha(4) = 1.0d+00
1186 !~~~> Gamma_i = \sum_j gamma_{i,j}
1187 ros_Gamma(1) = 0.5d+00
1188 ros_Gamma(2) = 1.5d+00
1189 ros_Gamma(3) = 0.0d+00
1190 ros_Gamma(4) = 0.0d+00
1192 END SUBROUTINE Rodas3
1194 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1195 SUBROUTINE Rodas4 (ros_S
,ros_A
,ros_C
,ros_M
,ros_E
,ros_Alpha
,&
1196 ros_Gamma
,ros_NewF
,ros_ELO
,ros_Name
)
1197 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1198 ! STIFFLY-STABLE ROSENBROCK METHOD OF ORDER 4, WITH 6 STAGES
1200 ! E. HAIRER AND G. WANNER, SOLVING ORDINARY DIFFERENTIAL
1201 ! EQUATIONS II. STIFF AND DIFFERENTIAL-ALGEBRAIC PROBLEMS.
1202 ! SPRINGER SERIES IN COMPUTATIONAL MATHEMATICS,
1203 ! SPRINGER-VERLAG (1996)
1204 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1208 INTEGER, PARAMETER :: S
= 6
1209 INTEGER, INTENT(OUT
) :: ros_S
1210 KPP_REAL
, DIMENSION(S
), INTENT(OUT
) :: ros_M
,ros_E
,ros_Alpha
,ros_Gamma
1211 KPP_REAL
, DIMENSION(S
*(S
-1)/2), INTENT(OUT
) :: ros_A
, ros_C
1212 KPP_REAL
, INTENT(OUT
) :: ros_ELO
1213 LOGICAL, DIMENSION(S
), INTENT(OUT
) :: ros_NewF
1214 CHARACTER(LEN
=12), INTENT(OUT
) :: ros_Name
1217 !~~~> Name of the method
1218 ros_Name
= 'RODAS-4'
1219 !~~~> Number of stages
1222 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
1223 ros_Alpha(1) = 0.000d0
1224 ros_Alpha(2) = 0.386d0
1225 ros_Alpha(3) = 0.210d0
1226 ros_Alpha(4) = 0.630d0
1227 ros_Alpha(5) = 1.000d0
1228 ros_Alpha(6) = 1.000d0
1230 !~~~> Gamma_i = \sum_j gamma_{i,j}
1231 ros_Gamma(1) = 0.2500000000000000d+00
1232 ros_Gamma(2) =-0.1043000000000000d+00
1233 ros_Gamma(3) = 0.1035000000000000d+00
1234 ros_Gamma(4) =-0.3620000000000023d-01
1235 ros_Gamma(5) = 0.0d0
1236 ros_Gamma(6) = 0.0d0
1238 !~~~> The coefficient matrices A and C are strictly lower triangular.
1239 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
1240 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
1241 ! The general mapping formula is: A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
1242 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
1244 ros_A(1) = 0.1544000000000000d+01
1245 ros_A(2) = 0.9466785280815826d+00
1246 ros_A(3) = 0.2557011698983284d+00
1247 ros_A(4) = 0.3314825187068521d+01
1248 ros_A(5) = 0.2896124015972201d+01
1249 ros_A(6) = 0.9986419139977817d+00
1250 ros_A(7) = 0.1221224509226641d+01
1251 ros_A(8) = 0.6019134481288629d+01
1252 ros_A(9) = 0.1253708332932087d+02
1253 ros_A(10) =-0.6878860361058950d+00
1254 ros_A(11) = ros_A(7)
1255 ros_A(12) = ros_A(8)
1256 ros_A(13) = ros_A(9)
1257 ros_A(14) = ros_A(10)
1260 ros_C(1) =-0.5668800000000000d+01
1261 ros_C(2) =-0.2430093356833875d+01
1262 ros_C(3) =-0.2063599157091915d+00
1263 ros_C(4) =-0.1073529058151375d+00
1264 ros_C(5) =-0.9594562251023355d+01
1265 ros_C(6) =-0.2047028614809616d+02
1266 ros_C(7) = 0.7496443313967647d+01
1267 ros_C(8) =-0.1024680431464352d+02
1268 ros_C(9) =-0.3399990352819905d+02
1269 ros_C(10) = 0.1170890893206160d+02
1270 ros_C(11) = 0.8083246795921522d+01
1271 ros_C(12) =-0.7981132988064893d+01
1272 ros_C(13) =-0.3152159432874371d+02
1273 ros_C(14) = 0.1631930543123136d+02
1274 ros_C(15) =-0.6058818238834054d+01
1276 !~~~> M_i = Coefficients for new step solution
1280 ros_M(4) = ros_A(10)
1284 !~~~> E_i = Coefficients for error estimator
1292 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
1293 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
1294 ros_NewF(1) = .TRUE
.
1295 ros_NewF(2) = .TRUE
.
1296 ros_NewF(3) = .TRUE
.
1297 ros_NewF(4) = .TRUE
.
1298 ros_NewF(5) = .TRUE
.
1299 ros_NewF(6) = .TRUE
.
1301 !~~~> ros_ELO = estimator of local order - the minimum between the
1302 ! main and the embedded scheme orders plus 1
1305 END SUBROUTINE Rodas4
1309 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1310 SUBROUTINE FunTemplate( T
, Y
, Ydot
)
1311 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1312 ! Template for the ODE function call.
1313 ! Updates the rate coefficients (and possibly the fixed species) at each call
1314 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1317 !~~~> Input variables
1319 !~~~> Output variables
1321 !~~~> Local variables
1327 CALL Update_RCONST()
1328 CALL Fun( Y
, FIX
, RCONST
, Ydot
)
1333 END SUBROUTINE FunTemplate
1336 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1337 SUBROUTINE JacTemplate( T
, Y
, Jcb
)
1338 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1339 ! Template for the ODE Jacobian call.
1340 ! Updates the rate coefficients (and possibly the fixed species) at each call
1341 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1344 !~~~> Input variables
1346 !~~~> Output variables
1347 KPP_REAL
Jcb(LU_NONZERO
)
1348 !~~~> Local variables
1354 CALL Update_RCONST()
1355 CALL Jac_SP( Y
, FIX
, RCONST
, Jcb
)
1360 END SUBROUTINE JacTemplate
1364 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1365 SUBROUTINE HessTemplate( T
, Y
, Hes
)
1366 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1367 ! Template for the ODE Hessian call.
1368 ! Updates the rate coefficients (and possibly the fixed species) at each call
1369 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1372 !~~~> Input variables
1374 !~~~> Output variables
1376 !~~~> Local variables
1382 CALL Update_RCONST()
1383 CALL Hessian( Y
, FIX
, RCONST
, Hes
)
1386 END SUBROUTINE HessTemplate
1388 END MODULE KPP_ROOT_Integrator