Line data Source code
1 :
2 : module mo_setrxt
3 :
4 : use shr_kind_mod, only : r8 => shr_kind_r8
5 :
6 : private
7 : public :: setrxt
8 : public :: setrxt_hrates
9 :
10 : contains
11 :
12 1489176 : subroutine setrxt( rate, temp, m, ncol )
13 :
14 : use ppgrid, only : pver, pcols
15 : use shr_kind_mod, only : r8 => shr_kind_r8
16 : use chem_mods, only : rxntot
17 : use mo_jpl, only : jpl
18 :
19 : implicit none
20 :
21 : !-------------------------------------------------------
22 : ! ... dummy arguments
23 : !-------------------------------------------------------
24 : integer, intent(in) :: ncol
25 : real(r8), intent(in) :: temp(pcols,pver)
26 : real(r8), intent(in) :: m(ncol,pver)
27 : real(r8), intent(inout) :: rate(ncol,pver,rxntot)
28 :
29 : !-------------------------------------------------------
30 : ! ... local variables
31 : !-------------------------------------------------------
32 : integer :: n
33 2978352 : real(r8) :: itemp(ncol,pver)
34 : real(r8) :: exp_fac(ncol,pver)
35 2978352 : real(r8) :: ko(ncol,pver)
36 2978352 : real(r8) :: kinf(ncol,pver)
37 :
38 2314006344 : rate(:,:,5) = 1.8e-12_r8
39 2314006344 : rate(:,:,15) = 1.157e-05_r8
40 2314006344 : itemp(:ncol,:) = 1._r8 / temp(:ncol,:)
41 1489176 : n = ncol*pver
42 2314006344 : rate(:,:,11) = 1.9e-13_r8 * exp( 520._r8 * itemp(:,:) )
43 2314006344 : rate(:,:,12) = 1.1e-11_r8 * exp( -280._r8 * itemp(:,:) )
44 :
45 2314006344 : itemp(:,:) = 300._r8 * itemp(:,:)
46 :
47 2314006344 : ko(:,:) = 2.9e-31_r8 * itemp(:,:)**4.1_r8
48 2314006344 : kinf(:,:) = 1.7e-12_r8 * itemp(:,:)**(-0.2_r8)
49 1489176 : call jpl( rate(1,1,13), m, 0.6_r8, ko, kinf, n )
50 :
51 1489176 : end subroutine setrxt
52 :
53 :
54 0 : subroutine setrxt_hrates( rate, temp, m, ncol, kbot )
55 :
56 : use ppgrid, only : pver, pcols
57 : use shr_kind_mod, only : r8 => shr_kind_r8
58 : use chem_mods, only : rxntot
59 : use mo_jpl, only : jpl
60 :
61 : implicit none
62 :
63 : !-------------------------------------------------------
64 : ! ... dummy arguments
65 : !-------------------------------------------------------
66 : integer, intent(in) :: ncol
67 : integer, intent(in) :: kbot
68 : real(r8), intent(in) :: temp(pcols,pver)
69 : real(r8), intent(in) :: m(ncol,pver)
70 : real(r8), intent(inout) :: rate(ncol,pver,rxntot)
71 :
72 : !-------------------------------------------------------
73 : ! ... local variables
74 : !-------------------------------------------------------
75 : integer :: n
76 : real(r8) :: itemp(ncol,kbot)
77 : real(r8) :: exp_fac(ncol,kbot)
78 : real(r8) :: ko(ncol,kbot)
79 : real(r8) :: kinf(ncol,kbot)
80 : real(r8) :: wrk(ncol,kbot)
81 :
82 :
83 0 : end subroutine setrxt_hrates
84 :
85 : end module mo_setrxt
|