Original WRF subgrid support version from John Michalakes without fire
[wrffire.git] / wrfv2_fire / chem / KPP / kpp / kpp-2.1 / int / WRF_conform / TEMP / ros2w.f90
blob99ed0897b2e694e108d058a806e2839cc5033ed9
1 MODULE KPP_ROOT_Integrator
3 USE KPP_ROOT_Parameters
4 USE KPP_ROOT_Precision
5 USE KPP_ROOT_JacobianSP
7 IMPLICIT NONE
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
25 ' ', & ! 0 (not used)
26 'Success ' /) ! 1
28 CONTAINS
30 SUBROUTINE KPP_ROOT_INTEGRATE( TIN, TOUT, &
31 FIX, VAR, RCONST, ATOL, RTOL, &
32 ICNTRL_U, RCNTRL_U, ISTATUS_U, RSTATUS_U, IERR_U )
34 USE KPP_ROOT_Parameters
35 !! USE KPP_ROOT_Global
36 IMPLICIT NONE
37 KPP_REAL, INTENT(INOUT), DIMENSION(NFIX) :: FIX
38 KPP_REAL, INTENT(INOUT), DIMENSION(NVAR) :: VAR
39 KPP_REAL, INTENT(IN), DIMENSION(NSPEC) :: ATOL, RTOL
40 KPP_REAL, INTENT(IN), DIMENSION(NREACT) :: RCONST
41 KPP_REAL, INTENT(IN) :: TIN ! Start Time
42 KPP_REAL, INTENT(IN) :: 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)
48 INTEGER, INTENT(OUT), OPTIONAL :: IERR_U
50 KPP_REAL :: STEPMIN
53 INTEGER :: N_stp, N_acc, N_rej, N_sng
54 SAVE N_stp, N_acc, N_rej, N_sng
55 INTEGER :: i, IERR
56 KPP_REAL :: RCNTRL(20), RSTATUS(20)
57 INTEGER :: ICNTRL(20), ISTATUS(20)
59 ICNTRL(:) = 0
60 RCNTRL(:) = 0.0_dp
61 ISTATUS(:) = 0
62 RSTATUS(:) = 0.0_dp
64 ! If optional parameters are given, and if they are >0,
65 ! then they overwrite default settings.
66 IF (PRESENT(ICNTRL_U)) THEN
67 WHERE(ICNTRL_U(:) > 0) ICNTRL(:) = ICNTRL_U(:)
68 END IF
69 IF (PRESENT(RCNTRL_U)) THEN
70 WHERE(RCNTRL_U(:) > 0) RCNTRL(:) = RCNTRL_U(:)
71 END IF
73 CALL KPP_ROOT_Rosenbrock(VAR, FIX, RCONST, TIN,TOUT, &
74 ATOL,RTOL, &
75 RCNTRL,ICNTRL,RSTATUS,ISTATUS,IERR)
77 STEPMIN = RCNTRL(ihexit)
78 ! if optional parameters are given for output they to return information
79 !!$ IF (PRESENT(ISTATUS_U)) ISTATUS_U(:) = ISTATUS(:)
80 !!$ IF (PRESENT(RSTATUS_U)) RSTATUS_U(:) = RSTATUS(:)
81 !!$ IF (PRESENT(IERR_U)) IERR_U = IERR
83 END SUBROUTINE KPP_ROOT_INTEGRATE
85 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
86 SUBROUTINE KPP_ROOT_Rosenbrock(Y, FIX, RCONST, Tstart,Tend, &
87 AbsTol,RelTol, &
88 RCNTRL,ICNTRL,RSTATUS,ISTATUS,IERR)
89 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
91 ! Solves the system y'=F(t,y) using a Rosenbrock method defined by:
93 ! G = 1/(H*gamma(1)) - Jac(t0,Y0)
94 ! T_i = t0 + Alpha(i)*H
95 ! Y_i = Y0 + \sum_{j=1}^{i-1} A(i,j)*K_j
96 ! G * K_i = Fun( T_i, Y_i ) + \sum_{j=1}^S C(i,j)/H * K_j +
97 ! gamma(i)*dF/dT(t0, Y0)
98 ! Y1 = Y0 + \sum_{j=1}^S M(j)*K_j
100 ! For details on Rosenbrock methods and their implementation consult:
101 ! E. Hairer and G. Wanner
102 ! "Solving ODEs II. Stiff and differential-algebraic problems".
103 ! Springer series in computational mathematics, Springer-Verlag, 1996.
104 ! The codes contained in the book inspired this implementation.
106 ! (C) Adrian Sandu, August 2004
107 ! Virginia Polytechnic Institute and State University
108 ! Contact: sandu@cs.vt.edu
109 ! This implementation is part of KPP - the Kinetic PreProcessor
110 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
112 !~~~> INPUT ARGUMENTS:
114 !- Y(NVAR) = vector of initial conditions (at T=Tstart)
115 !- [Tstart,Tend] = time range of integration
116 ! (if Tstart>Tend the integration is performed backwards in time)
117 !- RelTol, AbsTol = user precribed accuracy
118 !- SUBROUTINE Fun( T, Y, Ydot ) = ODE function,
119 ! returns Ydot = Y' = F(T,Y)
120 !- SUBROUTINE Jac( T, Y, Jcb ) = Jacobian of the ODE function,
121 ! returns Jcb = dFun/dY
122 !- ICNTRL(1:20) = integer inputs parameters
123 !- RCNTRL(1:20) = real inputs parameters
124 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
126 !~~~> OUTPUT ARGUMENTS:
128 !- Y(NVAR) -> vector of final states (at T->Tend)
129 !- ISTATUS(1:20) -> integer output parameters
130 !- RSTATUS(1:20) -> real output parameters
131 !- IERR -> job status upon return
132 ! success (positive value) or
133 ! failure (negative value)
134 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
136 !~~~> INPUT PARAMETERS:
138 ! Note: For input parameters equal to zero the default values of the
139 ! corresponding variables are used.
141 ! ICNTRL(1) = 1: F = F(y) Independent of T (AUTONOMOUS)
142 ! = 0: F = F(t,y) Depends on T (NON-AUTONOMOUS)
144 ! ICNTRL(2) = 0: AbsTol, RelTol are NVAR-dimensional vectors
145 ! = 1: AbsTol, RelTol are scalars
147 ! ICNTRL(3) -> selection of a particular Rosenbrock method
148 ! = 0 : default method is Rodas3
149 ! = 1 : method is Ros2
150 ! = 2 : method is Ros3
151 ! = 3 : method is Ros4
152 ! = 4 : method is Rodas3
153 ! = 5: method is Rodas4
155 ! ICNTRL(4) -> maximum number of integration steps
156 ! For ICNTRL(4)=0) the default value of 100000 is used
158 ! RCNTRL(1) -> Hmin, lower bound for the integration step size
159 ! It is strongly recommended to keep Hmin = ZERO
160 ! RCNTRL(2) -> Hmax, upper bound for the integration step size
161 ! RCNTRL(3) -> Hstart, starting value for the integration step size
163 ! RCNTRL(4) -> FacMin, lower bound on step decrease factor (default=0.2)
164 ! RCNTRL(5) -> FacMax, upper bound on step increase factor (default=6)
165 ! RCNTRL(6) -> FacRej, step decrease factor after multiple rejections
166 ! (default=0.1)
167 ! RCNTRL(7) -> FacSafe, by which the new step is slightly smaller
168 ! than the predicted value (default=0.9)
169 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
171 !~~~> OUTPUT PARAMETERS:
173 ! Note: each call to Rosenbrock adds the current no. of fcn calls
174 ! to previous value of ISTATUS(1), and similar for the other params.
175 ! Set ISTATUS(1:20) = 0 before call to avoid this accumulation.
177 ! ISTATUS(1) = No. of function calls
178 ! ISTATUS(2) = No. of jacobian calls
179 ! ISTATUS(3) = No. of steps
180 ! ISTATUS(4) = No. of accepted steps
181 ! ISTATUS(5) = No. of rejected steps (except at the beginning)
182 ! ISTATUS(6) = No. of LU decompositions
183 ! ISTATUS(7) = No. of forward/backward substitutions
184 ! ISTATUS(8) = No. of singular matrix decompositions
186 ! RSTATUS(1) -> Texit, the time corresponding to the
187 ! computed Y upon return
188 ! RSTATUS(2) -> Hexit, last accepted step before exit
189 ! For multiple restarts, use Hexit as Hstart in the following run
190 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
192 USE KPP_ROOT_Parameters
193 !! USE KPP_ROOT_LinearAlgebra
194 IMPLICIT NONE
196 !~~~> Arguments
197 KPP_REAL, INTENT(INOUT) :: Y(NVAR)
198 KPP_REAL, INTENT(IN), DIMENSION(NFIX) :: FIX
199 KPP_REAL, INTENT(IN), DIMENSION(NREACT) :: RCONST
200 KPP_REAL, INTENT(IN) :: Tstart,Tend
201 KPP_REAL, INTENT(IN) :: AbsTol(NVAR),RelTol(NVAR)
202 INTEGER, INTENT(IN) :: ICNTRL(20)
203 KPP_REAL, INTENT(IN) :: RCNTRL(20)
204 INTEGER, INTENT(INOUT) :: ISTATUS(20)
205 KPP_REAL, INTENT(INOUT) :: RSTATUS(20)
206 INTEGER, INTENT(OUT) :: IERR
207 !~~~> The method parameters
208 INTEGER, PARAMETER :: Smax = 6
209 INTEGER :: Method, ros_S
210 KPP_REAL, DIMENSION(Smax) :: ros_M, ros_E, ros_Alpha, ros_Gamma
211 KPP_REAL, DIMENSION(Smax*(Smax-1)/2) :: ros_A, ros_C
212 KPP_REAL :: ros_ELO
213 LOGICAL, DIMENSION(Smax) :: ros_NewF
214 CHARACTER(LEN=12) :: ros_Name
216 !~~~> Statistics on the work performed by the Rosenbrock method
217 INTEGER :: Nfun,Njac,Nstp,Nacc,Nrej,Ndec,Nsol,Nsng
220 !~~~> Local variables
221 KPP_REAL :: Roundoff, FacMin, FacMax, FacRej, FacSafe
222 KPP_REAL :: Hmin, Hmax, Hstart, Hexit
223 KPP_REAL :: Texit
224 INTEGER :: i, UplimTol, Max_no_steps
225 LOGICAL :: Autonomous, VectorTol
226 !~~~> Parameters
227 KPP_REAL, PARAMETER :: ZERO = 0.0_dp, ONE = 1.0_dp
228 KPP_REAL, PARAMETER :: DeltaMin = 1.0E-5_dp
230 !~~~> Initialize statistics
231 Nfun = ISTATUS(ifun)
232 Njac = ISTATUS(ijac)
233 Nstp = ISTATUS(istp)
234 Nacc = ISTATUS(iacc)
235 Nrej = ISTATUS(irej)
236 Ndec = ISTATUS(idec)
237 Nsol = ISTATUS(isol)
238 Nsng = ISTATUS(isng)
240 !~~~> Autonomous or time dependent ODE. Default is time dependent.
241 !!$ Autonomous = .NOT.(ICNTRL(1) == 0)
242 Autonomous=.FALSE.
247 !~~~> For Scalar tolerances (ICNTRL(2).NE.0) the code uses AbsTol(1) and RelTol(1)
248 ! For Vector tolerances (ICNTRL(2) == 0) the code uses AbsTol(1:NVAR) and RelTol(1:NVAR)
249 !!$ IF (ICNTRL(2) == 0) THEN
250 VectorTol = .TRUE.
251 UplimTol = NVAR
252 !!$ ELSE
253 !!$ VectorTol = .FALSE.
254 !!$ UplimTol = 1
255 !!$ END IF
257 !~~~> The particular Rosenbrock method chosen
258 !!$ IF (ICNTRL(3) == 0) THEN
259 !!$ Method = 4
260 !!$ ELSEIF ( (ICNTRL(3) >= 1).AND.(ICNTRL(3) <= 5) ) THEN
261 !!$ Method = ICNTRL(3)
262 !!$ ELSE
263 !!$ PRINT * , 'User-selected Rosenbrock method: ICNTRL(3)=', Method
264 !!$ CALL KPP_ROOT_ros_ErrorMsg(-2,Tstart,ZERO,IERR)
265 !!$ RETURN
266 !!$ END IF
268 !~~~> The maximum number of steps admitted
269 !!$ IF (ICNTRL(4) == 0) THEN
270 Max_no_steps = 10000
271 !!$ ELSEIF (ICNTRL(4) > 0) THEN
272 !!$ Max_no_steps=ICNTRL(4)
273 !!$ ELSE
274 !!$ PRINT * ,'User-selected max no. of steps: ICNTRL(4)=',ICNTRL(4)
275 !!$ CALL KPP_ROOT_ros_ErrorMsg(-1,Tstart,ZERO,IERR)
276 !!$ RETURN
277 !!$ END IF
279 !~~~> Unit roundoff (1+Roundoff>1)
280 Roundoff = KPP_ROOT_WLAMCH('E')
282 !~~~> Lower bound on the step size: (positive value)
283 !!$ IF (RCNTRL(1) == ZERO) THEN
284 Hmin = ZERO
285 !!$ ELSEIF (RCNTRL(1) > ZERO) THEN
286 !!$ Hmin = RCNTRL(1)
287 !!$ ELSE
288 !!$ PRINT * , 'User-selected Hmin: RCNTRL(1)=', RCNTRL(1)
289 !!$ CALL KPP_ROOT_ros_ErrorMsg(-3,Tstart,ZERO,IERR)
290 !!$ RETURN
291 !!$ END IF
292 !!$!~~~> Upper bound on the step size: (positive value)
293 !!$ IF (RCNTRL(2) == ZERO) THEN
294 Hmax = ABS(Tend-Tstart)
295 !!$ ELSEIF (RCNTRL(2) > ZERO) THEN
296 !!$ Hmax = MIN(ABS(RCNTRL(2)),ABS(Tend-Tstart))
297 !!$ ELSE
298 !!$ PRINT * , 'User-selected Hmax: RCNTRL(2)=', RCNTRL(2)
299 !!$ CALL KPP_ROOT_ros_ErrorMsg(-3,Tstart,ZERO,IERR)
300 !!$ RETURN
301 !!$ END IF
302 !!$!~~~> Starting step size: (positive value)
303 !!$ IF (RCNTRL(3) == ZERO) THEN
304 Hstart = MAX(Hmin,DeltaMin)
305 !!$ ELSEIF (RCNTRL(3) > ZERO) THEN
306 !!$ Hstart = MIN(ABS(RCNTRL(3)),ABS(Tend-Tstart))
307 !!$ ELSE
308 !!$ PRINT * , 'User-selected Hstart: RCNTRL(3)=', RCNTRL(3)
309 !!$ CALL KPP_ROOT_ros_ErrorMsg(-3,Tstart,ZERO,IERR)
310 !!$ RETURN
311 !!$ END IF
312 !~~~> Step size can be changed s.t. FacMin < Hnew/Hexit < FacMax
313 !!$ IF (RCNTRL(4) == ZERO) THEN
314 FacMin = 0.2_dp
315 !!$ ELSEIF (RCNTRL(4) > ZERO) THEN
316 !!$ FacMin = RCNTRL(4)
317 !!$ ELSE
318 !!$ PRINT * , 'User-selected FacMin: RCNTRL(4)=', RCNTRL(4)
319 !!$ CALL KPP_ROOT_ros_ErrorMsg(-4,Tstart,ZERO,IERR)
320 !!$ RETURN
321 !!$ END IF
322 !!$ IF (RCNTRL(5) == ZERO) THEN
323 FacMax = 6.0_dp
324 !!$ ELSEIF (RCNTRL(5) > ZERO) THEN
325 !!$ FacMax = RCNTRL(5)
326 !!$ ELSE
327 !!$ PRINT * , 'User-selected FacMax: RCNTRL(5)=', RCNTRL(5)
328 !!$ CALL KPP_ROOT_ros_ErrorMsg(-4,Tstart,ZERO,IERR)
329 !!$ RETURN
330 !!$ END IF
331 !~~~> FacRej: Factor to decrease step after 2 succesive rejections
332 !!$ IF (RCNTRL(6) == ZERO) THEN
333 FacRej = 0.1_dp
334 !!$ ELSEIF (RCNTRL(6) > ZERO) THEN
335 !!$ FacRej = RCNTRL(6)
336 !!$ ELSE
337 !!$ PRINT * , 'User-selected FacRej: RCNTRL(6)=', RCNTRL(6)
338 !!$ CALL KPP_ROOT_ros_ErrorMsg(-4,Tstart,ZERO,IERR)
339 !!$ RETURN
340 !!$ END IF
341 !~~~> FacSafe: Safety Factor in the computation of new step size
342 !!$ IF (RCNTRL(7) == ZERO) THEN
343 FacSafe = 0.9_dp
344 !!$ ELSEIF (RCNTRL(7) > ZERO) THEN
345 !!$ FacSafe = RCNTRL(7)
346 !!$ ELSE
347 !!$ PRINT * , 'User-selected FacSafe: RCNTRL(7)=', RCNTRL(7)
348 !!$ CALL KPP_ROOT_ros_ErrorMsg(-4,Tstart,ZERO,IERR)
349 !!$ RETURN
350 !!$ END IF
351 !!$!~~~> Check if tolerances are reasonable
352 !!$ DO i=1,UplimTol
353 !!$ IF ( (AbsTol(i) <= ZERO) .OR. (RelTol(i) <= 10.0_dp*Roundoff) &
354 !!$ .OR. (RelTol(i) >= 1.0_dp) ) THEN
355 !!$ PRINT * , ' AbsTol(',i,') = ',AbsTol(i)
356 !!$ PRINT * , ' RelTol(',i,') = ',RelTol(i)
357 !!$ CALL KPP_ROOT_ros_ErrorMsg(-5,Tstart,ZERO,IERR)
358 !!$ RETURN
359 !!$ END IF
360 !!$ END DO
363 !~~~> Initialize the particular Rosenbrock method
364 !!$ SELECT CASE (Method)
365 !!$ CASE (1)
366 CALL KPP_ROOT_Ros2(ros_S, ros_A, ros_C, ros_M, ros_E, &
367 ros_Alpha, ros_Gamma, ros_NewF, ros_ELO, ros_Name)
368 !!$ CASE (2)
369 !!$ CALL KPP_ROOT_Ros3(ros_S, ros_A, ros_C, ros_M, ros_E, &
370 !!$ ros_Alpha, ros_Gamma, ros_NewF, ros_ELO, ros_Name)
371 !!$ CASE (3)
372 !!$ CALL KPP_ROOT_Ros4(ros_S, ros_A, ros_C, ros_M, ros_E, &
373 !!$ ros_Alpha, ros_Gamma, ros_NewF, ros_ELO, ros_Name)
374 !!$ CASE (4)
375 !!$ CALL KPP_ROOT_Rodas3(ros_S, ros_A, ros_C, ros_M, ros_E, &
376 !!$ ros_Alpha, ros_Gamma, ros_NewF, ros_ELO, ros_Name)
377 !!$ CASE (5)
378 !!$ CALL KPP_ROOT_Rodas4(ros_S, ros_A, ros_C, ros_M, ros_E, &
379 !!$ ros_Alpha, ros_Gamma, ros_NewF, ros_ELO, ros_Name)
380 !!$ CASE DEFAULT
381 !!$ PRINT * , 'Unknown Rosenbrock method: ICNTRL(4)=', Method
382 !!$ CALL KPP_ROOT_ros_ErrorMsg(-2,Tstart,ZERO,IERR)
383 !!$ RETURN
384 !!$ END SELECT
386 !~~~> CALL Rosenbrock method
387 CALL KPP_ROOT_ros_Integrator(Y,Tstart,Tend,Texit, &
388 AbsTol, RelTol, &
389 ! Rosenbrock method coefficients
390 ros_S, ros_M, ros_E, ros_A, ros_C, &
391 ros_Alpha, ros_Gamma, ros_ELO, ros_NewF, &
392 ! Integration parameters
393 Autonomous, VectorTol, Max_no_steps, &
394 Roundoff, Hmin, Hmax, Hstart, Hexit, &
395 FacMin, FacMax, FacRej, FacSafe, &
396 ! Error indicator
397 IERR, &
398 ! Statistics on the work performed by the Rosenbrock method
399 Nfun,Njac,Nstp,Nacc,Nrej,Ndec,Nsol,Nsng,&
400 !~~~>
401 RCONST, FIX &
405 !~~~> Collect run statistics
406 ISTATUS(ifun) = Nfun
407 ISTATUS(ijac) = Njac
408 ISTATUS(istp) = Nstp
409 ISTATUS(iacc) = Nacc
410 ISTATUS(irej) = Nrej
411 ISTATUS(idec) = Ndec
412 ISTATUS(isol) = Nsol
413 ISTATUS(isng) = Nsng
414 !~~~> Last T and H
415 RSTATUS(itexit) = Texit
416 RSTATUS(ihexit) = Hexit
418 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
419 CONTAINS ! SUBROUTINES internal to Rosenbrock
420 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
422 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
423 SUBROUTINE KPP_ROOT_ros_ErrorMsg(Code,T,H,IERR)
424 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
425 ! Handles all error messages
426 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
427 USE KPP_ROOT_Precision
429 KPP_REAL, INTENT(IN) :: T, H
430 INTEGER, INTENT(IN) :: Code
431 INTEGER, INTENT(OUT) :: IERR
433 IERR = Code
434 PRINT * , &
435 'Forced exit from Rosenbrock due to the following error:'
436 IF ((Code>=-8).AND.(Code<=-1)) THEN
437 PRINT *, IERR_NAMES(Code)
438 ELSE
439 PRINT *, 'Unknown Error code: ', Code
440 ENDIF
442 PRINT *, "T=", T, "and H=", H
444 END SUBROUTINE KPP_ROOT_ros_ErrorMsg
446 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
447 SUBROUTINE KPP_ROOT_ros_Integrator (Y, Tstart, Tend, T, &
448 AbsTol, RelTol, &
449 !~~~> Rosenbrock method coefficients
450 ros_S, ros_M, ros_E, ros_A, ros_C, &
451 ros_Alpha, ros_Gamma, ros_ELO, ros_NewF, &
452 !~~~> Integration parameters
453 Autonomous, VectorTol, Max_no_steps, &
454 Roundoff, Hmin, Hmax, Hstart, Hexit, &
455 FacMin, FacMax, FacRej, FacSafe, &
456 !~~~> Error indicator
457 IERR, &
458 !~~~> Statistics on the work performed by the Rosenbrock method
459 Nfun,Njac,Nstp,Nacc,Nrej,Ndec,Nsol,Nsng, &
460 !~~~>
461 RCONST, FIX &
463 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
464 ! Template for the implementation of a generic Rosenbrock method
465 ! defined by ros_S (no of stages)
466 ! and its coefficients ros_{A,C,M,E,Alpha,Gamma}
467 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
469 IMPLICIT NONE
471 !~~~> Input: the initial condition at Tstart; Output: the solution at T
472 KPP_REAL, INTENT(INOUT) :: Y(NVAR)
473 !~~~> Input: integration interval
474 KPP_REAL, INTENT(IN) :: Tstart,Tend
475 !~~~> Output: time at which the solution is returned (T=Tend if success)
476 KPP_REAL, INTENT(OUT) :: T
477 !~~~> Input: tolerances
478 KPP_REAL, INTENT(IN) :: AbsTol(NVAR), RelTol(NVAR)
479 !~~~> Input: The Rosenbrock method parameters
480 INTEGER, INTENT(IN) :: ros_S
481 KPP_REAL, INTENT(IN) :: ros_M(ros_S), ros_E(ros_S), &
482 ros_Alpha(ros_S), ros_A(ros_S*(ros_S-1)/2), &
483 ros_Gamma(ros_S), ros_C(ros_S*(ros_S-1)/2), ros_ELO
484 LOGICAL, INTENT(IN) :: ros_NewF(ros_S)
485 !~~~> Input: integration parameters
486 LOGICAL, INTENT(IN) :: Autonomous, VectorTol
487 KPP_REAL, INTENT(IN) :: Hstart, Hmin, Hmax
488 INTEGER, INTENT(IN) :: Max_no_steps
489 KPP_REAL, INTENT(IN) :: Roundoff, FacMin, FacMax, FacRej, FacSafe
490 !~~~> Output: last accepted step
491 KPP_REAL, INTENT(OUT) :: Hexit
492 !~~~> Output: Error indicator
493 INTEGER, INTENT(OUT) :: IERR
494 !~~~> Input
495 KPP_REAL, INTENT(IN), DIMENSION(NFIX) :: FIX
496 !~~~> Input
497 KPP_REAL, INTENT(IN), DIMENSION(NREACT) :: RCONST
499 !~~~> Statistics on the work performed by the Rosenbrock method
500 INTEGER, INTENT(INOUT) :: Nfun,Njac,Nstp,Nacc,Nrej,Ndec,Nsol,Nsng
502 ! ~~~~ Local variables
503 KPP_REAL :: Ynew(NVAR), Fcn0(NVAR), Fcn(NVAR)
504 KPP_REAL :: K(NVAR*ros_S), dFdT(NVAR)
505 #ifdef FULL_ALGEBRA
506 KPP_REAL :: Jac0(NVAR,NVAR), Ghimj(NVAR,NVAR)
507 #else
508 KPP_REAL :: Jac0(LU_NONZERO), Ghimj(LU_NONZERO)
509 #endif
510 KPP_REAL :: H, Hnew, HC, HG, Fac, Tau
511 KPP_REAL :: Err, Yerr(NVAR)
512 INTEGER :: Pivot(NVAR), Direction, ioffset, j, istage
513 LOGICAL :: RejectLastH, RejectMoreH, Singular
514 !~~~> Local parameters
515 KPP_REAL, PARAMETER :: ZERO = 0.0_dp, ONE = 1.0_dp
516 KPP_REAL, PARAMETER :: DeltaMin = 1.0E-5_dp
517 !~~~> Locally called functions
518 ! KPP_REAL WLAMCH
519 ! EXTERNAL WLAMCH
520 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
523 !~~~> Initial preparations
524 T = Tstart
525 Hexit = 0.0_dp
526 H = MIN(Hstart,Hmax)
527 IF (ABS(H) <= 10.0_dp*Roundoff) H = DeltaMin
529 IF (Tend >= Tstart) THEN
530 Direction = +1
531 ELSE
532 Direction = -1
533 END IF
535 RejectLastH=.FALSE.
536 RejectMoreH=.FALSE.
538 !~~~> Time loop begins below
540 TimeLoop: DO WHILE ( (Direction > 0).AND.((T-Tend)+Roundoff <= ZERO) &
541 .OR. (Direction < 0).AND.((Tend-T)+Roundoff <= ZERO) )
543 IF ( Nstp > Max_no_steps ) THEN ! Too many steps
544 CALL KPP_ROOT_ros_ErrorMsg(-6,T,H,IERR)
545 RETURN
546 END IF
547 IF ( ((T+0.1_dp*H) == T).OR.(H <= Roundoff) ) THEN ! Step size too small
548 CALL KPP_ROOT_ros_ErrorMsg(-7,T,H,IERR)
549 RETURN
550 END IF
552 !~~~> Limit H if necessary to avoid going beyond Tend
553 Hexit = H
554 H = MIN(H,ABS(Tend-T))
556 !~~~> Compute the function at current time
557 CALL KPP_ROOT_FunTemplate(T,Y,Fcn0, RCONST, FIX, Nfun)
559 !~~~> Compute the function derivative with respect to T
560 IF (.NOT.Autonomous) THEN
561 CALL KPP_ROOT_ros_FunTimeDeriv ( T, Roundoff, Y, &
562 Fcn0, dFdT, RCONST, FIX, Nfun )
563 END IF
565 !~~~> Compute the Jacobian at current time
566 CALL KPP_ROOT_JacTemplate(T,Y,Jac0, FIX, Njac, RCONST)
568 !~~~> Repeat step calculation until current step accepted
569 UntilAccepted: DO
571 CALL KPP_ROOT_ros_PrepareMatrix(H,Direction,ros_Gamma(1), &
572 Jac0,Ghimj,Pivot,Singular, Ndec, Nsng )
573 IF (Singular) THEN ! More than 5 consecutive failed decompositions
574 CALL KPP_ROOT_ros_ErrorMsg(-8,T,H,IERR)
575 RETURN
576 END IF
578 !~~~> Compute the stages
579 Stage: DO istage = 1, ros_S
581 ! Current istage offset. Current istage vector is K(ioffset+1:ioffset+NVAR)
582 ioffset = NVAR*(istage-1)
584 ! For the 1st istage the function has been computed previously
585 IF ( istage == 1 ) THEN
586 CALL KPP_ROOT_WCOPY(NVAR,Fcn0,1,Fcn,1)
587 ! istage>1 and a new function evaluation is needed at the current istage
588 ELSEIF ( ros_NewF(istage) ) THEN
589 CALL KPP_ROOT_WCOPY(NVAR,Y,1,Ynew,1)
590 DO j = 1, istage-1
591 CALL KPP_ROOT_WAXPY(NVAR,ros_A((istage-1)*(istage-2)/2+j), &
592 K(NVAR*(j-1)+1),1,Ynew,1)
593 END DO
594 Tau = T + ros_Alpha(istage)*Direction*H
595 CALL KPP_ROOT_FunTemplate(Tau,Ynew,Fcn, RCONST, FIX, Nfun)
596 END IF ! if istage == 1 elseif ros_NewF(istage)
597 CALL KPP_ROOT_WCOPY(NVAR,Fcn,1,K(ioffset+1),1)
598 DO j = 1, istage-1
599 HC = ros_C((istage-1)*(istage-2)/2+j)/(Direction*H)
600 CALL KPP_ROOT_WAXPY(NVAR,HC,K(NVAR*(j-1)+1),1,K(ioffset+1),1)
601 END DO
602 IF ((.NOT. Autonomous).AND.(ros_Gamma(istage).NE.ZERO)) THEN
603 HG = Direction*H*ros_Gamma(istage)
604 CALL KPP_ROOT_WAXPY(NVAR,HG,dFdT,1,K(ioffset+1),1)
605 END IF
606 CALL KPP_ROOT_ros_Solve(Ghimj, Pivot, K(ioffset+1), Nsol)
608 END DO Stage
611 !~~~> Compute the new solution
612 CALL KPP_ROOT_WCOPY(NVAR,Y,1,Ynew,1)
613 DO j=1,ros_S
614 CALL KPP_ROOT_WAXPY(NVAR,ros_M(j),K(NVAR*(j-1)+1),1,Ynew,1)
615 END DO
617 !~~~> Compute the error estimation
618 CALL KPP_ROOT_WSCAL(NVAR,ZERO,Yerr,1)
619 DO j=1,ros_S
620 CALL KPP_ROOT_WAXPY(NVAR,ros_E(j),K(NVAR*(j-1)+1),1,Yerr,1)
621 END DO
622 Err = KPP_ROOT_ros_ErrorNorm ( Y, Ynew, Yerr, AbsTol, RelTol, VectorTol )
624 !~~~> New step size is bounded by FacMin <= Hnew/H <= FacMax
625 Fac = MIN(FacMax,MAX(FacMin,FacSafe/Err**(ONE/ros_ELO)))
626 Hnew = H*Fac
628 !~~~> Check the error magnitude and adjust step size
629 Nstp = Nstp+1
630 IF ( (Err <= ONE).OR.(H <= Hmin) ) THEN !~~~> Accept step
631 Nacc = Nacc+1
632 CALL KPP_ROOT_WCOPY(NVAR,Ynew,1,Y,1)
633 T = T + Direction*H
634 Hnew = MAX(Hmin,MIN(Hnew,Hmax))
635 IF (RejectLastH) THEN ! No step size increase after a rejected step
636 Hnew = MIN(Hnew,H)
637 END IF
638 RejectLastH = .FALSE.
639 RejectMoreH = .FALSE.
640 H = Hnew
641 EXIT UntilAccepted ! EXIT THE LOOP: WHILE STEP NOT ACCEPTED
642 ELSE !~~~> Reject step
643 IF (RejectMoreH) THEN
644 Hnew = H*FacRej
645 END IF
646 RejectMoreH = RejectLastH
647 RejectLastH = .TRUE.
648 H = Hnew
649 IF (Nacc >= 1) THEN
650 Nrej = Nrej+1
651 END IF
652 END IF ! Err <= 1
654 END DO UntilAccepted
656 END DO TimeLoop
658 !~~~> Succesful exit
659 IERR = 1 !~~~> The integration was successful
661 END SUBROUTINE KPP_ROOT_ros_Integrator
664 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
665 KPP_REAL FUNCTION KPP_ROOT_ros_ErrorNorm ( Y, Ynew, Yerr, &
666 AbsTol, RelTol, VectorTol )
667 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
668 !~~~> Computes the "scaled norm" of the error vector Yerr
669 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
670 IMPLICIT NONE
672 ! Input arguments
673 KPP_REAL, INTENT(IN) :: Y(NVAR), Ynew(NVAR), &
674 Yerr(NVAR), AbsTol(NVAR), RelTol(NVAR)
675 LOGICAL, INTENT(IN) :: VectorTol
676 ! Local variables
677 KPP_REAL :: Err, Scale, Ymax
678 INTEGER :: i
679 KPP_REAL, PARAMETER :: ZERO = 0.0_dp
681 Err = ZERO
682 DO i=1,NVAR
683 Ymax = MAX(ABS(Y(i)),ABS(Ynew(i)))
684 IF (VectorTol) THEN
685 Scale = AbsTol(i)+RelTol(i)*Ymax
686 ELSE
687 Scale = AbsTol(1)+RelTol(1)*Ymax
688 END IF
689 Err = Err+(Yerr(i)/Scale)**2
690 END DO
691 Err = SQRT(Err/NVAR)
693 KPP_ROOT_ros_ErrorNorm = Err
695 END FUNCTION KPP_ROOT_ros_ErrorNorm
698 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
699 SUBROUTINE KPP_ROOT_ros_FunTimeDeriv ( T, Roundoff, Y, &
700 Fcn0, dFdT, RCONST, FIX, Nfun )
701 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
702 !~~~> The time partial derivative of the function by finite differences
703 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
704 IMPLICIT NONE
706 !~~~> Input arguments
707 KPP_REAL, INTENT(IN) :: T, Roundoff, Y(NVAR), Fcn0(NVAR)
708 KPP_REAL, INTENT(IN) :: RCONST(NREACT), FIX(NFIX)
709 !~~~> Output arguments
710 KPP_REAL, INTENT(OUT) :: dFdT(NVAR)
711 !~~~> InOut args
712 INTEGER, INTENT(INOUT) ::Nfun
713 !~~~> Local variables
714 KPP_REAL :: Delta
715 KPP_REAL, PARAMETER :: ONE = 1.0_dp, DeltaMin = 1.0E-6_dp
717 Delta = SQRT(Roundoff)*MAX(DeltaMin,ABS(T))
718 CALL KPP_ROOT_FunTemplate(T+Delta,Y,dFdT, RCONST, FIX, Nfun)
719 CALL KPP_ROOT_WAXPY(NVAR,(-ONE),Fcn0,1,dFdT,1)
720 CALL KPP_ROOT_WSCAL(NVAR,(ONE/Delta),dFdT,1)
722 END SUBROUTINE KPP_ROOT_ros_FunTimeDeriv
725 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
726 SUBROUTINE KPP_ROOT_ros_PrepareMatrix ( H, Direction, gam, &
727 Jac0, Ghimj, Pivot, Singular, Ndec, Nsng )
728 ! --- --- --- --- --- --- --- --- --- --- --- --- ---
729 ! Prepares the LHS matrix for stage calculations
730 ! 1. Construct Ghimj = 1/(H*ham) - Jac0
731 ! "(Gamma H) Inverse Minus Jacobian"
732 ! 2. Repeat LU decomposition of Ghimj until successful.
733 ! -half the step size if LU decomposition fails and retry
734 ! -exit after 5 consecutive fails
735 ! --- --- --- --- --- --- --- --- --- --- --- --- ---
736 IMPLICIT NONE
738 !~~~> Input arguments
739 #ifdef FULL_ALGEBRA
740 KPP_REAL, INTENT(IN) :: Jac0(NVAR,NVAR)
741 #else
742 KPP_REAL, INTENT(IN) :: Jac0(LU_NONZERO)
743 #endif
744 KPP_REAL, INTENT(IN) :: gam
745 INTEGER, INTENT(IN) :: Direction
746 !~~~> Output arguments
747 #ifdef FULL_ALGEBRA
748 KPP_REAL, INTENT(OUT) :: Ghimj(NVAR,NVAR)
749 #else
750 KPP_REAL, INTENT(OUT) :: Ghimj(LU_NONZERO)
751 #endif
752 LOGICAL, INTENT(OUT) :: Singular
753 INTEGER, INTENT(OUT) :: Pivot(NVAR)
754 !~~~> Inout arguments
755 KPP_REAL, INTENT(INOUT) :: H ! step size is decreased when LU fails
756 INTEGER, INTENT(INOUT) :: Ndec, Nsng
757 !~~~> Local variables
758 INTEGER :: i, ising, Nconsecutive
759 KPP_REAL :: ghinv
760 KPP_REAL, PARAMETER :: ONE = 1.0_dp, HALF = 0.5_dp
762 Nconsecutive = 0
763 Singular = .TRUE.
765 DO WHILE (Singular)
767 !~~~> Construct Ghimj = 1/(H*gam) - Jac0
768 #ifdef FULL_ALGEBRA
769 CALL KPP_ROOT_WCOPY(NVAR*NVAR,Jac0,1,Ghimj,1)
770 CALL KPP_ROOT_WSCAL(NVAR*NVAR,(-ONE),Ghimj,1)
771 ghinv = ONE/(Direction*H*gam)
772 DO i=1,NVAR
773 Ghimj(i,i) = Ghimj(i,i)+ghinv
774 END DO
775 #else
776 CALL KPP_ROOT_WCOPY(LU_NONZERO,Jac0,1,Ghimj,1)
777 CALL KPP_ROOT_WSCAL(LU_NONZERO,(-ONE),Ghimj,1)
778 ghinv = ONE/(Direction*H*gam)
779 DO i=1,NVAR
780 Ghimj(LU_DIAG(i)) = Ghimj(LU_DIAG(i))+ghinv
781 END DO
782 #endif
783 !~~~> Compute LU decomposition
784 CALL KPP_ROOT_ros_Decomp( Ghimj, Pivot, ising, Ndec )
785 IF (ising == 0) THEN
786 !~~~> If successful done
787 Singular = .FALSE.
788 ELSE ! ising .ne. 0
789 !~~~> If unsuccessful half the step size; if 5 consecutive fails then return
790 Nsng = Nsng+1
791 Nconsecutive = Nconsecutive+1
792 Singular = .TRUE.
793 PRINT*,'Warning: LU Decomposition returned ising = ',ising
794 IF (Nconsecutive <= 5) THEN ! Less than 5 consecutive failed decompositions
795 H = H*HALF
796 ELSE ! More than 5 consecutive failed decompositions
797 RETURN
798 END IF ! Nconsecutive
799 END IF ! ising
801 END DO ! WHILE Singular
803 END SUBROUTINE KPP_ROOT_ros_PrepareMatrix
806 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
807 SUBROUTINE KPP_ROOT_ros_Decomp( A, Pivot, ising, Ndec )
808 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
809 ! Template for the LU decomposition
810 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
811 IMPLICIT NONE
812 !~~~> Inout variables
813 KPP_REAL, INTENT(INOUT) :: A(LU_NONZERO)
814 !~~~> Output variables
815 INTEGER, INTENT(OUT) :: Pivot(NVAR), ising
816 INTEGER, INTENT(INOUT) :: Ndec
818 #ifdef FULL_ALGEBRA
819 CALL KPP_ROOT_DGETRF( NVAR, NVAR, A, NVAR, Pivot, ising )
820 #else
821 CALL KPP_ROOT_KppDecomp ( A, ising )
822 Pivot(1) = 1
823 #endif
824 Ndec = Ndec + 1
826 END SUBROUTINE KPP_ROOT_ros_Decomp
829 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
830 SUBROUTINE KPP_ROOT_ros_Solve( A, Pivot, b, Nsol )
831 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
832 ! Template for the forward/backward substitution (using pre-computed LU decomposition)
833 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
834 IMPLICIT NONE
835 !~~~> Input variables
836 KPP_REAL, INTENT(IN) :: A(LU_NONZERO)
837 INTEGER, INTENT(IN) :: Pivot(NVAR)
838 !~~~~> InOut args
839 INTEGER, INTENT(INOUT) :: nsol
840 !~~~> InOut variables
841 KPP_REAL, INTENT(INOUT) :: b(NVAR)
844 #ifdef FULL_ALGEBRA
845 CALL KPP_ROOT_DGETRS( 'N', NVAR , 1, A, NVAR, Pivot, b, NVAR, 0 )
846 #else
847 CALL KPP_ROOT_KppSolve( A, b )
848 #endif
850 Nsol = Nsol+1
852 END SUBROUTINE KPP_ROOT_ros_Solve
856 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
857 SUBROUTINE KPP_ROOT_Ros2 (ros_S,ros_A,ros_C,ros_M,ros_E,ros_Alpha,&
858 ros_Gamma,ros_NewF,ros_ELO,ros_Name)
859 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
860 ! --- AN L-STABLE METHOD, 2 stages, order 2
861 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
863 IMPLICIT NONE
865 INTEGER, PARAMETER :: S=2
866 INTEGER, INTENT(OUT) :: ros_S
867 KPP_REAL, DIMENSION(S), INTENT(OUT) :: ros_M,ros_E,ros_Alpha,ros_Gamma
868 KPP_REAL, DIMENSION(S*(S-1)/2), INTENT(OUT) :: ros_A, ros_C
869 KPP_REAL, INTENT(OUT) :: ros_ELO
870 LOGICAL, DIMENSION(S), INTENT(OUT) :: ros_NewF
871 CHARACTER(LEN=12), INTENT(OUT) :: ros_Name
872 !cms DOUBLE PRECISION g
873 KPP_REAL :: g
875 g = 1.0_dp + 1.0_dp/SQRT(2.0_dp)
877 !~~~> Name of the method
878 ros_Name = 'ROS-2'
879 !~~~> Number of stages
880 ros_S = S
882 !~~~> The coefficient matrices A and C are strictly lower triangular.
883 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
884 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
885 ! The general mapping formula is:
886 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
887 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
889 ros_A(1) = (1.0_dp)/g
890 ros_C(1) = (-2.0_dp)/g
891 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
892 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
893 ros_NewF(1) = .TRUE.
894 ros_NewF(2) = .TRUE.
895 !~~~> M_i = Coefficients for new step solution
896 ros_M(1)= (3.0_dp)/(2.0_dp*g)
897 ros_M(2)= (1.0_dp)/(2.0_dp*g)
898 ! E_i = Coefficients for error estimator
899 ros_E(1) = 1.0_dp/(2.0_dp*g)
900 ros_E(2) = 1.0_dp/(2.0_dp*g)
901 !~~~> ros_ELO = estimator of local order - the minimum between the
902 ! main and the embedded scheme orders plus one
903 ros_ELO = 2.0_dp
904 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
905 ros_Alpha(1) = 0.0_dp
906 ros_Alpha(2) = 1.0_dp
907 !~~~> Gamma_i = \sum_j gamma_{i,j}
908 ros_Gamma(1) = g
909 ros_Gamma(2) =-g
911 END SUBROUTINE KPP_ROOT_Ros2
914 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
915 SUBROUTINE KPP_ROOT_Ros3 (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, 3 stages, order 3, 2 function evaluations
919 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
921 IMPLICIT NONE
923 INTEGER, PARAMETER :: S=3
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
931 !~~~> Name of the method
932 ros_Name = 'ROS-3'
933 !~~~> Number of stages
934 ros_S = S
936 !~~~> The coefficient matrices A and C are strictly lower triangular.
937 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
938 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
939 ! The general mapping formula is:
940 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
941 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
943 ros_A(1)= 1.0_dp
944 ros_A(2)= 1.0_dp
945 ros_A(3)= 0.0_dp
947 ros_C(1) = -0.10156171083877702091975600115545E+01_dp
948 ros_C(2) = 0.40759956452537699824805835358067E+01_dp
949 ros_C(3) = 0.92076794298330791242156818474003E+01_dp
950 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
951 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
952 ros_NewF(1) = .TRUE.
953 ros_NewF(2) = .TRUE.
954 ros_NewF(3) = .FALSE.
955 !~~~> M_i = Coefficients for new step solution
956 ros_M(1) = 0.1E+01_dp
957 ros_M(2) = 0.61697947043828245592553615689730E+01_dp
958 ros_M(3) = -0.42772256543218573326238373806514E+00_dp
959 ! E_i = Coefficients for error estimator
960 ros_E(1) = 0.5E+00_dp
961 ros_E(2) = -0.29079558716805469821718236208017E+01_dp
962 ros_E(3) = 0.22354069897811569627360909276199E+00_dp
963 !~~~> ros_ELO = estimator of local order - the minimum between the
964 ! main and the embedded scheme orders plus 1
965 ros_ELO = 3.0_dp
966 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
967 ros_Alpha(1)= 0.0E+00_dp
968 ros_Alpha(2)= 0.43586652150845899941601945119356E+00_dp
969 ros_Alpha(3)= 0.43586652150845899941601945119356E+00_dp
970 !~~~> Gamma_i = \sum_j gamma_{i,j}
971 ros_Gamma(1)= 0.43586652150845899941601945119356E+00_dp
972 ros_Gamma(2)= 0.24291996454816804366592249683314E+00_dp
973 ros_Gamma(3)= 0.21851380027664058511513169485832E+01_dp
975 END SUBROUTINE KPP_ROOT_Ros3
977 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
980 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
981 SUBROUTINE KPP_ROOT_Ros4 (ros_S,ros_A,ros_C,ros_M,ros_E,ros_Alpha,&
982 ros_Gamma,ros_NewF,ros_ELO,ros_Name)
983 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
984 ! L-STABLE ROSENBROCK METHOD OF ORDER 4, WITH 4 STAGES
985 ! L-STABLE EMBEDDED ROSENBROCK METHOD OF ORDER 3
987 ! E. HAIRER AND G. WANNER, SOLVING ORDINARY DIFFERENTIAL
988 ! EQUATIONS II. STIFF AND DIFFERENTIAL-ALGEBRAIC PROBLEMS.
989 ! SPRINGER SERIES IN COMPUTATIONAL MATHEMATICS,
990 ! SPRINGER-VERLAG (1990)
991 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
993 IMPLICIT NONE
995 INTEGER, PARAMETER :: S=4
996 INTEGER, INTENT(OUT) :: ros_S
997 KPP_REAL, DIMENSION(4), INTENT(OUT) :: ros_M,ros_E,ros_Alpha,ros_Gamma
998 KPP_REAL, DIMENSION(6), INTENT(OUT) :: ros_A, ros_C
999 KPP_REAL, INTENT(OUT) :: ros_ELO
1000 LOGICAL, DIMENSION(4), INTENT(OUT) :: ros_NewF
1001 CHARACTER(LEN=12), INTENT(OUT) :: ros_Name
1002 !cms DOUBLE PRECISION g
1003 KPP_REAL :: g
1006 !~~~> Name of the method
1007 ros_Name = 'ROS-4'
1008 !~~~> Number of stages
1009 ros_S = S
1011 !~~~> The coefficient matrices A and C are strictly lower triangular.
1012 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
1013 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
1014 ! The general mapping formula is:
1015 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
1016 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
1018 ros_A(1) = 0.2000000000000000E+01_dp
1019 ros_A(2) = 0.1867943637803922E+01_dp
1020 ros_A(3) = 0.2344449711399156E+00_dp
1021 ros_A(4) = ros_A(2)
1022 ros_A(5) = ros_A(3)
1023 ros_A(6) = 0.0_dp
1025 ros_C(1) =-0.7137615036412310E+01_dp
1026 ros_C(2) = 0.2580708087951457E+01_dp
1027 ros_C(3) = 0.6515950076447975E+00_dp
1028 ros_C(4) =-0.2137148994382534E+01_dp
1029 ros_C(5) =-0.3214669691237626E+00_dp
1030 ros_C(6) =-0.6949742501781779E+00_dp
1031 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
1032 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
1033 ros_NewF(1) = .TRUE.
1034 ros_NewF(2) = .TRUE.
1035 ros_NewF(3) = .TRUE.
1036 ros_NewF(4) = .FALSE.
1037 !~~~> M_i = Coefficients for new step solution
1038 ros_M(1) = 0.2255570073418735E+01_dp
1039 ros_M(2) = 0.2870493262186792E+00_dp
1040 ros_M(3) = 0.4353179431840180E+00_dp
1041 ros_M(4) = 0.1093502252409163E+01_dp
1042 !~~~> E_i = Coefficients for error estimator
1043 ros_E(1) =-0.2815431932141155E+00_dp
1044 ros_E(2) =-0.7276199124938920E-01_dp
1045 ros_E(3) =-0.1082196201495311E+00_dp
1046 ros_E(4) =-0.1093502252409163E+01_dp
1047 !~~~> ros_ELO = estimator of local order - the minimum between the
1048 ! main and the embedded scheme orders plus 1
1049 ros_ELO = 4.0_dp
1050 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
1051 ros_Alpha(1) = 0.0_dp
1052 ros_Alpha(2) = 0.1145640000000000E+01_dp
1053 ros_Alpha(3) = 0.6552168638155900E+00_dp
1054 ros_Alpha(4) = ros_Alpha(3)
1055 !~~~> Gamma_i = \sum_j gamma_{i,j}
1056 ros_Gamma(1) = 0.5728200000000000E+00_dp
1057 ros_Gamma(2) =-0.1769193891319233E+01_dp
1058 ros_Gamma(3) = 0.7592633437920482E+00_dp
1059 ros_Gamma(4) =-0.1049021087100450E+00_dp
1061 END SUBROUTINE KPP_ROOT_Ros4
1063 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1064 SUBROUTINE KPP_ROOT_Rodas3 (ros_S,ros_A,ros_C,ros_M,ros_E,ros_Alpha,&
1065 ros_Gamma,ros_NewF,ros_ELO,ros_Name)
1066 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1067 ! --- A STIFFLY-STABLE METHOD, 4 stages, order 3
1068 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1070 IMPLICIT NONE
1072 INTEGER, PARAMETER :: S=4
1073 INTEGER, INTENT(OUT) :: ros_S
1074 KPP_REAL, DIMENSION(S), INTENT(OUT) :: ros_M,ros_E,ros_Alpha,ros_Gamma
1075 KPP_REAL, DIMENSION(S*(S-1)/2), INTENT(OUT) :: ros_A, ros_C
1076 KPP_REAL, INTENT(OUT) :: ros_ELO
1077 LOGICAL, DIMENSION(S), INTENT(OUT) :: ros_NewF
1078 CHARACTER(LEN=12), INTENT(OUT) :: ros_Name
1079 !cms DOUBLE PRECISION g
1080 KPP_REAL :: g
1082 !~~~> Name of the method
1083 ros_Name = 'RODAS-3'
1084 !~~~> Number of stages
1085 ros_S = S
1087 !~~~> The coefficient matrices A and C are strictly lower triangular.
1088 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
1089 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
1090 ! The general mapping formula is:
1091 ! A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
1092 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
1094 ros_A(1) = 0.0E+00_dp
1095 ros_A(2) = 2.0E+00_dp
1096 ros_A(3) = 0.0E+00_dp
1097 ros_A(4) = 2.0E+00_dp
1098 ros_A(5) = 0.0E+00_dp
1099 ros_A(6) = 1.0E+00_dp
1101 ros_C(1) = 4.0E+00_dp
1102 ros_C(2) = 1.0E+00_dp
1103 ros_C(3) =-1.0E+00_dp
1104 ros_C(4) = 1.0E+00_dp
1105 ros_C(5) =-1.0E+00_dp
1106 ros_C(6) =-(8.0E+00_dp/3.0E+00_dp)
1108 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
1109 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
1110 ros_NewF(1) = .TRUE.
1111 ros_NewF(2) = .FALSE.
1112 ros_NewF(3) = .TRUE.
1113 ros_NewF(4) = .TRUE.
1114 !~~~> M_i = Coefficients for new step solution
1115 ros_M(1) = 2.0E+00_dp
1116 ros_M(2) = 0.0E+00_dp
1117 ros_M(3) = 1.0E+00_dp
1118 ros_M(4) = 1.0E+00_dp
1119 !~~~> E_i = Coefficients for error estimator
1120 ros_E(1) = 0.0E+00_dp
1121 ros_E(2) = 0.0E+00_dp
1122 ros_E(3) = 0.0E+00_dp
1123 ros_E(4) = 1.0E+00_dp
1124 !~~~> ros_ELO = estimator of local order - the minimum between the
1125 ! main and the embedded scheme orders plus 1
1126 ros_ELO = 3.0E+00_dp
1127 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
1128 ros_Alpha(1) = 0.0E+00_dp
1129 ros_Alpha(2) = 0.0E+00_dp
1130 ros_Alpha(3) = 1.0E+00_dp
1131 ros_Alpha(4) = 1.0E+00_dp
1132 !~~~> Gamma_i = \sum_j gamma_{i,j}
1133 ros_Gamma(1) = 0.5E+00_dp
1134 ros_Gamma(2) = 1.5E+00_dp
1135 ros_Gamma(3) = 0.0E+00_dp
1136 ros_Gamma(4) = 0.0E+00_dp
1138 END SUBROUTINE KPP_ROOT_Rodas3
1140 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1141 SUBROUTINE KPP_ROOT_Rodas4 (ros_S,ros_A,ros_C,ros_M,ros_E,ros_Alpha,&
1142 ros_Gamma,ros_NewF,ros_ELO,ros_Name)
1143 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1144 ! STIFFLY-STABLE ROSENBROCK METHOD OF ORDER 4, WITH 6 STAGES
1146 ! E. HAIRER AND G. WANNER, SOLVING ORDINARY DIFFERENTIAL
1147 ! EQUATIONS II. STIFF AND DIFFERENTIAL-ALGEBRAIC PROBLEMS.
1148 ! SPRINGER SERIES IN COMPUTATIONAL MATHEMATICS,
1149 ! SPRINGER-VERLAG (1996)
1150 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1152 IMPLICIT NONE
1154 INTEGER, PARAMETER :: S=6
1155 INTEGER, INTENT(OUT) :: ros_S
1156 KPP_REAL, DIMENSION(S), INTENT(OUT) :: ros_M,ros_E,ros_Alpha,ros_Gamma
1157 KPP_REAL, DIMENSION(S*(S-1)/2), INTENT(OUT) :: ros_A, ros_C
1158 KPP_REAL, INTENT(OUT) :: ros_ELO
1159 LOGICAL, DIMENSION(S), INTENT(OUT) :: ros_NewF
1160 CHARACTER(LEN=12), INTENT(OUT) :: ros_Name
1161 !cms DOUBLE PRECISION g
1162 KPP_REAL :: g
1164 !~~~> Name of the method
1165 ros_Name = 'RODAS-4'
1166 !~~~> Number of stages
1167 ros_S = 6
1169 !~~~> Y_stage_i ~ Y( T + H*Alpha_i )
1170 ros_Alpha(1) = 0.000_dp
1171 ros_Alpha(2) = 0.386_dp
1172 ros_Alpha(3) = 0.210_dp
1173 ros_Alpha(4) = 0.630_dp
1174 ros_Alpha(5) = 1.000_dp
1175 ros_Alpha(6) = 1.000_dp
1177 !~~~> Gamma_i = \sum_j gamma_{i,j}
1178 ros_Gamma(1) = 0.2500000000000000E+00_dp
1179 ros_Gamma(2) =-0.1043000000000000E+00_dp
1180 ros_Gamma(3) = 0.1035000000000000E+00_dp
1181 ros_Gamma(4) =-0.3620000000000023E-01_dp
1182 ros_Gamma(5) = 0.0_dp
1183 ros_Gamma(6) = 0.0_dp
1185 !~~~> The coefficient matrices A and C are strictly lower triangular.
1186 ! The lower triangular (subdiagonal) elements are stored in row-wise order:
1187 ! A(2,1) = ros_A(1), A(3,1)=ros_A(2), A(3,2)=ros_A(3), etc.
1188 ! The general mapping formula is: A(i,j) = ros_A( (i-1)*(i-2)/2 + j )
1189 ! C(i,j) = ros_C( (i-1)*(i-2)/2 + j )
1191 ros_A(1) = 0.1544000000000000E+01_dp
1192 ros_A(2) = 0.9466785280815826E+00_dp
1193 ros_A(3) = 0.2557011698983284E+00_dp
1194 ros_A(4) = 0.3314825187068521E+01_dp
1195 ros_A(5) = 0.2896124015972201E+01_dp
1196 ros_A(6) = 0.9986419139977817E+00_dp
1197 ros_A(7) = 0.1221224509226641E+01_dp
1198 ros_A(8) = 0.6019134481288629E+01_dp
1199 ros_A(9) = 0.1253708332932087E+02_dp
1200 ros_A(10) =-0.6878860361058950E+00_dp
1201 ros_A(11) = ros_A(7)
1202 ros_A(12) = ros_A(8)
1203 ros_A(13) = ros_A(9)
1204 ros_A(14) = ros_A(10)
1205 ros_A(15) = 1.0E+00_dp
1207 ros_C(1) =-0.5668800000000000E+01_dp
1208 ros_C(2) =-0.2430093356833875E+01_dp
1209 ros_C(3) =-0.2063599157091915E+00_dp
1210 ros_C(4) =-0.1073529058151375E+00_dp
1211 ros_C(5) =-0.9594562251023355E+01_dp
1212 ros_C(6) =-0.2047028614809616E+02_dp
1213 ros_C(7) = 0.7496443313967647E+01_dp
1214 ros_C(8) =-0.1024680431464352E+02_dp
1215 ros_C(9) =-0.3399990352819905E+02_dp
1216 ros_C(10) = 0.1170890893206160E+02_dp
1217 ros_C(11) = 0.8083246795921522E+01_dp
1218 ros_C(12) =-0.7981132988064893E+01_dp
1219 ros_C(13) =-0.3152159432874371E+02_dp
1220 ros_C(14) = 0.1631930543123136E+02_dp
1221 ros_C(15) =-0.6058818238834054E+01_dp
1223 !~~~> M_i = Coefficients for new step solution
1224 ros_M(1) = ros_A(7)
1225 ros_M(2) = ros_A(8)
1226 ros_M(3) = ros_A(9)
1227 ros_M(4) = ros_A(10)
1228 ros_M(5) = 1.0E+00_dp
1229 ros_M(6) = 1.0E+00_dp
1231 !~~~> E_i = Coefficients for error estimator
1232 ros_E(1) = 0.0E+00_dp
1233 ros_E(2) = 0.0E+00_dp
1234 ros_E(3) = 0.0E+00_dp
1235 ros_E(4) = 0.0E+00_dp
1236 ros_E(5) = 0.0E+00_dp
1237 ros_E(6) = 1.0E+00_dp
1239 !~~~> Does the stage i require a new function evaluation (ros_NewF(i)=TRUE)
1240 ! or does it re-use the function evaluation from stage i-1 (ros_NewF(i)=FALSE)
1241 ros_NewF(1) = .TRUE.
1242 ros_NewF(2) = .TRUE.
1243 ros_NewF(3) = .TRUE.
1244 ros_NewF(4) = .TRUE.
1245 ros_NewF(5) = .TRUE.
1246 ros_NewF(6) = .TRUE.
1248 !~~~> ros_ELO = estimator of local order - the minimum between the
1249 ! main and the embedded scheme orders plus 1
1250 ros_ELO = 4.0_dp
1252 END SUBROUTINE KPP_ROOT_Rodas4
1254 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1255 ! End of the set of internal Rosenbrock subroutines
1256 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1257 END SUBROUTINE KPP_ROOT_Rosenbrock
1258 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1261 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1262 SUBROUTINE KPP_ROOT_FunTemplate( T, Y, Ydot, RCONST, FIX, Nfun )
1263 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1264 ! Template for the ODE function call.
1265 ! Updates the rate coefficients (and possibly the fixed species) at each call
1266 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1267 USE KPP_ROOT_Parameters
1268 !! USE KPP_ROOT_Global
1269 !! USE KPP_ROOT_Function
1270 !! USE KPP_ROOT_Rates
1271 !~~~> Input variables
1272 KPP_REAL :: T, Y(NVAR)
1273 KPP_REAL :: RCONST(NREACT)
1274 KPP_REAL :: FIX(NFIX)
1275 !~~~> Output variables
1276 KPP_REAL :: Ydot(NVAR)
1277 INTEGER :: Nfun
1280 !~~~> Local variables
1281 !! KPP_REAL :: Told
1283 !! Told = TIME
1284 !! TIME = T
1285 !! CALL Update_SUN()
1286 !! CALL Update_RCONST()
1287 CALL KPP_ROOT_Fun( Y, FIX, RCONST, Ydot )
1288 !! TIME = Told
1290 Nfun = Nfun+1
1292 END SUBROUTINE KPP_ROOT_FunTemplate
1295 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1296 SUBROUTINE KPP_ROOT_JacTemplate( T, Y, Jcb, FIX, Njac, RCONST )
1297 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1298 ! Template for the ODE Jacobian call.
1299 ! Updates the rate coefficients (and possibly the fixed species) at each call
1300 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1301 USE KPP_ROOT_Parameters
1302 !!USE KPP_ROOT_Global
1303 USE KPP_ROOT_Jacobian
1304 !! USE KPP_ROOT_LinearAlgebra
1305 !! USE KPP_ROOT_Rates
1306 !~~~> Input variables
1307 KPP_REAL :: T, Y(NVAR)
1308 KPP_REAL :: FIX(NFIX)
1309 KPP_REAL :: RCONST(NREACT)
1311 INTEGER :: Njac
1313 !~~~> Output variables
1314 #ifdef FULL_ALGEBRA
1315 KPP_REAL :: JV(LU_NONZERO), Jcb(NVAR,NVAR)
1316 #else
1317 KPP_REAL :: Jcb(LU_NONZERO)
1318 #endif
1319 !~~~> Local variables
1320 KPP_REAL :: Told
1321 #ifdef FULL_ALGEBRA
1322 INTEGER :: i, j
1323 #endif
1325 !! Told = TIME
1326 !! TIME = T
1327 !! CALL Update_SUN()
1328 !! CALL Update_RCONST()
1329 #ifdef FULL_ALGEBRA
1330 CALL KPP_ROOT_Jac_SP(Y, FIX, RCONST, JV)
1331 DO j=1,NVAR
1332 DO i=1,NVAR
1333 Jcb(i,j) = 0.0d0
1334 END DO
1335 END DO
1336 DO i=1,LU_NONZERO
1337 Jcb(LU_IROW(i),LU_ICOL(i)) = JV(i)
1338 END DO
1339 #else
1340 CALL KPP_ROOT_Jac_SP( Y, FIX, RCONST, Jcb )
1341 #endif
1342 !! TIME = Told
1344 Njac = Njac+1
1346 END SUBROUTINE KPP_ROOT_JacTemplate
1348 !!!END MODULE KPP_ROOT_Integrator