1 MODULE KPP_ROOT_Integrator
3 USE KPP_ROOT_Parameters
5 USE KPP_ROOT_JacobianSP
11 INTEGER, PARAMETER :: ifun
=1, ijac
=2, istp
=3, iacc
=4, &
12 irej
=5, idec
=6, isol
=7, isng
=8, itexit
=1, ihexit
=2
15 ! description of the error numbers IERR
16 CHARACTER(LEN
=50), PARAMETER, DIMENSION(-8:1) :: IERR_NAMES
= (/ &
17 'Matrix is repeatedly singular ', & ! -8
18 'Step size too small ', & ! -7
19 'No of steps exceeds maximum bound ', & ! -6
20 'Improper tolerance values ', & ! -5
21 'FacMin/FacMax/FacRej must be positive ', & ! -4
22 'Hmin/Hmax/Hstart must be positive ', & ! -3
23 'Selected Rosenbrock method not implemented ', & ! -2
24 'Improper value for maximal no of steps ', & ! -1
30 SUBROUTINE KPP_ROOT_INTEGRATE( TIN
, TOUT
, &
31 FIX
, VAR
, RCONST
, ATOL
, RTOL
, IRR_WRK
, &
32 ICNTRL_U
, RCNTRL_U
, ISTATUS_U
, RSTATUS_U
, IERR_U
)
34 USE KPP_ROOT_Parameters
35 !! USE KPP_ROOT_Global
37 KPP_REAL
, INTENT(INOUT
), DIMENSION(NFIX
) :: FIX
38 KPP_REAL
, INTENT(INOUT
), DIMENSION(NVAR
) :: VAR
39 KPP_REAL
, INTENT(INOUT
) :: IRR_WRK(NREACT
)
40 KPP_REAL
, INTENT(IN
), DIMENSION(NSPEC
) :: ATOL
, RTOL
41 KPP_REAL
, INTENT(IN
), DIMENSION(NREACT
) :: RCONST
42 KPP_REAL
, INTENT(IN
) :: TIN
! Start Time
43 KPP_REAL
, INTENT(IN
) :: TOUT
! End Time
44 ! Optional input parameters and statistics
45 INTEGER, INTENT(IN
), OPTIONAL
:: ICNTRL_U(20)
46 KPP_REAL
, INTENT(IN
), OPTIONAL
:: RCNTRL_U(20)
47 INTEGER, INTENT(OUT
), OPTIONAL
:: ISTATUS_U(20)
48 KPP_REAL
, INTENT(OUT
), OPTIONAL
:: RSTATUS_U(20)
49 INTEGER, INTENT(OUT
), OPTIONAL
:: IERR_U
54 INTEGER :: N_stp
, N_acc
, N_rej
, N_sng
55 SAVE N_stp
, N_acc
, N_rej
, N_sng
57 KPP_REAL
:: RCNTRL(20), RSTATUS(20)
58 INTEGER :: ICNTRL(20), ISTATUS(20)
66 ! If optional parameters are given, and if they are >0,
67 ! then they overwrite default settings.
68 IF (PRESENT(ICNTRL_U
)) THEN
69 WHERE(ICNTRL_U(:) > 0) ICNTRL(:) = ICNTRL_U(:)
71 IF (PRESENT(RCNTRL_U
)) THEN
72 WHERE(RCNTRL_U(:) > 0) RCNTRL(:) = RCNTRL_U(:)
75 CALL KPP_ROOT_Rosenbrock(VAR
, FIX
, RCONST
, TIN
,TOUT
, &
77 RCNTRL
,ICNTRL
,RSTATUS
,ISTATUS
,IRR_WRK
,IERR
)
79 STEPMIN
= RCNTRL(ihexit
)
80 ! if optional parameters are given for output they to return information
81 IF (PRESENT(ISTATUS_U
)) ISTATUS_U(:) = ISTATUS(:)
82 IF (PRESENT(RSTATUS_U
)) RSTATUS_U(:) = RSTATUS(:)
83 IF (PRESENT(IERR_U
)) IERR_U
= IERR
85 END SUBROUTINE KPP_ROOT_INTEGRATE
87 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
88 SUBROUTINE KPP_ROOT_Rosenbrock(Y
, FIX
, RCONST
, Tstart
,Tend
, &
90 RCNTRL
,ICNTRL
,RSTATUS
,ISTATUS
,IRR_WRK
,IERR
)
91 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
93 ! Solves the system y'=F(t,y) using a Rosenbrock method defined by:
95 ! G = 1/(H*gamma(1)) - Jac(t0,Y0)
96 ! T_i = t0 + Alpha(i)*H
97 ! Y_i = Y0 + \sum_{j=1}^{i-1} A(i,j)*K_j
98 ! G * K_i = Fun( T_i, Y_i ) + \sum_{j=1}^S C(i,j)/H * K_j +
99 ! gamma(i)*dF/dT(t0, Y0)
100 ! Y1 = Y0 + \sum_{j=1}^S M(j)*K_j
102 ! For details on Rosenbrock methods and their implementation consult:
103 ! E. Hairer and G. Wanner
104 ! "Solving ODEs II. Stiff and differential-algebraic problems".
105 ! Springer series in computational mathematics, Springer-Verlag, 1996.
106 ! The codes contained in the book inspired this implementation.
108 ! (C) Adrian Sandu, August 2004
109 ! Virginia Polytechnic Institute and State University
110 ! Contact: sandu@cs.vt.edu
111 ! This implementation is part of KPP - the Kinetic PreProcessor
112 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
114 !~~~> INPUT ARGUMENTS:
116 !- Y(NVAR) = vector of initial conditions (at T=Tstart)
117 !- [Tstart,Tend] = time range of integration
118 ! (if Tstart>Tend the integration is performed backwards in time)
119 !- RelTol, AbsTol = user precribed accuracy
120 !- SUBROUTINE Fun( T, Y, Ydot ) = ODE function,
121 ! returns Ydot = Y' = F(T,Y)
122 !- SUBROUTINE Jac( T, Y, Jcb ) = Jacobian of the ODE function,
123 ! returns Jcb = dFun/dY
124 !- ICNTRL(1:20) = integer inputs parameters
125 !- RCNTRL(1:20) = real inputs parameters
126 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
128 !~~~> OUTPUT ARGUMENTS:
130 !- Y(NVAR) -> vector of final states (at T->Tend)
131 !- ISTATUS(1:20) -> integer output parameters
132 !- RSTATUS(1:20) -> real output parameters
133 !- IERR -> job status upon return
134 ! success (positive value) or
135 ! failure (negative value)
136 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
138 !~~~> INPUT PARAMETERS:
140 ! Note: For input parameters equal to zero the default values of the
141 ! corresponding variables are used.
143 ! ICNTRL(1) = 1: F = F(y) Independent of T (AUTONOMOUS)
144 ! = 0: F = F(t,y) Depends on T (NON-AUTONOMOUS)
146 ! ICNTRL(2) = 0: AbsTol, RelTol are NVAR-dimensional vectors
147 ! = 1: AbsTol, RelTol are scalars
149 ! ICNTRL(3) -> selection of a particular Rosenbrock method
150 ! = 0 : default method is Rodas3
151 ! = 1 : method is Ros2
152 ! = 2 : method is Ros3
153 ! = 3 : method is Ros4
154 ! = 4 : method is Rodas3
155 ! = 5: method is Rodas4
157 ! ICNTRL(4) -> maximum number of integration steps
158 ! For ICNTRL(4)=0) the default value of 100000 is used
160 ! RCNTRL(1) -> Hmin, lower bound for the integration step size
161 ! It is strongly recommended to keep Hmin = ZERO
162 ! RCNTRL(2) -> Hmax, upper bound for the integration step size
163 ! RCNTRL(3) -> Hstart, starting value for the integration step size
165 ! RCNTRL(4) -> FacMin, lower bound on step decrease factor (default=0.2)
166 ! RCNTRL(5) -> FacMax, upper bound on step increase factor (default=6)
167 ! RCNTRL(6) -> FacRej, step decrease factor after multiple rejections
169 ! RCNTRL(7) -> FacSafe, by which the new step is slightly smaller
170 ! than the predicted value (default=0.9)
171 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
173 !~~~> OUTPUT PARAMETERS:
175 ! Note: each call to Rosenbrock adds the current no. of fcn calls
176 ! to previous value of ISTATUS(1), and similar for the other params.
177 ! Set ISTATUS(1:20) = 0 before call to avoid this accumulation.
179 ! ISTATUS(1) = No. of function calls
180 ! ISTATUS(2) = No. of jacobian calls
181 ! ISTATUS(3) = No. of steps
182 ! ISTATUS(4) = No. of accepted steps
183 ! ISTATUS(5) = No. of rejected steps (except at the beginning)
184 ! ISTATUS(6) = No. of LU decompositions
185 ! ISTATUS(7) = No. of forward/backward substitutions
186 ! ISTATUS(8) = No. of singular matrix decompositions
188 ! RSTATUS(1) -> Texit, the time corresponding to the
189 ! computed Y upon return
190 ! RSTATUS(2) -> Hexit, last accepted step before exit
191 ! For multiple restarts, use Hexit as Hstart in the following run
192 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
194 USE KPP_ROOT_Parameters
195 !! USE KPP_ROOT_LinearAlgebra
199 KPP_REAL
, INTENT(INOUT
) :: Y(NVAR
)
200 KPP_REAL
, INTENT(INOUT
) :: IRR_WRK(NREACT
)
201 KPP_REAL
, INTENT(IN
), DIMENSION(NFIX
) :: FIX
202 KPP_REAL
, INTENT(IN
), DIMENSION(NREACT
) :: RCONST
203 KPP_REAL
, INTENT(IN
) :: Tstart
,Tend
204 KPP_REAL
, INTENT(IN
) :: AbsTol(NVAR
),RelTol(NVAR
)
205 INTEGER, INTENT(IN
) :: ICNTRL(20)
206 KPP_REAL
, INTENT(IN
) :: RCNTRL(20)
207 INTEGER, INTENT(INOUT
) :: ISTATUS(20)
208 KPP_REAL
, INTENT(INOUT
) :: RSTATUS(20)
209 INTEGER, INTENT(OUT
) :: IERR
210 !~~~> The method parameters
211 INTEGER, PARAMETER :: Smax
= 6
212 INTEGER :: Method
, ros_S
213 KPP_REAL
, DIMENSION(Smax
) :: ros_M
, ros_E
, ros_Alpha
, ros_Gamma
214 KPP_REAL
, DIMENSION(Smax
*(Smax
-1)/2) :: ros_A
, ros_C
216 LOGICAL, DIMENSION(Smax
) :: ros_NewF
217 CHARACTER(LEN
=12) :: ros_Name
219 !~~~> Statistics on the work performed by the Rosenbrock method
220 INTEGER :: Nfun
,Njac
,Nstp
,Nacc
,Nrej
,Ndec
,Nsol
,Nsng
223 !~~~> Local variables
224 KPP_REAL
:: Roundoff
, FacMin
, FacMax
, FacRej
, FacSafe
225 KPP_REAL
:: Hmin
, Hmax
, Hstart
, Hexit
227 INTEGER :: i
, UplimTol
, Max_no_steps
228 LOGICAL :: Autonomous
, VectorTol
230 KPP_REAL
, PARAMETER :: ZERO
= 0.0_dp
, ONE
= 1.0_dp
231 KPP_REAL
, PARAMETER :: DeltaMin
= 1.0E-5_dp
233 !~~~> Initialize statistics
243 !~~~> Autonomous or time dependent ODE. Default is time dependent.
244 Autonomous
= .NOT
.(ICNTRL(1) == 0)
250 !~~~> For Scalar tolerances (ICNTRL(2).NE.0) the code uses AbsTol(1) and RelTol(1)
251 ! For Vector tolerances (ICNTRL(2) == 0) the code uses AbsTol(1:NVAR) and RelTol(1:NVAR)
252 IF (ICNTRL(2) == 0) THEN
260 !~~~> The particular Rosenbrock method chosen
261 IF (ICNTRL(3) == 0) THEN
263 ELSEIF ( (ICNTRL(3) >= 1).AND
.(ICNTRL(3) <= 5) ) THEN
266 PRINT * , 'User-selected Rosenbrock method: ICNTRL(3)=', Method
267 CALL KPP_ROOT_ros_ErrorMsg(-2,Tstart
,ZERO
,IERR
)
271 !~~~> The maximum number of steps admitted
272 IF (ICNTRL(4) == 0) THEN
273 Max_no_steps
= 100000
274 ELSEIF (ICNTRL(4) > 0) THEN
275 Max_no_steps
=ICNTRL(4)
277 PRINT * ,'User-selected max no. of steps: ICNTRL(4)=',ICNTRL(4)
278 CALL KPP_ROOT_ros_ErrorMsg(-1,Tstart
,ZERO
,IERR
)
282 !~~~> Unit roundoff (1+Roundoff>1)
283 Roundoff
= KPP_ROOT_WLAMCH('E')
285 !~~~> Lower bound on the step size: (positive value)
286 IF (RCNTRL(1) == ZERO
) THEN
288 ELSEIF (RCNTRL(1) > ZERO
) THEN
291 PRINT * , 'User-selected Hmin: RCNTRL(1)=', RCNTRL(1)
292 CALL KPP_ROOT_ros_ErrorMsg(-3,Tstart
,ZERO
,IERR
)
295 !~~~> Upper bound on the step size: (positive value)
296 IF (RCNTRL(2) == ZERO
) THEN
297 Hmax
= ABS(Tend
-Tstart
)
298 ELSEIF (RCNTRL(2) > ZERO
) THEN
299 Hmax
= MIN(ABS(RCNTRL(2)),ABS(Tend
-Tstart
))
301 PRINT * , 'User-selected Hmax: RCNTRL(2)=', RCNTRL(2)
302 CALL KPP_ROOT_ros_ErrorMsg(-3,Tstart
,ZERO
,IERR
)
305 !~~~> Starting step size: (positive value)
306 IF (RCNTRL(3) == ZERO
) THEN
307 Hstart
= MAX(Hmin
,DeltaMin
)
308 ELSEIF (RCNTRL(3) > ZERO
) THEN
309 Hstart
= MIN(ABS(RCNTRL(3)),ABS(Tend
-Tstart
))
311 PRINT * , 'User-selected Hstart: RCNTRL(3)=', RCNTRL(3)
312 CALL KPP_ROOT_ros_ErrorMsg(-3,Tstart
,ZERO
,IERR
)
315 !~~~> Step size can be changed s.t. FacMin < Hnew/Hexit < FacMax
316 IF (RCNTRL(4) == ZERO
) THEN
318 ELSEIF (RCNTRL(4) > ZERO
) THEN
321 PRINT * , 'User-selected FacMin: RCNTRL(4)=', RCNTRL(4)
322 CALL KPP_ROOT_ros_ErrorMsg(-4,Tstart
,ZERO
,IERR
)
325 IF (RCNTRL(5) == ZERO
) THEN
327 ELSEIF (RCNTRL(5) > ZERO
) THEN
330 PRINT * , 'User-selected FacMax: RCNTRL(5)=', RCNTRL(5)
331 CALL KPP_ROOT_ros_ErrorMsg(-4,Tstart
,ZERO
,IERR
)
334 !~~~> FacRej: Factor to decrease step after 2 succesive rejections
335 IF (RCNTRL(6) == ZERO
) THEN
337 ELSEIF (RCNTRL(6) > ZERO
) THEN
340 PRINT * , 'User-selected FacRej: RCNTRL(6)=', RCNTRL(6)
341 CALL KPP_ROOT_ros_ErrorMsg(-4,Tstart
,ZERO
,IERR
)
344 !~~~> FacSafe: Safety Factor in the computation of new step size
345 IF (RCNTRL(7) == ZERO
) THEN
347 ELSEIF (RCNTRL(7) > ZERO
) THEN
350 PRINT * , 'User-selected FacSafe: RCNTRL(7)=', RCNTRL(7)
351 CALL KPP_ROOT_ros_ErrorMsg(-4,Tstart
,ZERO
,IERR
)
354 !~~~> Check if tolerances are reasonable
356 IF ( (AbsTol(i
) <= ZERO
) .OR
. (RelTol(i
) <= 10.0_dp
*Roundoff
) &
357 .OR
. (RelTol(i
) >= 1.0_dp
) ) THEN
358 PRINT * , ' AbsTol(',i
,') = ',AbsTol(i
)
359 PRINT * , ' RelTol(',i
,') = ',RelTol(i
)
360 CALL KPP_ROOT_ros_ErrorMsg(-5,Tstart
,ZERO
,IERR
)
366 !~~~> Initialize the particular Rosenbrock method
369 CALL KPP_ROOT_Ros2(ros_S
, ros_A
, ros_C
, ros_M
, ros_E
, &
370 ros_Alpha
, ros_Gamma
, ros_NewF
, ros_ELO
, ros_Name
)
372 CALL KPP_ROOT_Ros3(ros_S
, ros_A
, ros_C
, ros_M
, ros_E
, &
373 ros_Alpha
, ros_Gamma
, ros_NewF
, ros_ELO
, ros_Name
)
375 CALL KPP_ROOT_Ros4(ros_S
, ros_A
, ros_C
, ros_M
, ros_E
, &
376 ros_Alpha
, ros_Gamma
, ros_NewF
, ros_ELO
, ros_Name
)
378 CALL KPP_ROOT_Rodas3(ros_S
, ros_A
, ros_C
, ros_M
, ros_E
, &
379 ros_Alpha
, ros_Gamma
, ros_NewF
, ros_ELO
, ros_Name
)
381 CALL KPP_ROOT_Rodas4(ros_S
, ros_A
, ros_C
, ros_M
, ros_E
, &
382 ros_Alpha
, ros_Gamma
, ros_NewF
, ros_ELO
, ros_Name
)
384 PRINT * , 'Unknown Rosenbrock method: ICNTRL(4)=', Method
385 CALL KPP_ROOT_ros_ErrorMsg(-2,Tstart
,ZERO
,IERR
)
389 !~~~> CALL Rosenbrock method
390 CALL KPP_ROOT_ros_Integrator(Y
,Tstart
,Tend
,Texit
, &
392 ! Rosenbrock method coefficients
393 ros_S
, ros_M
, ros_E
, ros_A
, ros_C
, &
394 ros_Alpha
, ros_Gamma
, ros_ELO
, ros_NewF
, &
395 ! Integration parameters
396 Autonomous
, VectorTol
, Max_no_steps
, &
397 Roundoff
, Hmin
, Hmax
, Hstart
, Hexit
, &
398 FacMin
, FacMax
, FacRej
, FacSafe
, &
401 ! Statistics on the work performed by the Rosenbrock method
402 Nfun
,Njac
,Nstp
,Nacc
,Nrej
,Ndec
,Nsol
,Nsng
,&
408 !~~~> Collect run statistics
418 RSTATUS(itexit
) = Texit
419 RSTATUS(ihexit
) = Hexit
421 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
422 CONTAINS ! SUBROUTINES internal to Rosenbrock
423 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
425 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
426 SUBROUTINE KPP_ROOT_ros_ErrorMsg(Code
,T
,H
,IERR
)
427 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
428 ! Handles all error messages
429 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
430 USE KPP_ROOT_Precision
432 KPP_REAL
, INTENT(IN
) :: T
, H
433 INTEGER, INTENT(IN
) :: Code
434 INTEGER, INTENT(OUT
) :: IERR
438 'Forced exit from Rosenbrock due to the following error:'
439 IF ((Code
>=-8).AND
.(Code
<=-1)) THEN
440 PRINT *, IERR_NAMES(Code
)
442 PRINT *, 'Unknown Error code: ', Code
445 PRINT *, "T=", T
, "and H=", H
447 END SUBROUTINE KPP_ROOT_ros_ErrorMsg
449 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
450 SUBROUTINE KPP_ROOT_ros_Integrator (Y
, Tstart
, Tend
, T
, &
452 !~~~> Rosenbrock method coefficients
453 ros_S
, ros_M
, ros_E
, ros_A
, ros_C
, &
454 ros_Alpha
, ros_Gamma
, ros_ELO
, ros_NewF
, &
455 !~~~> Integration parameters
456 Autonomous
, VectorTol
, Max_no_steps
, &
457 Roundoff
, Hmin
, Hmax
, Hstart
, Hexit
, &
458 FacMin
, FacMax
, FacRej
, FacSafe
, &
459 !~~~> Error indicator
461 !~~~> Statistics on the work performed by the Rosenbrock method
462 Nfun
,Njac
,Nstp
,Nacc
,Nrej
,Ndec
,Nsol
,Nsng
, &
465 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
466 ! Template for the implementation of a generic Rosenbrock method
467 ! defined by ros_S (no of stages)
468 ! and its coefficients ros_{A,C,M,E,Alpha,Gamma}
469 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
473 !~~~> Input: the initial condition at Tstart; Output: the solution at T
474 KPP_REAL
, INTENT(INOUT
) :: Y(NVAR
)
475 !~~~> Output: the reaction rates
476 KPP_REAL
, INTENT(INOUT
) :: IRR_WRK(NREACT
)
477 !~~~> Input: integration interval
478 KPP_REAL
, INTENT(IN
) :: Tstart
,Tend
479 !~~~> Output: time at which the solution is returned (T=Tend if success)
480 KPP_REAL
, INTENT(OUT
) :: T
481 !~~~> Input: tolerances
482 KPP_REAL
, INTENT(IN
) :: AbsTol(NVAR
), RelTol(NVAR
)
483 !~~~> Input: The Rosenbrock method parameters
484 INTEGER, INTENT(IN
) :: ros_S
485 KPP_REAL
, INTENT(IN
) :: ros_M(ros_S
), ros_E(ros_S
), &
486 ros_Alpha(ros_S
), ros_A(ros_S
*(ros_S
-1)/2), &
487 ros_Gamma(ros_S
), ros_C(ros_S
*(ros_S
-1)/2), ros_ELO
488 LOGICAL, INTENT(IN
) :: ros_NewF(ros_S
)
489 !~~~> Input: integration parameters
490 LOGICAL, INTENT(IN
) :: Autonomous
, VectorTol
491 KPP_REAL
, INTENT(IN
) :: Hstart
, Hmin
, Hmax
492 INTEGER, INTENT(IN
) :: Max_no_steps
493 KPP_REAL
, INTENT(IN
) :: Roundoff
, FacMin
, FacMax
, FacRej
, FacSafe
494 !~~~> Output: last accepted step
495 KPP_REAL
, INTENT(OUT
) :: Hexit
496 !~~~> Output: Error indicator
497 INTEGER, INTENT(OUT
) :: IERR
499 KPP_REAL
, INTENT(IN
), DIMENSION(NFIX
) :: FIX
501 KPP_REAL
, INTENT(IN
), DIMENSION(NREACT
) :: RCONST
503 !~~~> Statistics on the work performed by the Rosenbrock method
504 INTEGER, INTENT(INOUT
) :: Nfun
,Njac
,Nstp
,Nacc
,Nrej
,Ndec
,Nsol
,Nsng
506 ! ~~~~ Local variables
507 KPP_REAL
:: Ynew(NVAR
), Fcn0(NVAR
), Fcn(NVAR
)
508 KPP_REAL
:: K(NVAR
*ros_S
), dFdT(NVAR
)
510 KPP_REAL
:: Jac0(NVAR
,NVAR
), Ghimj(NVAR
,NVAR
)
512 KPP_REAL
:: Jac0(LU_NONZERO
), Ghimj(LU_NONZERO
)
514 KPP_REAL
:: H
, Hnew
, HC
, HG
, Fac
, Tau
515 KPP_REAL
:: Err
, Yerr(NVAR
)
516 INTEGER :: Pivot(NVAR
), Direction
, ioffset
, j
, istage
517 LOGICAL :: RejectLastH
, RejectMoreH
, Singular
518 !~~~> Local parameters
519 KPP_REAL
, PARAMETER :: ZERO
= 0.0_dp
, ONE
= 1.0_dp
520 KPP_REAL
, PARAMETER :: DeltaMin
= 1.0E-5_dp
521 !~~~> Locally called functions
524 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
527 !~~~> Initial preparations
531 IF (ABS(H
) <= 10.0_dp
*Roundoff
) H
= DeltaMin
533 IF (Tend
>= Tstart
) THEN
542 !~~~> Time loop begins below
544 TimeLoop
: DO WHILE ( (Direction
> 0).AND
.((T
-Tend
)+Roundoff
<= ZERO
) &
545 .OR
. (Direction
< 0).AND
.((Tend
-T
)+Roundoff
<= ZERO
) )
547 IF ( Nstp
> Max_no_steps
) THEN ! Too many steps
548 CALL KPP_ROOT_ros_ErrorMsg(-6,T
,H
,IERR
)
551 IF ( ((T
+0.1_dp
*H
) == T
).OR
.(H
<= Roundoff
) ) THEN ! Step size too small
552 CALL KPP_ROOT_ros_ErrorMsg(-7,T
,H
,IERR
)
556 !~~~> Limit H if necessary to avoid going beyond Tend
558 H
= MIN(H
,ABS(Tend
-T
))
560 !~~~> Compute the function at current time
561 CALL KPP_ROOT_FunTemplate(T
,Y
,Fcn0
, RCONST
, FIX
, Nfun
)
562 IF( T
== Tstart
) THEN
563 CALL KPP_ROOT_IRRFun( Y
, FIX
, RCONST
, IRR_WRK
)
566 !~~~> Compute the function derivative with respect to T
567 IF (.NOT
.Autonomous
) THEN
568 CALL KPP_ROOT_ros_FunTimeDeriv ( T
, Roundoff
, Y
, &
569 Fcn0
, dFdT
, RCONST
, FIX
, Nfun
)
572 !~~~> Compute the Jacobian at current time
573 CALL KPP_ROOT_JacTemplate(T
,Y
,Jac0
, FIX
, Njac
, RCONST
)
575 !~~~> Repeat step calculation until current step accepted
578 CALL KPP_ROOT_ros_PrepareMatrix(H
,Direction
,ros_Gamma(1), &
579 Jac0
,Ghimj
,Pivot
,Singular
, Ndec
, Nsng
)
580 IF (Singular
) THEN ! More than 5 consecutive failed decompositions
581 CALL KPP_ROOT_ros_ErrorMsg(-8,T
,H
,IERR
)
585 !~~~> Compute the stages
586 Stage
: DO istage
= 1, ros_S
588 ! Current istage offset. Current istage vector is K(ioffset+1:ioffset+NVAR)
589 ioffset
= NVAR
*(istage
-1)
591 ! For the 1st istage the function has been computed previously
592 IF ( istage
== 1 ) THEN
593 CALL KPP_ROOT_WCOPY(NVAR
,Fcn0
,1,Fcn
,1)
594 ! istage>1 and a new function evaluation is needed at the current istage
595 ELSEIF ( ros_NewF(istage
) ) THEN
596 CALL KPP_ROOT_WCOPY(NVAR
,Y
,1,Ynew
,1)
598 CALL KPP_ROOT_WAXPY(NVAR
,ros_A((istage
-1)*(istage
-2)/2+j
), &
599 K(NVAR
*(j
-1)+1),1,Ynew
,1)
601 Tau
= T
+ ros_Alpha(istage
)*Direction
*H
602 CALL KPP_ROOT_FunTemplate(Tau
,Ynew
,Fcn
, RCONST
, FIX
, Nfun
)
603 END IF ! if istage == 1 elseif ros_NewF(istage)
604 CALL KPP_ROOT_WCOPY(NVAR
,Fcn
,1,K(ioffset
+1),1)
606 HC
= ros_C((istage
-1)*(istage
-2)/2+j
)/(Direction
*H
)
607 CALL KPP_ROOT_WAXPY(NVAR
,HC
,K(NVAR
*(j
-1)+1),1,K(ioffset
+1),1)
609 IF ((.NOT
. Autonomous
).AND
.(ros_Gamma(istage
).NE
.ZERO
)) THEN
610 HG
= Direction
*H
*ros_Gamma(istage
)
611 CALL KPP_ROOT_WAXPY(NVAR
,HG
,dFdT
,1,K(ioffset
+1),1)
613 CALL KPP_ROOT_ros_Solve(Ghimj
, Pivot
, K(ioffset
+1), Nsol
)
618 !~~~> Compute the new solution
619 CALL KPP_ROOT_WCOPY(NVAR
,Y
,1,Ynew
,1)
621 CALL KPP_ROOT_WAXPY(NVAR
,ros_M(j
),K(NVAR
*(j
-1)+1),1,Ynew
,1)
624 !~~~> Compute the error estimation
625 CALL KPP_ROOT_WSCAL(NVAR
,ZERO
,Yerr
,1)
627 CALL KPP_ROOT_WAXPY(NVAR
,ros_E(j
),K(NVAR
*(j
-1)+1),1,Yerr
,1)
629 Err
= KPP_ROOT_ros_ErrorNorm ( Y
, Ynew
, Yerr
, AbsTol
, RelTol
, VectorTol
)
631 !~~~> New step size is bounded by FacMin <= Hnew/H <= FacMax
632 Fac
= MIN(FacMax
,MAX(FacMin
,FacSafe
/Err
**(ONE
/ros_ELO
)))
635 !~~~> Check the error magnitude and adjust step size
637 IF ( (Err
<= ONE
).OR
.(H
<= Hmin
) ) THEN !~~~> Accept step
639 CALL KPP_ROOT_WCOPY(NVAR
,Ynew
,1,Y
,1)
641 Hnew
= MAX(Hmin
,MIN(Hnew
,Hmax
))
642 IF (RejectLastH
) THEN ! No step size increase after a rejected step
645 RejectLastH
= .FALSE
.
646 RejectMoreH
= .FALSE
.
648 EXIT UntilAccepted
! EXIT THE LOOP: WHILE STEP NOT ACCEPTED
649 ELSE !~~~> Reject step
650 IF (RejectMoreH
) THEN
653 RejectMoreH
= RejectLastH
666 IERR
= 1 !~~~> The integration was successful
668 END SUBROUTINE KPP_ROOT_ros_Integrator
671 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
672 KPP_REAL
FUNCTION KPP_ROOT_ros_ErrorNorm ( Y
, Ynew
, Yerr
, &
673 AbsTol
, RelTol
, VectorTol
)
674 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
675 !~~~> Computes the "scaled norm" of the error vector Yerr
676 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
680 KPP_REAL
, INTENT(IN
) :: Y(NVAR
), Ynew(NVAR
), &
681 Yerr(NVAR
), AbsTol(NVAR
), RelTol(NVAR
)
682 LOGICAL, INTENT(IN
) :: VectorTol
684 KPP_REAL
:: Err
, Scale
, Ymax
686 KPP_REAL
, PARAMETER :: ZERO
= 0.0_dp
690 Ymax
= MAX(ABS(Y(i
)),ABS(Ynew(i
)))
692 Scale
= AbsTol(i
)+RelTol(i
)*Ymax
694 Scale
= AbsTol(1)+RelTol(1)*Ymax
696 Err
= Err
+(Yerr(i
)/Scale
)**2
700 KPP_ROOT_ros_ErrorNorm
= Err
702 END FUNCTION KPP_ROOT_ros_ErrorNorm
705 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
706 SUBROUTINE KPP_ROOT_ros_FunTimeDeriv ( T
, Roundoff
, Y
, &
707 Fcn0
, dFdT
, RCONST
, FIX
, Nfun
)
708 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
709 !~~~> The time partial derivative of the function by finite differences
710 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
713 !~~~> Input arguments
714 KPP_REAL
, INTENT(IN
) :: T
, Roundoff
, Y(NVAR
), Fcn0(NVAR
)
715 KPP_REAL
, INTENT(IN
) :: RCONST(NREACT
), FIX(NFIX
)
716 !~~~> Output arguments
717 KPP_REAL
, INTENT(OUT
) :: dFdT(NVAR
)
719 INTEGER, INTENT(INOUT
) ::Nfun
720 !~~~> Local variables
722 KPP_REAL
, PARAMETER :: ONE
= 1.0_dp
, DeltaMin
= 1.0E-6_dp
724 Delta
= SQRT(Roundoff
)*MAX(DeltaMin
,ABS(T
))
725 CALL KPP_ROOT_FunTemplate(T
+Delta
,Y
,dFdT
, RCONST
, FIX
, Nfun
)
726 CALL KPP_ROOT_WAXPY(NVAR
,(-ONE
),Fcn0
,1,dFdT
,1)
727 CALL KPP_ROOT_WSCAL(NVAR
,(ONE
/Delta
),dFdT
,1)
729 END SUBROUTINE KPP_ROOT_ros_FunTimeDeriv
732 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
733 SUBROUTINE KPP_ROOT_ros_PrepareMatrix ( H
, Direction
, gam
, &
734 Jac0
, Ghimj
, Pivot
, Singular
, Ndec
, Nsng
)
735 ! --- --- --- --- --- --- --- --- --- --- --- --- ---
736 ! Prepares the LHS matrix for stage calculations
737 ! 1. Construct Ghimj = 1/(H*ham) - Jac0
738 ! "(Gamma H) Inverse Minus Jacobian"
739 ! 2. Repeat LU decomposition of Ghimj until successful.
740 ! -half the step size if LU decomposition fails and retry
741 ! -exit after 5 consecutive fails
742 ! --- --- --- --- --- --- --- --- --- --- --- --- ---
745 !~~~> Input arguments
747 KPP_REAL
, INTENT(IN
) :: Jac0(NVAR
,NVAR
)
749 KPP_REAL
, INTENT(IN
) :: Jac0(LU_NONZERO
)
751 KPP_REAL
, INTENT(IN
) :: gam
752 INTEGER, INTENT(IN
) :: Direction
753 !~~~> Output arguments
755 KPP_REAL
, INTENT(OUT
) :: Ghimj(NVAR
,NVAR
)
757 KPP_REAL
, INTENT(OUT
) :: Ghimj(LU_NONZERO
)
759 LOGICAL, INTENT(OUT
) :: Singular
760 INTEGER, INTENT(OUT
) :: Pivot(NVAR
)
761 !~~~> Inout arguments
762 KPP_REAL
, INTENT(INOUT
) :: H
! step size is decreased when LU fails
763 INTEGER, INTENT(INOUT
) :: Ndec
, Nsng
764 !~~~> Local variables
765 INTEGER :: i
, ising
, Nconsecutive
767 KPP_REAL
, PARAMETER :: ONE
= 1.0_dp
, HALF
= 0.5_dp
774 !~~~> Construct Ghimj = 1/(H*gam) - Jac0
776 CALL KPP_ROOT_WCOPY(NVAR
*NVAR
,Jac0
,1,Ghimj
,1)
777 CALL KPP_ROOT_WSCAL(NVAR
*NVAR
,(-ONE
),Ghimj
,1)
778 ghinv
= ONE
/(Direction
*H
*gam
)
780 Ghimj(i
,i
) = Ghimj(i
,i
)+ghinv
783 CALL KPP_ROOT_WCOPY(LU_NONZERO
,Jac0
,1,Ghimj
,1)
784 CALL KPP_ROOT_WSCAL(LU_NONZERO
,(-ONE
),Ghimj
,1)
785 ghinv
= ONE
/(Direction
*H
*gam
)
787 Ghimj(LU_DIAG(i
)) = Ghimj(LU_DIAG(i
))+ghinv
790 !~~~> Compute LU decomposition
791 CALL KPP_ROOT_ros_Decomp( Ghimj
, Pivot
, ising
, Ndec
)
793 !~~~> If successful done
796 !~~~> If unsuccessful half the step size; if 5 consecutive fails then return
798 Nconsecutive
= Nconsecutive
+1
800 PRINT*,'Warning: LU Decomposition returned ising = ',ising
801 IF (Nconsecutive
<= 5) THEN ! Less than 5 consecutive failed decompositions
803 ELSE ! More than 5 consecutive failed decompositions
805 END IF ! Nconsecutive
808 END DO ! WHILE Singular
810 END SUBROUTINE KPP_ROOT_ros_PrepareMatrix
813 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
814 SUBROUTINE KPP_ROOT_ros_Decomp( A
, Pivot
, ising
, Ndec
)
815 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
816 ! Template for the LU decomposition
817 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
819 !~~~> Inout variables
820 KPP_REAL
, INTENT(INOUT
) :: A(LU_NONZERO
)
821 !~~~> Output variables
822 INTEGER, INTENT(OUT
) :: Pivot(NVAR
), ising
823 INTEGER, INTENT(INOUT
) :: Ndec
826 CALL KPP_ROOT_DGETRF( NVAR
, NVAR
, A
, NVAR
, Pivot
, ising
)
828 CALL KPP_ROOT_KppDecomp ( A
, ising
)
833 END SUBROUTINE KPP_ROOT_ros_Decomp
836 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
837 SUBROUTINE KPP_ROOT_ros_Solve( A
, Pivot
, b
, Nsol
)
838 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
839 ! Template for the forward/backward substitution (using pre-computed LU decomposition)
840 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
842 !~~~> Input variables
843 KPP_REAL
, INTENT(IN
) :: A(LU_NONZERO
)
844 INTEGER, INTENT(IN
) :: Pivot(NVAR
)
846 INTEGER, INTENT(INOUT
) :: nsol
847 !~~~> InOut variables
848 KPP_REAL
, INTENT(INOUT
) :: b(NVAR
)
852 CALL KPP_ROOT_DGETRS( 'N', NVAR
, 1, A
, NVAR
, Pivot
, b
, NVAR
, 0 )
854 CALL KPP_ROOT_KppSolve( A
, b
)
859 END SUBROUTINE KPP_ROOT_ros_Solve
863 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
864 SUBROUTINE KPP_ROOT_Ros2 (ros_S
,ros_A
,ros_C
,ros_M
,ros_E
,ros_Alpha
,&
865 ros_Gamma
,ros_NewF
,ros_ELO
,ros_Name
)
866 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
867 ! --- AN L-STABLE METHOD, 2 stages, order 2
868 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
872 INTEGER, PARAMETER :: S
=2
873 INTEGER, INTENT(OUT
) :: ros_S
874 KPP_REAL
, DIMENSION(S
), INTENT(OUT
) :: ros_M
,ros_E
,ros_Alpha
,ros_Gamma
875 KPP_REAL
, DIMENSION(S
*(S
-1)/2), INTENT(OUT
) :: ros_A
, ros_C
876 KPP_REAL
, INTENT(OUT
) :: ros_ELO
877 LOGICAL, DIMENSION(S
), INTENT(OUT
) :: ros_NewF
878 CHARACTER(LEN
=12), INTENT(OUT
) :: ros_Name
879 !cms DOUBLE PRECISION g
882 g
= 1.0_dp
+ 1.0_dp
/SQRT(2.0_dp
)
884 !~~~> Name of the method
886 !~~~> Number of stages
889 !~~~> The coefficient matrices A and C are strictly lower triangular.
890 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
891 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
892 ! The general mapping formula is:
893 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
894 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
896 ros_A(1) = (1.0_dp
)/g
897 ros_C(1) = (-2.0_dp
)/g
898 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
899 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
902 !~~~> M_i = Coefficients for new step solution
903 ros_M(1)= (3.0_dp
)/(2.0_dp
*g
)
904 ros_M(2)= (1.0_dp
)/(2.0_dp
*g
)
905 ! E_i = Coefficients for error estimator
906 ros_E(1) = 1.0_dp
/(2.0_dp
*g
)
907 ros_E(2) = 1.0_dp
/(2.0_dp
*g
)
908 !~~~> ros_ELO = estimator of local order - the minimum between the
909 ! main and the embedded scheme orders plus one
911 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
912 ros_Alpha(1) = 0.0_dp
913 ros_Alpha(2) = 1.0_dp
914 !~~~> Gamma_i = \sum_j gamma_{i,j}
918 END SUBROUTINE KPP_ROOT_Ros2
921 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
922 SUBROUTINE KPP_ROOT_Ros3 (ros_S
,ros_A
,ros_C
,ros_M
,ros_E
,ros_Alpha
,&
923 ros_Gamma
,ros_NewF
,ros_ELO
,ros_Name
)
924 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
925 ! --- AN L-STABLE METHOD, 3 stages, order 3, 2 function evaluations
926 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
930 INTEGER, PARAMETER :: S
=3
931 INTEGER, INTENT(OUT
) :: ros_S
932 KPP_REAL
, DIMENSION(S
), INTENT(OUT
) :: ros_M
,ros_E
,ros_Alpha
,ros_Gamma
933 KPP_REAL
, DIMENSION(S
*(S
-1)/2), INTENT(OUT
) :: ros_A
, ros_C
934 KPP_REAL
, INTENT(OUT
) :: ros_ELO
935 LOGICAL, DIMENSION(S
), INTENT(OUT
) :: ros_NewF
936 CHARACTER(LEN
=12), INTENT(OUT
) :: ros_Name
938 !~~~> Name of the method
940 !~~~> Number of stages
943 !~~~> The coefficient matrices A and C are strictly lower triangular.
944 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
945 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
946 ! The general mapping formula is:
947 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
948 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
954 ros_C(1) = -0.10156171083877702091975600115545E+01_dp
955 ros_C(2) = 0.40759956452537699824805835358067E+01_dp
956 ros_C(3) = 0.92076794298330791242156818474003E+01_dp
957 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
958 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
961 ros_NewF(3) = .FALSE
.
962 !~~~> M_i = Coefficients for new step solution
963 ros_M(1) = 0.1E+01_dp
964 ros_M(2) = 0.61697947043828245592553615689730E+01_dp
965 ros_M(3) = -0.42772256543218573326238373806514E+00_dp
966 ! E_i = Coefficients for error estimator
967 ros_E(1) = 0.5E+00_dp
968 ros_E(2) = -0.29079558716805469821718236208017E+01_dp
969 ros_E(3) = 0.22354069897811569627360909276199E+00_dp
970 !~~~> ros_ELO = estimator of local order - the minimum between the
971 ! main and the embedded scheme orders plus 1
973 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
974 ros_Alpha(1)= 0.0E+00_dp
975 ros_Alpha(2)= 0.43586652150845899941601945119356E+00_dp
976 ros_Alpha(3)= 0.43586652150845899941601945119356E+00_dp
977 !~~~> Gamma_i = \sum_j gamma_{i,j}
978 ros_Gamma(1)= 0.43586652150845899941601945119356E+00_dp
979 ros_Gamma(2)= 0.24291996454816804366592249683314E+00_dp
980 ros_Gamma(3)= 0.21851380027664058511513169485832E+01_dp
982 END SUBROUTINE KPP_ROOT_Ros3
984 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
987 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
988 SUBROUTINE KPP_ROOT_Ros4 (ros_S
,ros_A
,ros_C
,ros_M
,ros_E
,ros_Alpha
,&
989 ros_Gamma
,ros_NewF
,ros_ELO
,ros_Name
)
990 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
991 ! L-STABLE ROSENBROCK METHOD OF ORDER 4, WITH 4 STAGES
992 ! L-STABLE EMBEDDED ROSENBROCK METHOD OF ORDER 3
994 ! E. HAIRER AND G. WANNER, SOLVING ORDINARY DIFFERENTIAL
995 ! EQUATIONS II. STIFF AND DIFFERENTIAL-ALGEBRAIC PROBLEMS.
996 ! SPRINGER SERIES IN COMPUTATIONAL MATHEMATICS,
997 ! SPRINGER-VERLAG (1990)
998 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1002 INTEGER, PARAMETER :: S
=4
1003 INTEGER, INTENT(OUT
) :: ros_S
1004 KPP_REAL
, DIMENSION(4), INTENT(OUT
) :: ros_M
,ros_E
,ros_Alpha
,ros_Gamma
1005 KPP_REAL
, DIMENSION(6), INTENT(OUT
) :: ros_A
, ros_C
1006 KPP_REAL
, INTENT(OUT
) :: ros_ELO
1007 LOGICAL, DIMENSION(4), INTENT(OUT
) :: ros_NewF
1008 CHARACTER(LEN
=12), INTENT(OUT
) :: ros_Name
1009 !cms DOUBLE PRECISION g
1013 !~~~> Name of the method
1015 !~~~> Number of stages
1018 !~~~> The coefficient matrices A and C are strictly lower triangular.
1019 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
1020 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
1021 ! The general mapping formula is:
1022 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
1023 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
1025 ros_A(1) = 0.2000000000000000E+01_dp
1026 ros_A(2) = 0.1867943637803922E+01_dp
1027 ros_A(3) = 0.2344449711399156E+00_dp
1032 ros_C(1) =-0.7137615036412310E+01_dp
1033 ros_C(2) = 0.2580708087951457E+01_dp
1034 ros_C(3) = 0.6515950076447975E+00_dp
1035 ros_C(4) =-0.2137148994382534E+01_dp
1036 ros_C(5) =-0.3214669691237626E+00_dp
1037 ros_C(6) =-0.6949742501781779E+00_dp
1038 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
1039 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
1040 ros_NewF(1) = .TRUE
.
1041 ros_NewF(2) = .TRUE
.
1042 ros_NewF(3) = .TRUE
.
1043 ros_NewF(4) = .FALSE
.
1044 !~~~> M_i = Coefficients for new step solution
1045 ros_M(1) = 0.2255570073418735E+01_dp
1046 ros_M(2) = 0.2870493262186792E+00_dp
1047 ros_M(3) = 0.4353179431840180E+00_dp
1048 ros_M(4) = 0.1093502252409163E+01_dp
1049 !~~~> E_i = Coefficients for error estimator
1050 ros_E(1) =-0.2815431932141155E+00_dp
1051 ros_E(2) =-0.7276199124938920E-01_dp
1052 ros_E(3) =-0.1082196201495311E+00_dp
1053 ros_E(4) =-0.1093502252409163E+01_dp
1054 !~~~> ros_ELO = estimator of local order - the minimum between the
1055 ! main and the embedded scheme orders plus 1
1057 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
1058 ros_Alpha(1) = 0.0_dp
1059 ros_Alpha(2) = 0.1145640000000000E+01_dp
1060 ros_Alpha(3) = 0.6552168638155900E+00_dp
1061 ros_Alpha(4) = ros_Alpha(3)
1062 !~~~> Gamma_i = \sum_j gamma_{i,j}
1063 ros_Gamma(1) = 0.5728200000000000E+00_dp
1064 ros_Gamma(2) =-0.1769193891319233E+01_dp
1065 ros_Gamma(3) = 0.7592633437920482E+00_dp
1066 ros_Gamma(4) =-0.1049021087100450E+00_dp
1068 END SUBROUTINE KPP_ROOT_Ros4
1070 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1071 SUBROUTINE KPP_ROOT_Rodas3 (ros_S
,ros_A
,ros_C
,ros_M
,ros_E
,ros_Alpha
,&
1072 ros_Gamma
,ros_NewF
,ros_ELO
,ros_Name
)
1073 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1074 ! --- A STIFFLY-STABLE METHOD, 4 stages, order 3
1075 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1079 INTEGER, PARAMETER :: S
=4
1080 INTEGER, INTENT(OUT
) :: ros_S
1081 KPP_REAL
, DIMENSION(S
), INTENT(OUT
) :: ros_M
,ros_E
,ros_Alpha
,ros_Gamma
1082 KPP_REAL
, DIMENSION(S
*(S
-1)/2), INTENT(OUT
) :: ros_A
, ros_C
1083 KPP_REAL
, INTENT(OUT
) :: ros_ELO
1084 LOGICAL, DIMENSION(S
), INTENT(OUT
) :: ros_NewF
1085 CHARACTER(LEN
=12), INTENT(OUT
) :: ros_Name
1086 !cms DOUBLE PRECISION g
1089 !~~~> Name of the method
1090 ros_Name
= 'RODAS-3'
1091 !~~~> Number of stages
1094 !~~~> The coefficient matrices A and C are strictly lower triangular.
1095 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
1096 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
1097 ! The general mapping formula is:
1098 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
1099 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
1101 ros_A(1) = 0.0E+00_dp
1102 ros_A(2) = 2.0E+00_dp
1103 ros_A(3) = 0.0E+00_dp
1104 ros_A(4) = 2.0E+00_dp
1105 ros_A(5) = 0.0E+00_dp
1106 ros_A(6) = 1.0E+00_dp
1108 ros_C(1) = 4.0E+00_dp
1109 ros_C(2) = 1.0E+00_dp
1110 ros_C(3) =-1.0E+00_dp
1111 ros_C(4) = 1.0E+00_dp
1112 ros_C(5) =-1.0E+00_dp
1113 ros_C(6) =-(8.0E+00_dp
/3.0E+00_dp
)
1115 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
1116 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
1117 ros_NewF(1) = .TRUE
.
1118 ros_NewF(2) = .FALSE
.
1119 ros_NewF(3) = .TRUE
.
1120 ros_NewF(4) = .TRUE
.
1121 !~~~> M_i = Coefficients for new step solution
1122 ros_M(1) = 2.0E+00_dp
1123 ros_M(2) = 0.0E+00_dp
1124 ros_M(3) = 1.0E+00_dp
1125 ros_M(4) = 1.0E+00_dp
1126 !~~~> E_i = Coefficients for error estimator
1127 ros_E(1) = 0.0E+00_dp
1128 ros_E(2) = 0.0E+00_dp
1129 ros_E(3) = 0.0E+00_dp
1130 ros_E(4) = 1.0E+00_dp
1131 !~~~> ros_ELO = estimator of local order - the minimum between the
1132 ! main and the embedded scheme orders plus 1
1133 ros_ELO
= 3.0E+00_dp
1134 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
1135 ros_Alpha(1) = 0.0E+00_dp
1136 ros_Alpha(2) = 0.0E+00_dp
1137 ros_Alpha(3) = 1.0E+00_dp
1138 ros_Alpha(4) = 1.0E+00_dp
1139 !~~~> Gamma_i = \sum_j gamma_{i,j}
1140 ros_Gamma(1) = 0.5E+00_dp
1141 ros_Gamma(2) = 1.5E+00_dp
1142 ros_Gamma(3) = 0.0E+00_dp
1143 ros_Gamma(4) = 0.0E+00_dp
1145 END SUBROUTINE KPP_ROOT_Rodas3
1147 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1148 SUBROUTINE KPP_ROOT_Rodas4 (ros_S
,ros_A
,ros_C
,ros_M
,ros_E
,ros_Alpha
,&
1149 ros_Gamma
,ros_NewF
,ros_ELO
,ros_Name
)
1150 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1151 ! STIFFLY-STABLE ROSENBROCK METHOD OF ORDER 4, WITH 6 STAGES
1153 ! E. HAIRER AND G. WANNER, SOLVING ORDINARY DIFFERENTIAL
1154 ! EQUATIONS II. STIFF AND DIFFERENTIAL-ALGEBRAIC PROBLEMS.
1155 ! SPRINGER SERIES IN COMPUTATIONAL MATHEMATICS,
1156 ! SPRINGER-VERLAG (1996)
1157 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1161 INTEGER, PARAMETER :: S
=6
1162 INTEGER, INTENT(OUT
) :: ros_S
1163 KPP_REAL
, DIMENSION(S
), INTENT(OUT
) :: ros_M
,ros_E
,ros_Alpha
,ros_Gamma
1164 KPP_REAL
, DIMENSION(S
*(S
-1)/2), INTENT(OUT
) :: ros_A
, ros_C
1165 KPP_REAL
, INTENT(OUT
) :: ros_ELO
1166 LOGICAL, DIMENSION(S
), INTENT(OUT
) :: ros_NewF
1167 CHARACTER(LEN
=12), INTENT(OUT
) :: ros_Name
1168 !cms DOUBLE PRECISION g
1171 !~~~> Name of the method
1172 ros_Name
= 'RODAS-4'
1173 !~~~> Number of stages
1176 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
1177 ros_Alpha(1) = 0.000_dp
1178 ros_Alpha(2) = 0.386_dp
1179 ros_Alpha(3) = 0.210_dp
1180 ros_Alpha(4) = 0.630_dp
1181 ros_Alpha(5) = 1.000_dp
1182 ros_Alpha(6) = 1.000_dp
1184 !~~~> Gamma_i = \sum_j gamma_{i,j}
1185 ros_Gamma(1) = 0.2500000000000000E+00_dp
1186 ros_Gamma(2) =-0.1043000000000000E+00_dp
1187 ros_Gamma(3) = 0.1035000000000000E+00_dp
1188 ros_Gamma(4) =-0.3620000000000023E-01_dp
1189 ros_Gamma(5) = 0.0_dp
1190 ros_Gamma(6) = 0.0_dp
1192 !~~~> The coefficient matrices A and C are strictly lower triangular.
1193 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
1194 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
1195 ! The general mapping formula is: A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
1196 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
1198 ros_A(1) = 0.1544000000000000E+01_dp
1199 ros_A(2) = 0.9466785280815826E+00_dp
1200 ros_A(3) = 0.2557011698983284E+00_dp
1201 ros_A(4) = 0.3314825187068521E+01_dp
1202 ros_A(5) = 0.2896124015972201E+01_dp
1203 ros_A(6) = 0.9986419139977817E+00_dp
1204 ros_A(7) = 0.1221224509226641E+01_dp
1205 ros_A(8) = 0.6019134481288629E+01_dp
1206 ros_A(9) = 0.1253708332932087E+02_dp
1207 ros_A(10) =-0.6878860361058950E+00_dp
1208 ros_A(11) = ros_A(7)
1209 ros_A(12) = ros_A(8)
1210 ros_A(13) = ros_A(9)
1211 ros_A(14) = ros_A(10)
1212 ros_A(15) = 1.0E+00_dp
1214 ros_C(1) =-0.5668800000000000E+01_dp
1215 ros_C(2) =-0.2430093356833875E+01_dp
1216 ros_C(3) =-0.2063599157091915E+00_dp
1217 ros_C(4) =-0.1073529058151375E+00_dp
1218 ros_C(5) =-0.9594562251023355E+01_dp
1219 ros_C(6) =-0.2047028614809616E+02_dp
1220 ros_C(7) = 0.7496443313967647E+01_dp
1221 ros_C(8) =-0.1024680431464352E+02_dp
1222 ros_C(9) =-0.3399990352819905E+02_dp
1223 ros_C(10) = 0.1170890893206160E+02_dp
1224 ros_C(11) = 0.8083246795921522E+01_dp
1225 ros_C(12) =-0.7981132988064893E+01_dp
1226 ros_C(13) =-0.3152159432874371E+02_dp
1227 ros_C(14) = 0.1631930543123136E+02_dp
1228 ros_C(15) =-0.6058818238834054E+01_dp
1230 !~~~> M_i = Coefficients for new step solution
1234 ros_M(4) = ros_A(10)
1235 ros_M(5) = 1.0E+00_dp
1236 ros_M(6) = 1.0E+00_dp
1238 !~~~> E_i = Coefficients for error estimator
1239 ros_E(1) = 0.0E+00_dp
1240 ros_E(2) = 0.0E+00_dp
1241 ros_E(3) = 0.0E+00_dp
1242 ros_E(4) = 0.0E+00_dp
1243 ros_E(5) = 0.0E+00_dp
1244 ros_E(6) = 1.0E+00_dp
1246 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
1247 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
1248 ros_NewF(1) = .TRUE
.
1249 ros_NewF(2) = .TRUE
.
1250 ros_NewF(3) = .TRUE
.
1251 ros_NewF(4) = .TRUE
.
1252 ros_NewF(5) = .TRUE
.
1253 ros_NewF(6) = .TRUE
.
1255 !~~~> ros_ELO = estimator of local order - the minimum between the
1256 ! main and the embedded scheme orders plus 1
1259 END SUBROUTINE KPP_ROOT_Rodas4
1261 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1262 ! End of the set of internal Rosenbrock subroutines
1263 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1264 END SUBROUTINE KPP_ROOT_Rosenbrock
1265 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1268 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1269 SUBROUTINE KPP_ROOT_FunTemplate( T
, Y
, Ydot
, RCONST
, FIX
, Nfun
)
1270 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1271 ! Template for the ODE function call.
1272 ! Updates the rate coefficients (and possibly the fixed species) at each call
1273 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1274 USE KPP_ROOT_Parameters
1275 !! USE KPP_ROOT_Global
1276 !! USE KPP_ROOT_Function
1277 !! USE KPP_ROOT_Rates
1278 !~~~> Input variables
1279 KPP_REAL
:: T
, Y(NVAR
)
1280 KPP_REAL
:: RCONST(NREACT
)
1281 KPP_REAL
:: FIX(NFIX
)
1282 !~~~> Output variables
1283 KPP_REAL
:: Ydot(NVAR
)
1287 !~~~> Local variables
1292 !! CALL Update_SUN()
1293 !! CALL Update_RCONST()
1294 CALL KPP_ROOT_Fun( Y
, FIX
, RCONST
, Ydot
)
1299 END SUBROUTINE KPP_ROOT_FunTemplate
1302 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1303 SUBROUTINE KPP_ROOT_JacTemplate( T
, Y
, Jcb
, FIX
, Njac
, RCONST
)
1304 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1305 ! Template for the ODE Jacobian call.
1306 ! Updates the rate coefficients (and possibly the fixed species) at each call
1307 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1308 USE KPP_ROOT_Parameters
1309 !!USE KPP_ROOT_Global
1310 USE KPP_ROOT_Jacobian
1311 !! USE KPP_ROOT_LinearAlgebra
1312 !! USE KPP_ROOT_Rates
1313 !~~~> Input variables
1314 KPP_REAL
:: T
, Y(NVAR
)
1315 KPP_REAL
:: FIX(NFIX
)
1316 KPP_REAL
:: RCONST(NREACT
)
1320 !~~~> Output variables
1322 KPP_REAL
:: JV(LU_NONZERO
), Jcb(NVAR
,NVAR
)
1324 KPP_REAL
:: Jcb(LU_NONZERO
)
1326 !~~~> Local variables
1334 !! CALL Update_SUN()
1335 !! CALL Update_RCONST()
1337 CALL KPP_ROOT_Jac_SP(Y
, FIX
, RCONST
, JV
)
1344 Jcb(LU_IROW(i
),LU_ICOL(i
)) = JV(i
)
1347 CALL KPP_ROOT_Jac_SP( Y
, FIX
, RCONST
, Jcb
)
1353 END SUBROUTINE KPP_ROOT_JacTemplate
1355 !!!END MODULE KPP_ROOT_Integrator