00001 module seq_drydep_mod
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019 use shr_sys_mod, only : shr_sys_abort
00020 use shr_log_mod, only : s_loglev => shr_log_Level
00021 use shr_kind_mod, only : r8 => shr_kind_r8, CS => SHR_KIND_CS, CX => SHR_KIND_CX
00022 use shr_const_mod, only : SHR_CONST_G, SHR_CONST_RDAIR, &
00023 SHR_CONST_CPDAIR, SHR_CONST_MWWV
00024
00025 implicit none
00026 save
00027
00028 private
00029
00030
00031
00032 public :: seq_drydep_read
00033 public :: seq_drydep_init
00034 public :: seq_drydep_setHCoeff
00035
00036
00037
00038 integer, private, parameter :: maxspc = 100
00039 integer, private, parameter :: n_species_table = 53
00040 integer, private, parameter :: NSeas = 5
00041 integer, private, parameter :: NLUse = 11
00042
00043
00044
00045
00046 character(16),public,parameter :: DD_XATM = 'xactive_atm'
00047 character(16),public,parameter :: DD_XLND = 'xactive_lnd'
00048 character(16),public,parameter :: DD_TABL = 'table'
00049 character(16),public :: drydep_method = DD_XLND
00050
00051 real(r8), public, parameter :: ph = 1.e-5_r8
00052
00053 logical, public :: lnd_drydep
00054 integer, public :: n_drydep = 0
00055 character(len=32), public, dimension(maxspc) :: drydep_list = ''
00056
00057 character(len=CS), public :: drydep_fields_token = ''
00058
00059 real(r8), public, allocatable, dimension(:) :: foxd
00060 real(r8), public, allocatable, dimension(:) :: drat
00061 integer, public, allocatable, dimension(:) :: mapping
00062
00063 integer, public :: h2_ndx, ch4_ndx, co_ndx, pan_ndx, mpan_ndx, so2_ndx, o3_ndx, o3a_ndx
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099 real(r8), parameter, public :: h2_a(NLUse) =
00100 (/ 0.000_r8, 0.000_r8, 0.270_r8, 0.000_r8, 0.000_r8,
00101 0.000_r8, 0.000_r8, 0.000_r8, 0.000_r8, 0.000_r8, 0.000_r8/)
00102
00103 real(r8), parameter, public :: h2_b(NLUse) =
00104 (/ 0.000_r8,-41.390_r8, -0.472_r8,-41.900_r8,-41.900_r8,
00105 -41.900_r8, 0.000_r8, 0.000_r8, 0.000_r8,-41.390_r8, 0.000_r8/)
00106
00107 real(r8), parameter, public :: h2_c(NLUse) =
00108 (/ 0.000_r8, 16.850_r8, 1.235_r8, 19.700_r8, 19.700_r8,
00109 19.700_r8, 0.000_r8, 0.000_r8, 0.000_r8, 17.700_r8, 1.000_r8/)
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121 real(r8), public, dimension(NSeas,NLUse) :: ri, rlu, rac, rgss, rgso, rcls, rclo
00122
00123 data ri (1,1:NLUse) &
00124 /1.e36_r8, 60._r8, 120._r8, 70._r8, 130._r8, 100._r8,1.e36_r8,1.e36_r8, 80._r8, 100._r8, 150._r8/
00125 data rlu (1,1:NLUse) &
00126 /1.e36_r8,2000._r8,2000._r8,2000._r8,2000._r8,2000._r8,1.e36_r8,1.e36_r8,2500._r8,2000._r8,4000._r8/
00127 data rac (1,1:NLUse) &
00128 / 100._r8, 200._r8, 100._r8,2000._r8,2000._r8,2000._r8, 0._r8, 0._r8, 300._r8, 150._r8, 200._r8/
00129 data rgss(1,1:NLUse) &
00130 / 400._r8, 150._r8, 350._r8, 500._r8, 500._r8, 100._r8, 0._r8,1000._r8, 0._r8, 220._r8, 400._r8/
00131 data rgso(1,1:NLUse) &
00132 / 300._r8, 150._r8, 200._r8, 200._r8, 200._r8, 300._r8,2000._r8, 400._r8,1000._r8, 180._r8, 200._r8/
00133 data rcls(1,1:NLUse) &
00134 /1.e36_r8,2000._r8,2000._r8,2000._r8,2000._r8,2000._r8,1.e36_r8,1.e36_r8,2500._r8,2000._r8,4000._r8/
00135 data rclo(1,1:NLUse) &
00136 /1.e36_r8,1000._r8,1000._r8,1000._r8,1000._r8,1000._r8,1.e36_r8,1.e36_r8,1000._r8,1000._r8,1000._r8/
00137
00138 data ri (2,1:NLUse) &
00139 /1.e36_r8,1.e36_r8,1.e36_r8,1.e36_r8, 250._r8, 500._r8,1.e36_r8,1.e36_r8,1.e36_r8,1.e36_r8,1.e36_r8/
00140 data rlu (2,1:NLUse) &
00141 /1.e36_r8,9000._r8,9000._r8,9000._r8,4000._r8,8000._r8,1.e36_r8,1.e36_r8,9000._r8,9000._r8,9000._r8/
00142 data rac (2,1:NLUse) &
00143 / 100._r8, 150._r8, 100._r8,1500._r8,2000._r8,1700._r8, 0._r8, 0._r8, 200._r8, 120._r8, 140._r8/
00144 data rgss(2,1:NLUse) &
00145 / 400._r8, 200._r8, 350._r8, 500._r8, 500._r8, 100._r8, 0._r8,1000._r8, 0._r8, 300._r8, 400._r8/
00146 data rgso(2,1:NLUse) &
00147 / 300._r8, 150._r8, 200._r8, 200._r8, 200._r8, 300._r8,2000._r8, 400._r8, 800._r8, 180._r8, 200._r8/
00148 data rcls(2,1:NLUse) &
00149 /1.e36_r8,9000._r8,9000._r8,9000._r8,2000._r8,4000._r8,1.e36_r8,1.e36_r8,9000._r8,9000._r8,9000._r8/
00150 data rclo(2,1:NLUse) &
00151 /1.e36_r8, 400._r8, 400._r8, 400._r8,1000._r8, 600._r8,1.e36_r8,1.e36_r8, 400._r8, 400._r8, 400._r8/
00152
00153 data ri (3,1:NLUse) &
00154 /1.e36_r8,1.e36_r8,1.e36_r8,1.e36_r8, 250._r8, 500._r8,1.e36_r8,1.e36_r8,1.e36_r8,1.e36_r8,1.e36_r8/
00155 data rlu (3,1:NLUse) &
00156 /1.e36_r8,1.e36_r8,9000._r8,9000._r8,4000._r8,8000._r8,1.e36_r8,1.e36_r8,9000._r8,9000._r8,9000._r8/
00157 data rac (3,1:NLUse) &
00158 / 100._r8, 10._r8, 100._r8,1000._r8,2000._r8,1500._r8, 0._r8, 0._r8, 100._r8, 50._r8, 120._r8/
00159 data rgss(3,1:NLUse) &
00160 / 400._r8, 150._r8, 350._r8, 500._r8, 500._r8, 200._r8, 0._r8,1000._r8, 0._r8, 200._r8, 400._r8/
00161 data rgso(3,1:NLUse) &
00162 / 300._r8, 150._r8, 200._r8, 200._r8, 200._r8, 300._r8,2000._r8, 400._r8,1000._r8, 180._r8, 200._r8/
00163 data rcls(3,1:NLUse) &
00164 /1.e36_r8,1.e36_r8,9000._r8,9000._r8,3000._r8,6000._r8,1.e36_r8,1.e36_r8,9000._r8,9000._r8,9000._r8/
00165 data rclo(3,1:NLUse) &
00166 /1.e36_r8,1000._r8, 400._r8, 400._r8,1000._r8, 600._r8,1.e36_r8,1.e36_r8, 800._r8, 600._r8, 600._r8/
00167
00168 data ri (4,1:NLUse) &
00169 /1.e36_r8,1.e36_r8,1.e36_r8,1.e36_r8, 400._r8, 800._r8,1.e36_r8,1.e36_r8,1.e36_r8,1.e36_r8,1.e36_r8/
00170 data rlu (4,1:NLUse) &
00171 /1.e36_r8,1.e36_r8,1.e36_r8,1.e36_r8,6000._r8,9000._r8,1.e36_r8,1.e36_r8,9000._r8,9000._r8,9000._r8/
00172 data rac (4,1:NLUse) &
00173 / 100._r8, 10._r8, 10._r8,1000._r8,2000._r8,1500._r8, 0._r8, 0._r8, 50._r8, 10._r8, 50._r8/
00174 data rgss(4,1:NLUse) &
00175 / 100._r8, 100._r8, 100._r8, 100._r8, 100._r8, 100._r8, 0._r8,1000._r8, 100._r8, 100._r8, 50._r8/
00176 data rgso(4,1:NLUse) &
00177 / 600._r8,3500._r8,3500._r8,3500._r8,3500._r8,3500._r8,2000._r8, 400._r8,3500._r8,3500._r8,3500._r8/
00178 data rcls(4,1:NLUse) &
00179 /1.e36_r8,1.e36_r8,1.e36_r8,9000._r8, 200._r8, 400._r8,1.e36_r8,1.e36_r8,9000._r8,1.e36_r8,9000._r8/
00180 data rclo(4,1:NLUse) &
00181 /1.e36_r8,1000._r8,1000._r8, 400._r8,1500._r8, 600._r8,1.e36_r8,1.e36_r8, 800._r8,1000._r8, 800._r8/
00182
00183 data ri (5,1:NLUse) &
00184 /1.e36_r8, 120._r8, 240._r8, 140._r8, 250._r8, 190._r8,1.e36_r8,1.e36_r8, 160._r8, 200._r8, 300._r8/
00185 data rlu (5,1:NLUse) &
00186 /1.e36_r8,4000._r8,4000._r8,4000._r8,2000._r8,3000._r8,1.e36_r8,1.e36_r8,4000._r8,4000._r8,8000._r8/
00187 data rac (5,1:NLUse) &
00188 / 100._r8, 50._r8, 80._r8,1200._r8,2000._r8,1500._r8, 0._r8, 0._r8, 200._r8, 60._r8, 120._r8/
00189 data rgss(5,1:NLUse) &
00190 / 500._r8, 150._r8, 350._r8, 500._r8, 500._r8, 200._r8, 0._r8,1000._r8, 0._r8, 250._r8, 400._r8/
00191 data rgso(5,1:NLUse) &
00192 / 300._r8, 150._r8, 200._r8, 200._r8, 200._r8, 300._r8,2000._r8, 400._r8,1000._r8, 180._r8, 200._r8/
00193 data rcls(5,1:NLUse) &
00194 /1.e36_r8,4000._r8,4000._r8,4000._r8,2000._r8,3000._r8,1.e36_r8,1.e36_r8,4000._r8,4000._r8,8000._r8/
00195 data rclo(5,1:NLUse) &
00196 /1.e36_r8,1000._r8, 500._r8, 500._r8,1500._r8, 700._r8,1.e36_r8,1.e36_r8, 600._r8, 800._r8, 800._r8/
00197
00198
00199
00200
00201 real(r8), public, dimension(NSeas,NLUse) :: z0
00202
00203 data z0 (1,1:NLUse) &
00204 /1.000_r8,0.250_r8,0.050_r8,1.000_r8,1.000_r8,1.000_r8,0.0006_r8,0.002_r8,0.150_r8,0.100_r8,0.100_r8/
00205 data z0 (2,1:NLUse) &
00206 /1.000_r8,0.100_r8,0.050_r8,1.000_r8,1.000_r8,1.000_r8,0.0006_r8,0.002_r8,0.100_r8,0.080_r8,0.080_r8/
00207 data z0 (3,1:NLUse) &
00208 /1.000_r8,0.005_r8,0.050_r8,1.000_r8,1.000_r8,1.000_r8,0.0006_r8,0.002_r8,0.100_r8,0.020_r8,0.060_r8/
00209 data z0 (4,1:NLUse) &
00210 /1.000_r8,0.001_r8,0.001_r8,1.000_r8,1.000_r8,1.000_r8,0.0006_r8,0.002_r8,0.001_r8,0.001_r8,0.040_r8/
00211 data z0 (5,1:NLUse) &
00212 /1.000_r8,0.030_r8,0.020_r8,1.000_r8,1.000_r8,1.000_r8,0.0006_r8,0.002_r8,0.010_r8,0.030_r8,0.060_r8/
00213
00214
00215
00216
00217
00218
00219
00220
00221
00222
00223
00224
00225
00226 real(r8), public, parameter :: dfoxd(n_species_table) =
00227 (/ 1._r8
00228 ,1._r8
00229 ,1._r8
00230 ,.1_r8
00231 ,1.e-36_r8
00232 ,1.e-36_r8
00233 ,1._r8
00234 ,.1_r8
00235 ,1.e-36_r8
00236 ,0._r8
00237 ,0._r8
00238 ,.1_r8
00239 ,1.e-36_r8
00240 ,1.e-36_r8
00241 ,1.e-36_r8
00242 ,.1_r8
00243 ,1._r8
00244 ,1.e-36_r8
00245 ,.1_r8
00246 ,1._r8
00247 ,1.e-36_r8
00248 ,.1_r8
00249 ,.1_r8
00250 ,.1_r8
00251 ,.1_r8
00252 ,1.e-36_r8
00253 ,1.e-36_r8
00254 ,.1_r8
00255 ,1.e-36_r8
00256 ,.1_r8
00257 ,1.e-36_r8
00258 ,.1_r8
00259 ,.1_r8
00260 ,1.e-36_r8
00261 ,1.e-36_r8
00262 ,1.e-36_r8
00263 ,1.e-36_r8
00264 ,.1_r8
00265 ,1.e-36_r8
00266 ,.1_r8
00267 ,1.e-36_r8
00268 ,.1_r8
00269 ,.1_r8
00270 ,.1_r8
00271 ,1.e-36_r8
00272 ,1.e-36_r8
00273 ,1.e-36_r8
00274 ,1.e-36_r8
00275 ,1.e-36_r8
00276 ,.1_r8
00277 ,.1_r8
00278 ,.1_r8
00279 ,1.e-36_r8
00280 /)
00281
00282
00283
00284
00285 Interface seq_drydep_setHCoeff
00286 Module Procedure set_hcoeff_scalar
00287 Module Procedure set_hcoeff_vector
00288 End Interface
00289
00290 real(r8), private, parameter :: inf = O'0777600000000000000000'
00291 real(r8), private, parameter :: small_value = 1.e-36_r8
00292
00293
00294
00295
00296
00297
00298 character(len=20), private, parameter :: species_name_table(n_species_table) =
00299 (/ 'OX '
00300 ,'H2O2 '
00301 ,'OH '
00302 ,'HO2 '
00303 ,'CO '
00304 ,'CH4 '
00305 ,'CH3O2 '
00306 ,'CH3OOH '
00307 ,'CH2O '
00308 ,'CHOOH '
00309 ,'NO '
00310 ,'NO2 '
00311 ,'HNO3 '
00312 ,'CO2 '
00313 ,'NH3 '
00314 ,'N2O5 '
00315 ,'NO3 '
00316 ,'CH3OH '
00317 ,'HO2NO2 '
00318 ,'O1D '
00319 ,'C2H6 '
00320 ,'C2H5O2 '
00321 ,'PO2 '
00322 ,'MACRO2 '
00323 ,'ISOPO2 '
00324 ,'C4H10 '
00325 ,'CH3CHO '
00326 ,'C2H5OOH '
00327 ,'C3H6 '
00328 ,'POOH '
00329 ,'C2H4 '
00330 ,'PAN '
00331 ,'CH3COOOH'
00332 ,'C10H16 '
00333 ,'CHOCHO '
00334 ,'CH3COCHO'
00335 ,'GLYALD '
00336 ,'CH3CO3 '
00337 ,'C3H8 '
00338 ,'C3H7O2 '
00339 ,'CH3COCH3'
00340 ,'C3H7OOH '
00341 ,'RO2 '
00342 ,'ROOH '
00343 ,'Rn '
00344 ,'ISOP '
00345 ,'MVK '
00346 ,'MACR '
00347 ,'C2H5OH '
00348 ,'ONITR '
00349 ,'ONIT '
00350 ,'ISOPNO3 '
00351 ,'HYDRALD '
00352 /)
00353
00354
00355 real(r8), private, parameter :: dheff(n_species_table*6) =
00356 (/1.15e-02_r8, 2560._r8,0._r8 , 0._r8,0._r8 , 0._r8
00357 ,8.33e+04_r8, 7379._r8,2.2e-12_r8,-3730._r8,0._r8 , 0._r8
00358 ,3.00e+01_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00359 ,2.00e+03_r8, 6600._r8,3.5e-05_r8, 0._r8,0._r8 , 0._r8
00360 ,1.00e-03_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00361 ,1.70e-03_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00362 ,7.47e+00_r8, 5241._r8,0._r8 , 0._r8,0._r8 , 0._r8
00363 ,3.11e+02_r8, 5241._r8,0._r8 , 0._r8,0._r8 , 0._r8
00364 ,6.30e+03_r8, 6425._r8,0._r8 , 0._r8,0._r8 , 0._r8
00365 ,5.53e+03_r8, 5700._r8,1.8e-04_r8,-1510._r8,0._r8 , 0._r8
00366 ,1.90e-03_r8, 1480._r8,0._r8 , 0._r8,0._r8 , 0._r8
00367 ,6.40e-03_r8, 2500._r8,0._r8 , 0._r8,0._r8 , 0._r8
00368 ,0._r8 , 0._r8,2.6e+06_r8, 8700._r8,0._r8 , 0._r8
00369 ,3.40e-02_r8, 2420._r8,4.5e-07_r8,-1000._r8,3.6e-11_r8,-1760._r8
00370 ,7.40e+01_r8, 3400._r8,1.7e-05_r8, -450._r8,1.0e-14_r8,-6716._r8
00371 ,2.14e+00_r8, 3362._r8,0._r8 , 0._r8,0._r8 , 0._r8
00372 ,0.65e+00_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00373 ,2.20e+02_r8, 4934._r8,0._r8 , 0._r8,0._r8 , 0._r8
00374 ,0._r8 , 0._r8,3.2e+01_r8, 0._r8,0._r8 , 0._r8
00375 ,1.00e-16_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00376 ,1.70e-03_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00377 ,7.47e+00_r8, 5241._r8,0._r8 , 0._r8,0._r8 , 0._r8
00378 ,7.47e+00_r8, 5241._r8,0._r8 , 0._r8,0._r8 , 0._r8
00379 ,7.47e+00_r8, 5241._r8,0._r8 , 0._r8,0._r8 , 0._r8
00380 ,7.47e+00_r8, 5241._r8,0._r8 , 0._r8,0._r8 , 0._r8
00381 ,1.70e-03_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00382 ,1.14e+01_r8, 6267._r8,0._r8 , 0._r8,0._r8 , 0._r8
00383 ,3.36e+02_r8, 5995._r8,0._r8 , 0._r8,0._r8 , 0._r8
00384 ,1.70e-03_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00385 ,2.20e+02_r8, 5653._r8,0._r8 , 0._r8,0._r8 , 0._r8
00386 ,1.70e-03_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00387 ,5.00e+00_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00388 ,8.37e+02_r8, 5308._r8,1.8e-04_r8,-1510._r8,0._r8 , 0._r8
00389 ,1.70e-03_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00390 ,3.00e+05_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00391 ,3.71e+03_r8, 7541._r8,0._r8 , 0._r8,0._r8 , 0._r8
00392 ,4.14e+04_r8, 4630._r8,0._r8 , 0._r8,0._r8 , 0._r8
00393 ,7.47e+00_r8, 5241._r8,0._r8 , 0._r8,0._r8 , 0._r8
00394 ,1.45e-03_r8, 2700._r8,0._r8 , 0._r8,0._r8 , 0._r8
00395 ,3.00e+06_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00396 ,2.70e+01_r8, 5300._r8,0._r8 , 0._r8,0._r8 , 0._r8
00397 ,3.36e+02_r8, 5995._r8,0._r8 , 0._r8,0._r8 , 0._r8
00398 ,7.47e+00_r8, 5241._r8,0._r8 , 0._r8,0._r8 , 0._r8
00399 ,3.36e+02_r8, 5995._r8,0._r8 , 0._r8,0._r8 , 0._r8
00400 ,0.00e+00_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00401 ,1.70e-03_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00402 ,1.70e-03_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00403 ,1.70e-03_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00404 ,2.00e+02_r8, 6500._r8,0._r8 , 0._r8,0._r8 , 0._r8
00405 ,7.51e+03_r8, 6485._r8,0._r8 , 0._r8,0._r8 , 0._r8
00406 ,1.00e+03_r8, 6000._r8,0._r8 , 0._r8,0._r8 , 0._r8
00407 ,1.00e+01_r8, 0._r8,0._r8 , 0._r8,0._r8 , 0._r8
00408 ,7.00e+01_r8, 6000._r8,0._r8 , 0._r8,0._r8 , 0._r8
00409 /)
00410
00411 real(r8), private, parameter :: wh2o = SHR_CONST_MWWV
00412 real(r8), private, parameter :: mol_wgts(n_species_table) =
00413 (/ 47.9981995_r8, 34.0135994_r8, 17.0067997_r8, 33.0061989_r8, 28.0104008_r8,
00414 16.0405998_r8, 47.0320015_r8, 48.0393982_r8, 30.0251999_r8, 46.0246010_r8,
00415 30.0061398_r8, 46.0055389_r8, 63.0123405_r8, 44.0098000_r8, 17.0289402_r8,
00416 108.010483_r8, 62.0049400_r8, 32.0400009_r8, 79.0117416_r8, 15.9994001_r8,
00417 30.0664005_r8, 61.0578003_r8, 91.0830002_r8, 119.093399_r8, 117.119797_r8,
00418 58.1180000_r8, 44.0509987_r8, 62.0652008_r8, 42.0774002_r8, 92.0904007_r8,
00419 28.0515995_r8, 121.047943_r8, 76.0497971_r8, 136.228394_r8, 58.0355988_r8,
00420 72.0614014_r8, 60.0503998_r8, 75.0423965_r8, 44.0922012_r8, 75.0836029_r8,
00421 58.0768013_r8, 76.0910034_r8, 31.9988003_r8, 33.0061989_r8, 222.000000_r8,
00422 68.1141968_r8, 70.0877991_r8, 70.0877991_r8, 46.0657997_r8, 147.125946_r8,
00423 119.074341_r8, 162.117935_r8, 100.112999_r8 /)
00424
00425
00426 CONTAINS
00427
00428
00429
00430
00431 subroutine seq_drydep_read(NLFilename, seq_drydep_fields)
00432
00433
00434
00435
00436
00437
00438
00439
00440
00441 use shr_file_mod,only : shr_file_getUnit, shr_file_freeUnit
00442 use shr_log_mod, only : s_logunit => shr_log_Unit
00443
00444 implicit none
00445
00446 character(len=*), intent(in) :: NLFilename
00447 character(len=*), intent(out) :: seq_drydep_fields
00448
00449
00450 integer :: i
00451 integer :: unitn
00452 integer :: ierr
00453 logical :: exists
00454 character(len=8) :: token
00455
00456
00457 character(*),parameter :: subName = '(seq_drydep_read) '
00458 character(*),parameter :: F00 = "('(seq_drydep_read) ',8a)"
00459 character(*),parameter :: FI1 = "('(seq_drydep_init) ',a,I2)"
00460
00461 namelist /drydep_inparm/ drydep_list, drydep_method
00462
00463
00464
00465
00466
00467
00468
00469 if ( len_trim(NLFilename) == 0 )then
00470 call shr_sys_abort( subName//'ERROR: nlfilename not set' )
00471 end if
00472 inquire( file=trim(NLFileName), exist=exists)
00473 if ( exists ) then
00474 unitn = shr_file_getUnit()
00475 open( unitn, file=trim(NLFilename), status='old' )
00476 if ( s_loglev > 0 ) write(s_logunit,F00) &
00477 'Read in drydep_inparm namelist from: ', trim(NLFilename)
00478 ierr = 1
00479 do while ( ierr /= 0 )
00480 read(unitn, drydep_inparm, iostat=ierr)
00481 if (ierr < 0) then
00482 call shr_sys_abort( subName//'ERROR: encountered end-of-file on namelist read' )
00483 endif
00484 end do
00485 close( unitn )
00486 call shr_file_freeUnit( unitn )
00487 end if
00488
00489 n_drydep = 0
00490
00491
00492 seq_drydep_fields = ' '
00493 do i=1,maxspc
00494 if ( len_trim(drydep_list(i))==0 ) exit
00495 write(token,333) i
00496 seq_drydep_fields = trim(seq_drydep_fields)//':'//trim(token)
00497 if ( i == 1 ) then
00498 drydep_fields_token = token
00499 endif
00500 n_drydep = n_drydep+1
00501 enddo
00502
00503
00504 lnd_drydep = n_drydep>0 .and. drydep_method == DD_XLND
00505
00506 if ( s_loglev > 0 ) then
00507 write(s_logunit,*) 'seq_drydep_read: drydep_method: ', trim(drydep_method)
00508 if ( n_drydep == 0 )then
00509 write(s_logunit,F00) 'No dry deposition fields will be transfered'
00510 else
00511 write(s_logunit,FI1) 'Number of dry deposition fields transfered is ', &
00512 n_drydep
00513 end if
00514 end if
00515
00516 if ( trim(drydep_method)/=trim(DD_XATM) .and. &
00517 trim(drydep_method)/=trim(DD_XLND) .and. &
00518 trim(drydep_method)/=trim(DD_TABL) ) then
00519 if ( s_loglev > 0 ) then
00520 write(s_logunit,*) 'seq_drydep_read: drydep_method : ', trim(drydep_method)
00521 write(s_logunit,*) 'seq_drydep_read: drydep_method must be set to : ', &
00522 DD_XATM,', ', DD_XLND,', or ', DD_TABL
00523 end if
00524 call shr_sys_abort('seq_drydep_read: incorrect dry deposition method specification')
00525 endif
00526
00527
00528 333 format ('dd',i3.3)
00529
00530 end subroutine seq_drydep_read
00531
00532
00533
00534 subroutine seq_drydep_init( )
00535
00536
00537
00538
00539
00540
00541
00542
00543
00544
00545
00546 use shr_log_mod, only : s_logunit => shr_log_Unit
00547
00548 implicit none
00549
00550
00551 integer :: i, l
00552 character(len=32) :: test_name
00553
00554 character(*),parameter :: subName = '(seq_drydep_init) '
00555 character(*),parameter :: F00 = "('(seq_drydep_init) ',8a)"
00556
00557
00558
00559
00560
00561 if ( n_drydep > 0 ) then
00562
00563 allocate( foxd(n_drydep) )
00564 allocate( drat(n_drydep) )
00565 allocate( mapping(n_drydep) )
00566
00567 foxd = inf
00568 drat = inf
00569
00570 mapping(:) = 0
00571
00572 end if
00573
00574 h2_ndx=-1; ch4_ndx=-1; co_ndx=-1; mpan_ndx = -1; pan_ndx = -1; so2_ndx=-1; o3_ndx=-1
00575
00576
00577 do i=1,n_drydep
00578 if ( len_trim(drydep_list(i))==0 ) exit
00579
00580 test_name = drydep_list(i)
00581
00582 if( trim(test_name) == 'O3' ) then
00583 test_name = 'OX'
00584 end if
00585
00586
00587 do l = 1,n_species_table
00588 if( trim( test_name ) == trim( species_name_table(l) ) ) then
00589 mapping(i) = l
00590 exit
00591 end if
00592 end do
00593
00594
00595 if( mapping(i) < 1 ) then
00596 select case( trim(test_name) )
00597 case( 'H2' )
00598 test_name = 'CO'
00599 case( 'HYAC', 'CH3COOH' )
00600 test_name = 'CH2O'
00601 case( 'O3S', 'O3INERT', 'MPAN' )
00602 test_name = 'OX'
00603 case( 'ISOPOOH', 'MACROOH', 'Pb', 'XOOH', 'H2SO4' )
00604 test_name = 'HNO3'
00605 case( 'ALKOOH', 'MEKOOH', 'TOLOOH', 'BENOOH', 'XYLOOH', 'TERPOOH','SOGM','SOGI','SOGT','SOGB','SOGX' )
00606 test_name = 'CH3OOH'
00607 case( 'SOA', 'SO2', 'SO4', 'CB1', 'CB2', 'OC1', 'OC2', 'NH3', 'NH4', 'SA1', 'SA2', 'SA3', 'SA4' )
00608 test_name = 'OX'
00609 case( 'SOAM', 'SOAI', 'SOAT', 'SOAB', 'SOAX' )
00610 test_name = 'OX'
00611
00612 case( 'O3A', 'XMPAN' )
00613 test_name = 'OX'
00614 case( 'XPAN' )
00615 test_name = 'PAN'
00616 case( 'XNO' )
00617 test_name = 'NO'
00618 case( 'XNO2' )
00619 test_name = 'NO2'
00620 case( 'XHNO3' )
00621 test_name = 'HNO3'
00622 case( 'XONIT' )
00623 test_name = 'ONIT'
00624 case( 'XONITR' )
00625 test_name = 'ONITR'
00626 case( 'XHO2NO2')
00627 test_name = 'HO2NO2'
00628 case( 'XNH4NO3' )
00629 test_name = 'HNO3'
00630
00631 case( 'NH4NO3' )
00632 test_name = 'HNO3'
00633 case default
00634 test_name = 'blank'
00635 end select
00636
00637
00638 if( trim(test_name) /= 'blank' ) then
00639 do l = 1,n_species_table
00640 if( trim( test_name ) == trim( species_name_table(l) ) ) then
00641 mapping(i) = l
00642 exit
00643 end if
00644 end do
00645 else
00646 if ( s_loglev > 0 ) write(s_logunit,F00) trim(drydep_list(i)), &
00647 ' not in tables; will have dep vel = 0'
00648 end if
00649 end if
00650
00651
00652 if ( trim(drydep_list(i)) == 'H2' ) h2_ndx = i
00653 if ( trim(drydep_list(i)) == 'CO' ) co_ndx = i
00654 if ( trim(drydep_list(i)) == 'CH4' ) ch4_ndx = i
00655 if ( trim(drydep_list(i)) == 'MPAN' ) mpan_ndx = i
00656 if ( trim(drydep_list(i)) == 'PAN' ) pan_ndx = i
00657 if ( trim(drydep_list(i)) == 'SO2' ) so2_ndx = i
00658 if ( trim(drydep_list(i)) == 'OX' .or. trim(drydep_list(i)) == 'O3' ) o3_ndx = i
00659 if ( trim(drydep_list(i)) == 'O3A' ) o3a_ndx = i
00660
00661 if( mapping(i) > 0) then
00662 l = mapping(i)
00663 foxd(i) = dfoxd(l)
00664 drat(i) = sqrt(mol_wgts(l)/wh2o)
00665 endif
00666
00667 enddo
00668
00669 where( rgss < 1._r8 )
00670 rgss = 1._r8
00671 endwhere
00672
00673 where( rac < small_value)
00674 rac = small_value
00675 endwhere
00676
00677 end subroutine seq_drydep_init
00678
00679
00680
00681 subroutine set_hcoeff_scalar( sfc_temp, heff )
00682
00683
00684
00685
00686
00687
00688
00689
00690
00691
00692 implicit none
00693
00694 real(r8), intent(in) :: sfc_temp
00695 real(r8), intent(out) :: heff(n_drydep)
00696
00697
00698 real(r8) :: sfc_temp_tmp(1)
00699
00700 sfc_temp_tmp(:) = sfc_temp
00701 call set_hcoeff_vector( 1, sfc_temp_tmp, heff(:n_drydep) )
00702
00703 end subroutine set_hcoeff_scalar
00704
00705
00706
00707 subroutine set_hcoeff_vector( ncol, sfc_temp, heff )
00708
00709
00710
00711
00712
00713
00714
00715
00716 use shr_log_mod, only : s_logunit => shr_log_Unit
00717
00718 implicit none
00719
00720 integer, intent(in) :: ncol
00721 real(r8), intent(in) :: sfc_temp(ncol)
00722 real(r8), intent(out) :: heff(ncol,n_drydep)
00723
00724
00725 real(r8), parameter :: t0 = 298._r8
00726 real(r8), parameter :: ph_inv = 1._r8/ph
00727 integer :: m, l, id
00728 real(r8) :: e298
00729 real(r8) :: dhr
00730 real(r8) :: dk1s(ncol)
00731 real(r8) :: dk2s(ncol)
00732 real(r8) :: wrk(ncol)
00733
00734
00735 character(*),parameter :: subName = '(seq_drydep_set_hcoeff) '
00736 character(*),parameter :: F00 = "('(seq_drydep_set_hcoeff) ',8a)"
00737
00738
00739
00740
00741
00742 wrk(:) = (t0 - sfc_temp(:))/(t0*sfc_temp(:))
00743 do m = 1,n_drydep
00744 l = mapping(m)
00745 id = 6*(l - 1)
00746 e298 = dheff(id+1)
00747 dhr = dheff(id+2)
00748 heff(:,m) = e298*exp( dhr*wrk(:) )
00749
00750 if( dheff(id+3) /= 0._r8 .and. dheff(id+5) == 0._r8 ) then
00751 e298 = dheff(id+3)
00752 dhr = dheff(id+4)
00753 dk1s(:) = e298*exp( dhr*wrk(:) )
00754 where( heff(:,m) /= 0._r8 )
00755 heff(:,m) = heff(:,m)*(1._r8 + dk1s(:)*ph_inv)
00756 elsewhere
00757 heff(:,m) = dk1s(:)*ph_inv
00758 endwhere
00759 end if
00760
00761 if( dheff(id+5) /= 0._r8 ) then
00762 if( trim( drydep_list(m) ) == 'CO2' .or. trim( drydep_list(m) ) == 'NH3' ) then
00763 e298 = dheff(id+3)
00764 dhr = dheff(id+4)
00765 dk1s(:) = e298*exp( dhr*wrk(:) )
00766 e298 = dheff(id+5)
00767 dhr = dheff(id+6)
00768 dk2s(:) = e298*exp( dhr*wrk(:) )
00769
00770 if( trim(drydep_list(m)) == 'CO2' ) then
00771 heff(:,m) = heff(:,m)*(1._r8 + dk1s(:)*ph_inv)*(1._r8 + dk2s(:)*ph_inv)
00772
00773 else if( trim( drydep_list(m) ) == 'NH3' ) then
00774 heff(:,m) = heff(:,m)*(1._r8 + dk1s(:)*ph/dk2s(:))
00775
00776 else
00777 write(s_logunit,F00) 'Bad species ',drydep_list(m)
00778 call shr_sys_abort( subName//'ERROR: in assigning coefficients' )
00779 end if
00780 end if
00781 end if
00782 end do
00783
00784 end subroutine set_hcoeff_vector
00785
00786
00787
00788 end module seq_drydep_mod