Line data Source code
1 : !---------------------------------------------------------------------
2 : ! ... compute chemical potential heating
3 : !---------------------------------------------------------------------
4 :
5 : module mo_cph
6 :
7 : use shr_kind_mod,only : r8 => shr_kind_r8
8 : use chem_mods, only : ncph=>enthalpy_cnt, exotherm=>cph_enthalpy, cph_rid
9 :
10 : implicit none
11 :
12 : save
13 :
14 : !==============================================================
15 : !... Doug Kinnison, dkin@ucar.edu
16 : !
17 : !... Enthalpy Data are taken from Atkinson et al.,
18 : ! Evaluated kinetic and photochemical data for atmospheric
19 : ! chemistry: Volume I, Atmos. Chem. Phys., 4, 1461-1738.
20 : !... Heats of formation at 0K.
21 : !... Units: kJ mol-1
22 : !
23 : !... Exception to the Atkinson et al. reference (@0K unless noted)
24 : ! (4), (5), (8), (9), (10, (11), (14), (15), (27), (28) at 298K
25 : ! (7) h + o2 -> oh + o2 is multiplied by 0.6 (Mlynczak) to represent
26 : ! AG loss of excited OH.
27 : ! (25) n2d + o2 -> no + o1d taken from Roble, UMLT, Johnson and Killeen.
28 : ! (26) n2d + o -> n + o taken from Roble, UMLT, Johnson and Killeen.
29 : ! (30-41) Taken from Roble, UMLT, Johnson and Killeen Ed., Geophys. Mono. 87
30 : !==============================================================
31 : !... Enthalpy Data are specified in preprocessor input mechanism file
32 : !... F Vitt -- 29 Oct 2013
33 : !==============================================================
34 :
35 : private
36 : public :: cph, init_cph
37 :
38 : logical :: has_cph
39 : character(len=24) :: fldnames(ncph)
40 : logical, parameter :: debug = .false.
41 :
42 : contains
43 :
44 : !-------------------------------------------------------------------------
45 : !-------------------------------------------------------------------------
46 1536 : subroutine init_cph
47 :
48 : use mo_chem_utls, only : get_rxt_ndx, get_spc_ndx
49 : use cam_history, only : addfld, add_default
50 : use chem_mods, only : rxt_tag_lst, rxt_tag_map, rxt_tag_cnt
51 : use cam_abortutils, only : endrun
52 :
53 : implicit none
54 :
55 : character(len=64) :: longname
56 : integer :: i, n, tagndx
57 :
58 1536 : has_cph = ncph > 0
59 :
60 : if (.not.has_cph) return
61 :
62 : if ( any(exotherm(:) == 0._r8) ) then
63 : call endrun('init_cph: Enthalpies for chemical heating must be specified in mechanism file')
64 : endif
65 :
66 : do i = 1,ncph
67 :
68 : findtagndx: do n = 1,rxt_tag_cnt
69 : if ( rxt_tag_map(n) == cph_rid(i) ) then
70 : tagndx = n
71 : exit findtagndx
72 : endif
73 : enddo findtagndx
74 :
75 : if (debug) then
76 : if ( i< 10 ) then
77 : write(fldnames(i),fmt='("CPH",i1)') i
78 : else if (i<100) then
79 : write(fldnames(i),fmt='("CPH",i2)') i
80 : else if (i<1000) then
81 : write(fldnames(i),fmt='("CPH",i3)') i
82 : endif
83 : else
84 : fldnames(i) = 'CPH_'//trim(rxt_tag_lst(tagndx))
85 : endif
86 :
87 : write( longname, fmt='(f12.6)') exotherm(i)
88 : longname = trim(adjustl(longname))//' kcal/mol chem pot heating rate for rxtn '//trim(rxt_tag_lst(tagndx))
89 : call addfld( fldnames(i), (/ 'lev' /), 'I', 'K/s', trim(longname) )
90 : if (debug) then
91 : call add_default( fldnames(i), 10, ' ' )
92 : endif
93 : enddo
94 :
95 : call addfld( 'QCP', (/ 'lev' /), 'I', 'K/s', 'chem pot heating rate' )
96 :
97 1536 : end subroutine init_cph
98 :
99 : !-------------------------------------------------------------------------
100 : !-------------------------------------------------------------------------
101 0 : subroutine cph( cph_tot, vmr, rxt, cp, mbar, kbot, ncol, lchnk )
102 :
103 : !-----------------------------------------------------------------------
104 : ! ... forms the chemical potential heating rates
105 : !-----------------------------------------------------------------------
106 :
107 1536 : use chem_mods, only : gas_pcnst, rxntot
108 : use ppgrid, only : pver
109 : use cam_history, only : outfld
110 : use mo_rxt_rates_conv, only : set_rates
111 :
112 : implicit none
113 :
114 : !-----------------------------------------------------------------------
115 : ! ... dummy arguments
116 : !-----------------------------------------------------------------------
117 : integer, intent(in) :: ncol ! columns in chunck
118 : integer, intent(in) :: lchnk ! chunk index
119 : integer, intent(in) :: kbot ! bottom vert index
120 : real(r8), intent(in) :: rxt(ncol,pver,rxntot) ! rxt rates (1/cm^3/s)
121 : real(r8), intent(in) :: cp(ncol,pver) ! specific heat capacity (J/K/kg)
122 : real(r8), intent(in) :: mbar(ncol,pver) ! atm mean mass (g/mole)
123 : real(r8), intent(in) :: vmr(ncol,pver,gas_pcnst) ! concentrations (mol/mol)
124 : real(r8), intent(out) :: cph_tot(ncol,pver) ! total heating (K/s)
125 :
126 : !-----------------------------------------------------------------------
127 : ! ... local variables
128 : !-----------------------------------------------------------------------
129 : integer :: i, k
130 0 : real(r8) :: tmp(ncol,pver)
131 0 : real(r8) :: cph_rate(ncol,pver,ncph)
132 0 : real(r8) :: rxt_rates(ncol,pver,rxntot)
133 :
134 0 : if (.not.has_cph) return
135 :
136 : ! get the reaction rates from rate constants ...
137 0 : rxt_rates(:ncol,:,:) = rxt(:ncol,:,:)
138 0 : call set_rates( rxt_rates, vmr, ncol )
139 :
140 : ! compute corresponding chem heating rates ...
141 : cph_rate(:,:,:) = 0._r8
142 0 : tmp(:ncol,:) = 1._r8 / (1.e-6_r8*cp(:ncol,:)*mbar(:ncol,:))
143 : do i = 1,ncph
144 : cph_rate(:ncol,:,i) = tmp(:ncol,:) * rxt_rates(:ncol,:,cph_rid(i)) * exotherm(i)
145 : enddo
146 :
147 : ! compute total heating rate ...
148 0 : cph_tot(:,:) = 0._r8
149 0 : do k = 1,kbot
150 0 : do i = 1,ncol
151 0 : cph_tot(i,k) = sum( cph_rate(i,k,:) )
152 : end do
153 : end do
154 :
155 : ! output diagnostics
156 : do i = 1,ncph
157 : if ( exotherm(i)>0._r8) then
158 : call outfld( fldnames(i), cph_rate(:,:,i), ncol, lchnk )
159 : endif
160 : enddo
161 :
162 0 : call outfld( 'QCP', cph_tot(:,:), ncol, lchnk )
163 :
164 0 : end subroutine cph
165 :
166 : end module mo_cph
|