Line data Source code
1 : module nlte_aliarms
2 :
3 : !
4 : ! provides calculation of non-LTE heating rates by ALI-ARMS non-LTE code
5 : !
6 : use ppgrid, only: pcols, pver
7 : use shr_kind_mod, only: r8 => shr_kind_r8
8 : use cam_logfile, only: iulog
9 : use spmd_utils, only: masterproc
10 : use iso_c_binding, only: c_float, c_int
11 :
12 : implicit none
13 : private
14 : save
15 :
16 : ! Public interfaces
17 : public :: nlte_aliarms_init
18 : public :: nlte_aliarms_calc
19 :
20 :
21 : integer(c_int) :: pver_c = -1 ! pver for the ALI_ARMS C code (limited by max_pressure_lw)
22 :
23 : real(r8) :: o1_mw_inv = -huge(1.0_r8) ! O molecular weight (inverse)
24 : real(r8) :: o2_mw_inv = -huge(1.0_r8) ! O2 molecular weight (inverse)
25 : real(r8) :: co2_mw_inv = -huge(1.0_r8) ! CO2 molecular weight (inverse)
26 : real(r8) :: n2_mw_inv = -huge(1.0_r8) ! N2 molecular weight (inverse)
27 : contains
28 :
29 : !-----------------------------------------------------------------
30 0 : subroutine nlte_aliarms_init(max_pressure_lw,co2_mw,n2_mw,o1_mw,o2_mw)
31 : !-----------------------------------------------------------------
32 : !
33 : !
34 : !-----------------------------------------------------------------
35 :
36 : use cam_history, only: addfld
37 : use ref_pres, only: pref_mid
38 :
39 : real(r8), intent(in) :: max_pressure_lw ! Pa
40 : real(r8), intent(in) :: o1_mw ! O molecular weight
41 : real(r8), intent(in) :: o2_mw ! O2 molecular weight
42 : real(r8), intent(in) :: co2_mw ! CO2 molecular weight
43 : real(r8), intent(in) :: n2_mw ! N2 molecular weight
44 :
45 : integer :: iver
46 :
47 0 : if (masterproc) then
48 0 : write(iulog,*) 'init: ALI-ARMS non-LTE code'
49 : end if
50 :
51 0 : call addfld ('ALIARMS_Q',(/ 'lev' /), 'A','K/s','Non-LTE LW CO2 heating rate')
52 :
53 0 : co2_mw_inv = 1._r8/co2_mw
54 0 : o1_mw_inv = 1._r8/o1_mw
55 0 : o2_mw_inv = 1._r8/o2_mw
56 0 : n2_mw_inv = 1._r8/n2_mw
57 :
58 0 : pver_c=0
59 0 : do iver = 1,pver
60 0 : if (pref_mid(iver) < max_pressure_lw) then
61 0 : pver_c=pver_c+1
62 : else
63 : exit ! Have gone past the maximum pressure
64 : end if
65 : end do
66 :
67 0 : end subroutine nlte_aliarms_init
68 :
69 : !-----------------------------------------------------------------
70 0 : subroutine nlte_aliarms_calc (lchnk,ncol,state_zm,pmid,t,xo2mmr,xommr,xn2mmr,xco2mmr,cool)
71 : !-----------------------------------------------------------------
72 : !
73 : !
74 : !-----------------------------------------------------------------
75 :
76 0 : use air_composition, only: mbarv
77 : use cam_history, only: outfld
78 : use shr_infnan_mod, only: is_nan => shr_infnan_isnan
79 : use shr_kind_mod, only: SHR_KIND_CM
80 : use cam_abortutils, only: endrun
81 :
82 : ! Input variables
83 : integer, intent(in) :: ncol ! number of atmospheric columns
84 : integer, intent(in) :: lchnk ! chunk identifier
85 :
86 : real(r8), intent(in) :: state_zm(pcols,pver) ! model height (m)
87 : real(r8), intent(in) :: pmid(pcols,pver) ! model pressure at mid-point (Pa)
88 : real(r8), intent(in) :: t(pcols,pver) ! Neutral temperature (K)
89 :
90 : real(r8), intent(in) :: xco2mmr(pcols,pver) ! CO2 mass mixing ratio profile
91 : real(r8), intent(in) :: xn2mmr(pcols,pver) ! N2 mass mixing ratio profile
92 : real(r8), intent(in) :: xommr(pcols,pver) ! O mass mixing ratio profile
93 : real(r8), intent(in) :: xo2mmr(pcols,pver) ! O2 mass mixing ratio profile
94 :
95 : ! Output variables
96 : real(r8), intent(out) :: cool(pcols,pver) ! CO2 NLTE cooling rate (K/s)
97 :
98 : ! local variables
99 :
100 0 : real(c_float), dimension(pver_c) :: p, tn, zkm
101 0 : real(c_float), dimension(pver_c) :: co2_vmr, o_vmr, n2_vmr, o2_vmr
102 0 : real(c_float), dimension(pver_c) :: ali_cool
103 :
104 : integer :: icol, iver, i, j
105 :
106 : character (len=SHR_KIND_CM) :: errstring
107 :
108 : ! Interface to ali C routine
109 : ! Note that ali uses single precision C floats, so conversions from r8 to c_float are made before calling ali
110 : interface
111 : subroutine ali_(zkm, p, tn, co2_vmr, o_vmr, n2_vmr, o2_vmr, ali_cool, pver_c) bind(c,name='ali_')
112 : use iso_c_binding, only: c_float, c_int
113 : real(c_float), dimension(*) :: p, tn, zkm ! (in) input pressure(Pa), temperature(K) and height (km)
114 : real(c_float), dimension(*) :: co2_vmr, o_vmr, n2_vmr, o2_vmr ! (in) volume mixing ratios
115 : real(c_float), dimension(*) :: ali_cool ! (out) cooling rate (K/s)
116 : integer(c_int) :: pver_c
117 : end subroutine ali_
118 : end interface
119 :
120 0 : cool(:,:) = 0.0_r8
121 :
122 0 : do icol=1,ncol
123 :
124 0 : ali_cool(:) = 0.0_c_float
125 0 : zkm(:) = 0.0_c_float
126 0 : tn(:) = 0.0_c_float
127 0 : co2_vmr(:) = 0.0_c_float
128 0 : o_vmr(:) = 0.0_c_float
129 0 : n2_vmr(:) = 0.0_c_float
130 0 : o2_vmr(:) = 0.0_c_float
131 :
132 0 : p(:pver_c) = pmid(icol,:pver_c)*1.0e-5_c_float ! convert pmid in Pa to bars
133 0 : zkm(:pver_c) = state_zm(icol,:pver_c)*1.e-3_c_float
134 0 : tn(:pver_c) = t(icol,:pver_c)
135 :
136 : ! Convert to VMR from mmr
137 0 : co2_vmr(:pver_c) = mbarv(icol,:pver_c ,lchnk) * xco2mmr(icol,:pver_c) * co2_mw_inv
138 0 : o_vmr(:pver_c) = mbarv(icol,:pver_c ,lchnk) * xommr(icol,:pver_c) * o1_mw_inv
139 0 : n2_vmr(:pver_c) = mbarv(icol,:pver_c ,lchnk) * xn2mmr(icol,:pver_c) * n2_mw_inv
140 0 : o2_vmr(:pver_c) = mbarv(icol,:pver_c ,lchnk) * xo2mmr(icol,:pver_c) * o2_mw_inv
141 :
142 0 : call ali(zkm, p, tn, co2_vmr, o_vmr, n2_vmr, o2_vmr, ali_cool, pver_c)
143 :
144 0 : cool(icol,:pver_c) = ali_cool(:pver_c)
145 :
146 : enddo
147 :
148 : ! Check the rates
149 0 : do j=1,pver
150 0 : do i=1,ncol
151 0 : if (is_nan(cool(i,j))) then
152 0 : write(errstring,*) 'nlte_aliarms_calc: Nan in qrlaliarms for chunk=', lchnk, ' column index=',i,' vertical index=',j
153 0 : call endrun (errstring)
154 : end if
155 : end do
156 : end do
157 :
158 : ! Check cool for any way out-of-bounds values
159 0 : if (any(cool(:ncol,:pver_c) > 0.2_r8) .or. any(cool(:ncol,:pver_c)<-1.0_r8)) then
160 0 : write(errstring,*) 'nlte_aliarms_calc: Cooling rate (cool) is greater than .2 or less than -1 K/s for chunk ', lchnk
161 0 : call endrun (errstring)
162 : end if
163 :
164 0 : call outfld ('ALIARMS_Q', cool, pcols, lchnk)
165 :
166 0 : end subroutine nlte_aliarms_calc
167 :
168 : end module nlte_aliarms
169 :
|