1 | #undef DEBUG |
---|
2 | #define WRF_PORT |
---|
3 | |
---|
4 | module cldwat |
---|
5 | !----------------------------------------------------------------------- |
---|
6 | ! |
---|
7 | ! Purpose: Prognostic cloud water data and methods. |
---|
8 | ! |
---|
9 | ! Public interfaces: |
---|
10 | ! |
---|
11 | ! inimc -- Initialize constants |
---|
12 | ! pcond -- Calculate prognostic condensate |
---|
13 | ! |
---|
14 | ! cldwat_fice -- calculate fraction of condensate in ice phase (radiation partitioning) |
---|
15 | ! |
---|
16 | ! Author: P. Rasch, with Modifications by Minghua Zhang |
---|
17 | ! January 2010, modified by J. Kay to add precip fluxes for COSP simulator |
---|
18 | ! Ported to WRF by William.Gustafson@pnl.gov, Nov. 2009 |
---|
19 | ! updated to CESM1_0_1, Nov. 2010 |
---|
20 | ! |
---|
21 | !----------------------------------------------------------------------- |
---|
22 | use shr_kind_mod, only: r8 => shr_kind_r8 |
---|
23 | use physconst, only: tmelt |
---|
24 | #ifndef WRF_PORT |
---|
25 | use spmd_utils, only: masterproc |
---|
26 | use ppgrid, only: pcols, pver, pverp |
---|
27 | use wv_saturation, only: estblf, hlatv, tmin, hlatf, rgasv, pcf, & |
---|
28 | cp, epsqs, ttrice |
---|
29 | use abortutils, only: endrun |
---|
30 | use cam_logfile, only: iulog |
---|
31 | #else |
---|
32 | use module_cam_support, only: masterproc, & |
---|
33 | pcols, pver, pverp, & |
---|
34 | endrun, & |
---|
35 | iulog |
---|
36 | #endif |
---|
37 | |
---|
38 | implicit none |
---|
39 | |
---|
40 | !----------------------------------------------------------------------- |
---|
41 | ! PUBLIC: Make default data and interfaces private |
---|
42 | !----------------------------------------------------------------------- |
---|
43 | private |
---|
44 | save |
---|
45 | #ifndef WRF_PORT |
---|
46 | public inimc, pcond ! Public interfaces |
---|
47 | #endif |
---|
48 | public cldwat_fice ! Public interfaces |
---|
49 | #ifndef WRF_PORT |
---|
50 | public cldwat_readnl |
---|
51 | #endif |
---|
52 | integer, public:: ktop ! Level above 10 hPa |
---|
53 | |
---|
54 | #ifndef WRF_PORT |
---|
55 | real(r8),public :: icritc ! threshold for autoconversion of cold ice |
---|
56 | real(r8),public :: icritw ! threshold for autoconversion of warm ice |
---|
57 | !!$ real(r8),public,parameter:: conke = 1.e-6 ! tunable constant for evaporation of precip |
---|
58 | !!$ real(r8),public,parameter:: conke = 2.e-6 ! tunable constant for evaporation of precip |
---|
59 | real(r8),public :: conke ! tunable constant for evaporation of precip |
---|
60 | #else |
---|
61 | ! Currently, the WRF_PORT bipasses the namelist initialization of these |
---|
62 | ! tunable parameters. We are hard-coding them to the default values for |
---|
63 | ! the fv 0.23x0.31 grid. One might choose to implement an option to set |
---|
64 | ! these via the WRF Registry in the future. |
---|
65 | real(r8),public :: icritw = 2.0e-4 |
---|
66 | real(r8),public :: icritc = 45.0e-6 |
---|
67 | real(r8),public :: conke = 5.0e-6 |
---|
68 | |
---|
69 | #endif |
---|
70 | |
---|
71 | !----------------------------------------------------------------------- |
---|
72 | ! PRIVATE: Everything else is private to this module |
---|
73 | !----------------------------------------------------------------------- |
---|
74 | real(r8), private, parameter :: tmax_fice = tmelt - 10._r8 ! max temperature for cloud ice formation |
---|
75 | !!$ real(r8), private, parameter :: tmax_fice = tmelt ! max temperature for cloud ice formation |
---|
76 | !!$ real(r8), private, parameter :: tmin_fice = tmax_fice - 20.! min temperature for cloud ice formation |
---|
77 | real(r8), private, parameter :: tmin_fice = tmax_fice - 30._r8 ! min temperature for cloud ice formation |
---|
78 | ! pjr |
---|
79 | real(r8), private, parameter :: tmax_fsnow = tmelt ! max temperature for transition to convective snow |
---|
80 | real(r8), private, parameter :: tmin_fsnow = tmelt-5._r8 ! min temperature for transition to convective snow |
---|
81 | real(r8), private:: rhonot ! air density at surface |
---|
82 | real(r8), private:: t0 ! Freezing temperature |
---|
83 | real(r8), private:: cldmin ! assumed minimum cloud amount |
---|
84 | real(r8), private:: small ! small number compared to unity |
---|
85 | real(r8), private:: c ! constant for graupel like snow cm**(1-d)/s |
---|
86 | real(r8), private:: d ! constant for graupel like snow |
---|
87 | real(r8), private:: esi ! collection efficient for ice by snow |
---|
88 | real(r8), private:: esw ! collection efficient for water by snow |
---|
89 | real(r8), private:: nos ! particles snow / cm**4 |
---|
90 | real(r8), private:: pi ! Mathematical constant |
---|
91 | real(r8), private:: gravit ! Gravitational acceleration at surface |
---|
92 | real(r8), private:: rh2o |
---|
93 | real(r8), private:: prhonos |
---|
94 | real(r8), private:: thrpd ! numerical three added to d |
---|
95 | real(r8), private:: gam3pd ! gamma function on (3+d) |
---|
96 | real(r8), private:: gam4pd ! gamma function on (4+d) |
---|
97 | real(r8), private:: rhoi ! ice density |
---|
98 | real(r8), private:: rhos ! snow density |
---|
99 | real(r8), private:: rhow ! water density |
---|
100 | real(r8), private:: mcon01 ! constants used in cloud microphysics |
---|
101 | real(r8), private:: mcon02 ! constants used in cloud microphysics |
---|
102 | real(r8), private:: mcon03 ! constants used in cloud microphysics |
---|
103 | real(r8), private:: mcon04 ! constants used in cloud microphysics |
---|
104 | real(r8), private:: mcon05 ! constants used in cloud microphysics |
---|
105 | real(r8), private:: mcon06 ! constants used in cloud microphysics |
---|
106 | real(r8), private:: mcon07 ! constants used in cloud microphysics |
---|
107 | real(r8), private:: mcon08 ! constants used in cloud microphysics |
---|
108 | |
---|
109 | integer, private :: k1mb ! index of the eta level near 1 mb |
---|
110 | |
---|
111 | ! Parameters used in findmcnew |
---|
112 | real(r8) :: r3lcrit ! critical radius at which autoconversion become efficient |
---|
113 | real(r8) :: capnsi ! sea ice cloud particles / cm3 |
---|
114 | real(r8) :: capnc ! cold and oceanic cloud particles / cm3 |
---|
115 | real(r8) :: capnw ! warm continental cloud particles / cm3 |
---|
116 | real(r8) :: kconst ! const for terminal velocity (stokes regime) |
---|
117 | real(r8) :: effc ! collection efficiency |
---|
118 | real(r8) :: alpha ! ratio of 3rd moment radius to 2nd |
---|
119 | real(r8) :: capc ! constant for autoconversion |
---|
120 | real(r8) :: convfw ! constant used for fall velocity calculation |
---|
121 | real(r8) :: cracw ! constant used for rain accreting water |
---|
122 | real(r8) :: critpr ! critical precip rate collection efficiency changes |
---|
123 | real(r8) :: ciautb ! coefficient of autoconversion of ice (1/s) |
---|
124 | |
---|
125 | #ifdef DEBUG |
---|
126 | integer, private,parameter :: nlook = 1 ! Number of points to examine |
---|
127 | integer, private :: ilook(nlook) ! Longitude index to examine |
---|
128 | integer, private :: latlook(nlook) ! Latitude index to examine |
---|
129 | integer, private :: lchnklook(nlook) ! Chunk index to examine |
---|
130 | integer, private :: icollook(nlook) ! Column index to examine |
---|
131 | #endif |
---|
132 | ! Private data |
---|
133 | real(r8), parameter :: unset_r8 = huge(1.0_r8) |
---|
134 | |
---|
135 | contains |
---|
136 | !=============================================================================== |
---|
137 | subroutine cldwat_readnl(nlfile) |
---|
138 | #ifndef WRF_PORT |
---|
139 | use namelist_utils, only: find_group_name |
---|
140 | use units, only: getunit, freeunit |
---|
141 | use mpishorthand |
---|
142 | #endif |
---|
143 | |
---|
144 | character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input |
---|
145 | #ifndef WRF_PORT |
---|
146 | |
---|
147 | |
---|
148 | |
---|
149 | ! Namelist variables |
---|
150 | real(r8) :: cldwat_icritw = unset_r8 ! icritw = threshold for autoconversion of warm ice |
---|
151 | real(r8) :: cldwat_icritc = unset_r8 ! icritc = threshold for autoconversion of cold ice |
---|
152 | real(r8) :: cldwat_conke = unset_r8 ! conke = tunable constant for evaporation of precip |
---|
153 | |
---|
154 | ! Local variables |
---|
155 | integer :: unitn, ierr |
---|
156 | character(len=*), parameter :: subname = 'cldwat_readnl' |
---|
157 | |
---|
158 | namelist /cldwat_nl/ cldwat_icritw, cldwat_icritc, cldwat_conke |
---|
159 | |
---|
160 | !----------------------------------------------------------------------------- |
---|
161 | |
---|
162 | if (masterproc) then |
---|
163 | unitn = getunit() |
---|
164 | open( unitn, file=trim(nlfile), status='old' ) |
---|
165 | call find_group_name(unitn, 'cldwat_nl', status=ierr) |
---|
166 | if (ierr == 0) then |
---|
167 | read(unitn, cldwat_nl, iostat=ierr) |
---|
168 | if (ierr /= 0) then |
---|
169 | call endrun(subname // ':: ERROR reading namelist') |
---|
170 | end if |
---|
171 | end if |
---|
172 | close(unitn) |
---|
173 | call freeunit(unitn) |
---|
174 | |
---|
175 | ! set local variables |
---|
176 | icritw = cldwat_icritw |
---|
177 | icritc = cldwat_icritc |
---|
178 | conke = cldwat_conke |
---|
179 | |
---|
180 | end if |
---|
181 | |
---|
182 | |
---|
183 | |
---|
184 | #ifdef SPMD |
---|
185 | ! Broadcast namelist variables |
---|
186 | call mpibcast(icritw, 1, mpir8, 0, mpicom) |
---|
187 | call mpibcast(icritc, 1, mpir8, 0, mpicom) |
---|
188 | call mpibcast(conke, 1, mpir8, 0, mpicom) |
---|
189 | #endif |
---|
190 | |
---|
191 | !endif for WRF_PORT: |
---|
192 | #endif |
---|
193 | |
---|
194 | end subroutine cldwat_readnl |
---|
195 | |
---|
196 | !================================================================================================ |
---|
197 | subroutine cldwat_fice(ncol, t, fice, fsnow) |
---|
198 | ! |
---|
199 | ! Compute the fraction of the total cloud water which is in ice phase. |
---|
200 | ! The fraction depends on temperature only. |
---|
201 | ! This is the form that was used for radiation, the code came from cldefr originally |
---|
202 | ! |
---|
203 | ! Author: B. A. Boville Sept 10, 2002 |
---|
204 | ! modified: PJR 3/13/03 (added fsnow to ascribe snow production for convection ) |
---|
205 | !----------------------------------------------------------------------- |
---|
206 | implicit none |
---|
207 | |
---|
208 | ! Arguments |
---|
209 | integer, intent(in) :: ncol ! number of active columns |
---|
210 | real(r8), intent(in) :: t(pcols,pver) ! temperature |
---|
211 | |
---|
212 | real(r8), intent(out) :: fice(pcols,pver) ! Fractional ice content within cloud |
---|
213 | real(r8), intent(out) :: fsnow(pcols,pver) ! Fractional snow content for convection |
---|
214 | |
---|
215 | ! Local variables |
---|
216 | integer :: i,k ! loop indexes |
---|
217 | |
---|
218 | !----------------------------------------------------------------------- |
---|
219 | |
---|
220 | ! Define fractional amount of cloud that is ice |
---|
221 | do k=1,pver |
---|
222 | do i=1,ncol |
---|
223 | |
---|
224 | ! If warmer than tmax then water phase |
---|
225 | if (t(i,k) > tmax_fice) then |
---|
226 | fice(i,k) = 0.0_r8 |
---|
227 | |
---|
228 | ! If colder than tmin then ice phase |
---|
229 | else if (t(i,k) < tmin_fice) then |
---|
230 | fice(i,k) = 1.0_r8 |
---|
231 | |
---|
232 | ! Otherwise mixed phase, with ice fraction decreasing linearly from tmin to tmax |
---|
233 | else |
---|
234 | fice(i,k) =(tmax_fice - t(i,k)) / (tmax_fice - tmin_fice) |
---|
235 | end if |
---|
236 | |
---|
237 | ! snow fraction partitioning |
---|
238 | |
---|
239 | ! If warmer than tmax then water phase |
---|
240 | if (t(i,k) > tmax_fsnow) then |
---|
241 | fsnow(i,k) = 0.0_r8 |
---|
242 | |
---|
243 | ! If colder than tmin then ice phase |
---|
244 | else if (t(i,k) < tmin_fsnow) then |
---|
245 | fsnow(i,k) = 1.0_r8 |
---|
246 | |
---|
247 | ! Otherwise mixed phase, with ice fraction decreasing linearly from tmin to tmax |
---|
248 | else |
---|
249 | fsnow(i,k) =(tmax_fsnow - t(i,k)) / (tmax_fsnow - tmin_fsnow) |
---|
250 | end if |
---|
251 | |
---|
252 | end do |
---|
253 | end do |
---|
254 | |
---|
255 | return |
---|
256 | end subroutine cldwat_fice |
---|
257 | |
---|
258 | #ifndef WRF_PORT |
---|
259 | subroutine inimc( tmeltx, rhonotx, gravitx, rh2ox) |
---|
260 | !----------------------------------------------------------------------- |
---|
261 | ! |
---|
262 | ! Purpose: |
---|
263 | ! initialize constants for the prognostic condensate |
---|
264 | ! |
---|
265 | ! Author: P. Rasch, April 1997 |
---|
266 | ! |
---|
267 | !----------------------------------------------------------------------- |
---|
268 | use pmgrid, only: plev, plevp |
---|
269 | use dycore, only: dycore_is, get_resolution |
---|
270 | use hycoef, only: hypm |
---|
271 | |
---|
272 | integer k |
---|
273 | real(r8), intent(in) :: tmeltx |
---|
274 | real(r8), intent(in) :: rhonotx |
---|
275 | real(r8), intent(in) :: gravitx |
---|
276 | real(r8), intent(in) :: rh2ox |
---|
277 | |
---|
278 | #ifdef UNICOSMP |
---|
279 | real(r8) signgam ! variable required by cray gamma function |
---|
280 | external gamma |
---|
281 | #endif |
---|
282 | rhonot = rhonotx ! air density at surface (gm/cm3) |
---|
283 | gravit = gravitx |
---|
284 | rh2o = rh2ox |
---|
285 | rhos = .1_r8 ! assumed snow density (gm/cm3) |
---|
286 | rhow = 1._r8 ! water density |
---|
287 | rhoi = 1._r8 ! ice density |
---|
288 | esi = 1.0_r8 ! collection efficient for ice by snow |
---|
289 | esw = 0.1_r8 ! collection efficient for water by snow |
---|
290 | t0 = tmeltx ! approximate freezing temp |
---|
291 | cldmin = 0.02_r8 ! assumed minimum cloud amount |
---|
292 | small = 1.e-22_r8 ! a small number compared to unity |
---|
293 | c = 152.93_r8 ! constant for graupel like snow cm**(1-d)/s |
---|
294 | d = 0.25_r8 ! constant for graupel like snow |
---|
295 | nos = 3.e-2_r8 ! particles snow / cm**4 |
---|
296 | pi = 4._r8*atan(1.0_r8) |
---|
297 | prhonos = pi*rhos*nos |
---|
298 | thrpd = 3._r8 + d |
---|
299 | if (d==0.25_r8) then |
---|
300 | gam3pd = 2.549256966718531_r8 ! only right for d = 0.25 |
---|
301 | gam4pd = 8.285085141835282_r8 |
---|
302 | else |
---|
303 | #ifdef UNICOSMP |
---|
304 | call gamma(3._r8+d, signgam, gam3pd) |
---|
305 | gam3pd = sign(exp(gam3pd),signgam) |
---|
306 | call gamma(4._r8+d, signgam, gam4pd) |
---|
307 | gam4pd = sign(exp(gam4pd),signgam) |
---|
308 | write(iulog,*) ' d, gamma(3+d), gamma(4+d) =', gam3pd, gam4pd |
---|
309 | #else |
---|
310 | write(iulog,*) ' can only use d ne 0.25 on a cray ' |
---|
311 | stop |
---|
312 | #endif |
---|
313 | endif |
---|
314 | mcon01 = pi*nos*c*gam3pd/4._r8 |
---|
315 | mcon02 = 1._r8/(c*gam4pd*sqrt(rhonot)/(6*prhonos**(d/4._r8))) |
---|
316 | mcon03 = -(0.5_r8+d/4._r8) |
---|
317 | mcon04 = 4._r8/(4._r8+d) |
---|
318 | mcon05 = (3+d)/(4+d) |
---|
319 | mcon06 = (3+d)/4._r8 |
---|
320 | mcon07 = mcon01*sqrt(rhonot)*mcon02**mcon05/prhonos**mcon06 |
---|
321 | mcon08 = -0.5_r8/(4._r8+d) |
---|
322 | |
---|
323 | ! find the level about 1mb, we wont do the microphysics above this level |
---|
324 | k1mb = 1 |
---|
325 | do k=1,pver-1 |
---|
326 | if (hypm(k) < 1.e2_r8 .and. hypm(k+1) >= 1.e2_r8) then |
---|
327 | if (1.e2_r8-hypm(k) < hypm(k+1)-1.e2_r8) then |
---|
328 | k1mb = k |
---|
329 | else |
---|
330 | k1mb = k + 1 |
---|
331 | end if |
---|
332 | goto 20 |
---|
333 | end if |
---|
334 | end do |
---|
335 | if (masterproc) then |
---|
336 | write(iulog,*)'inimc: model levels bracketing 1 mb not found' |
---|
337 | end if |
---|
338 | ! call endrun |
---|
339 | k1mb = 1 |
---|
340 | 20 if( masterproc ) write(iulog,*)'inimc: model level nearest 1 mb is',k1mb,'which is',hypm(k1mb),'pascals' |
---|
341 | |
---|
342 | if( masterproc ) write(iulog,*) 'cloud water initialization by inimc complete ' |
---|
343 | |
---|
344 | ! Initialize parameters used by findmcnew |
---|
345 | capnw = 400._r8 ! warm continental cloud particles / cm3 |
---|
346 | capnc = 150._r8 ! cold and oceanic cloud particles / cm3 |
---|
347 | ! capnsi = 40._r8 ! sea ice cloud particles density / cm3 |
---|
348 | capnsi = 75._r8 ! sea ice cloud particles density / cm3 |
---|
349 | |
---|
350 | kconst = 1.18e6_r8 ! const for terminal velocity |
---|
351 | |
---|
352 | ! effc = 1._r8 ! autoconv collection efficiency following boucher 96 |
---|
353 | ! effc = .55*0.05_r8 ! autoconv collection efficiency following baker 93 |
---|
354 | effc = 0.55_r8 ! autoconv collection efficiency following tripoli and cotton |
---|
355 | ! effc = 0._r8 ! turn off warm-cloud autoconv |
---|
356 | alpha = 1.1_r8**4 |
---|
357 | capc = pi**(-.333_r8)*kconst*effc *(0.75_r8)**(1.333_r8)*alpha ! constant for autoconversion |
---|
358 | |
---|
359 | !r3lcrit: critical radius where liq conversion begins (10.0 micron) |
---|
360 | r3lcrit = 10.0e-6_r8 |
---|
361 | |
---|
362 | ! critical precip rate at which we assume the collector drops can change the |
---|
363 | ! drop size enough to enhance the auto-conversion process (mm/day) |
---|
364 | critpr = 0.5_r8 |
---|
365 | convfw = 1.94_r8*2.13_r8*sqrt(rhow*1000._r8*9.81_r8*2.7e-4_r8) |
---|
366 | |
---|
367 | ! liquid microphysics |
---|
368 | ! cracw = 6_r8 ! beheng |
---|
369 | cracw = .884_r8*sqrt(9.81_r8/(rhow*1000._r8*2.7e-4_r8)) ! tripoli and cotton |
---|
370 | |
---|
371 | ! ice microphysics |
---|
372 | ciautb = 5.e-4_r8 |
---|
373 | |
---|
374 | if ( masterproc ) then |
---|
375 | write(iulog,*)'tuning parameters cldwat: icritw',icritw,'icritc',icritc,'conke',conke |
---|
376 | write(iulog,*)'tuning parameters cldwat: capnw',capnw,'capnc',capnc,'capnsi',capnsi,'kconst',kconst |
---|
377 | write(iulog,*)'tuning parameters cldwat: effc',effc,'alpha',alpha,'capc',capc,'r3lcrit',r3lcrit |
---|
378 | write(iulog,*)'tuning parameters cldwat: critpr',critpr,'convfw',convfw,'cracw',cracw,'ciautb',ciautb |
---|
379 | endif |
---|
380 | |
---|
381 | return |
---|
382 | end subroutine inimc |
---|
383 | |
---|
384 | subroutine pcond (lchnk ,ncol , & |
---|
385 | tn ,ttend ,qn ,qtend ,omega , & |
---|
386 | cwat ,p ,pdel ,cldn ,fice , fsnow, & |
---|
387 | cme ,prodprec,prodsnow,evapprec,evapsnow,evapheat, prfzheat, & |
---|
388 | meltheat,precip ,snowab ,deltat ,fwaut , & |
---|
389 | fsaut ,fracw ,fsacw ,fsaci ,lctend , & |
---|
390 | rhdfda ,rhu00 ,seaicef, zi ,ice2pr, liq2pr, & |
---|
391 | liq2snow, snowh, rkflxprc, rkflxsnw) |
---|
392 | !----------------------------------------------------------------------- |
---|
393 | ! |
---|
394 | ! Purpose: |
---|
395 | ! The public interface to the cloud water parameterization |
---|
396 | ! returns tendencies to water vapor, temperature and cloud water variables |
---|
397 | ! |
---|
398 | ! For basic method |
---|
399 | ! See: Rasch, P. J, and J. E. Kristjansson, A Comparison of the CCM3 |
---|
400 | ! model climate using diagnosed and |
---|
401 | ! predicted condensate parameterizations, 1998, J. Clim., 11, |
---|
402 | ! pp1587---1614. |
---|
403 | ! |
---|
404 | ! For important modifications to improve the method of determining |
---|
405 | ! condensation/evaporation see Zhang et al (2001, in preparation) |
---|
406 | ! |
---|
407 | ! Authors: M. Zhang, W. Lin, P. Rasch and J.E. Kristjansson |
---|
408 | ! B. A. Boville (latent heat of fusion) |
---|
409 | !----------------------------------------------------------------------- |
---|
410 | use wv_saturation, only: vqsatd |
---|
411 | use cam_control_mod, only: nlvdry |
---|
412 | ! |
---|
413 | !--------------------------------------------------------------------- |
---|
414 | ! |
---|
415 | ! Input Arguments |
---|
416 | ! |
---|
417 | integer, intent(in) :: lchnk ! chunk identifier |
---|
418 | integer, intent(in) :: ncol ! number of atmospheric columns |
---|
419 | |
---|
420 | real(r8), intent(in) :: fice(pcols,pver) ! fraction of cwat that is ice |
---|
421 | real(r8), intent(in) :: fsnow(pcols,pver) ! fraction of rain that freezes to snow |
---|
422 | real(r8), intent(in) :: cldn(pcols,pver) ! new value of cloud fraction (fraction) |
---|
423 | real(r8), intent(in) :: cwat(pcols,pver) ! cloud water (kg/kg) |
---|
424 | real(r8), intent(in) :: omega(pcols,pver) ! vert pressure vel (Pa/s) |
---|
425 | real(r8), intent(in) :: p(pcols,pver) ! pressure (K) |
---|
426 | real(r8), intent(in) :: pdel(pcols,pver) ! pressure thickness (Pa) |
---|
427 | real(r8), intent(in) :: qn(pcols,pver) ! new water vapor (kg/kg) |
---|
428 | real(r8), intent(in) :: qtend(pcols,pver) ! mixing ratio tend (kg/kg/s) |
---|
429 | real(r8), intent(in) :: tn(pcols,pver) ! new temperature (K) |
---|
430 | real(r8), intent(in) :: ttend(pcols,pver) ! temp tendencies (K/s) |
---|
431 | real(r8), intent(in) :: deltat ! time step to advance solution over |
---|
432 | real(r8), intent(in) :: lctend(pcols,pver) ! cloud liquid water tendencies ====wlin |
---|
433 | real(r8), intent(in) :: rhdfda(pcols,pver) ! dG(a)/da, rh=G(a), when rh>u00 ====wlin |
---|
434 | real(r8), intent(in) :: rhu00 (pcols,pver) ! Rhlim for cloud ====wlin |
---|
435 | real(r8), intent(in) :: seaicef(pcols) ! sea ice fraction (fraction) |
---|
436 | real(r8), intent(in) :: zi(pcols,pverp) ! layer interfaces (m) |
---|
437 | real(r8), intent(in) :: snowh(pcols) ! Snow depth over land, water equivalent (m) |
---|
438 | ! |
---|
439 | ! Output Arguments |
---|
440 | ! |
---|
441 | real(r8), intent(out) :: cme (pcols,pver) ! rate of cond-evap of condensate (1/s) |
---|
442 | real(r8), intent(out) :: prodprec(pcols,pver) ! rate of conversion of condensate to precip (1/s) |
---|
443 | real(r8), intent(out) :: evapprec(pcols,pver) ! rate of evaporation of falling precip (1/s) |
---|
444 | real(r8), intent(out) :: evapsnow(pcols,pver) ! rate of evaporation of falling snow (1/s) |
---|
445 | real(r8), intent(out) :: evapheat(pcols,pver) ! heating rate due to evaporation of precip (W/kg) |
---|
446 | real(r8), intent(out) :: prfzheat(pcols,pver) ! heating rate due to freezing of precip (W/kg) |
---|
447 | real(r8), intent(out) :: meltheat(pcols,pver) ! heating rate due to snow melt (W/kg) |
---|
448 | real(r8), intent(out) :: precip(pcols) ! rate of precipitation (kg / (m**2 * s)) |
---|
449 | real(r8), intent(out) :: snowab(pcols) ! rate of snow (kg / (m**2 * s)) |
---|
450 | real(r8), intent(out) :: ice2pr(pcols,pver) ! rate of conversion of ice to precip |
---|
451 | real(r8), intent(out) :: liq2pr(pcols,pver) ! rate of conversion of liquid to precip |
---|
452 | real(r8), intent(out) :: liq2snow(pcols,pver) ! rate of conversion of liquid to snow |
---|
453 | real(r8), intent(out) :: rkflxprc(pcols,pverp) ! grid-box mean RK flux_large_scale_cloud_rain at interfaces (kg m^-2 s^-1) |
---|
454 | real(r8), intent(out) :: rkflxsnw(pcols,pverp) ! grid-box mean RK flux_large_scale_cloud_snow at interfaces (kg m^-2 s^-1) |
---|
455 | |
---|
456 | real(r8) nice2pr ! rate of conversion of ice to snow |
---|
457 | real(r8) nliq2pr ! rate of conversion of liquid to precip |
---|
458 | real(r8) nliq2snow ! rate of conversion of liquid to snow |
---|
459 | real(r8) prodsnow(pcols,pver) ! rate of production of snow |
---|
460 | |
---|
461 | ! |
---|
462 | ! Local workspace |
---|
463 | ! |
---|
464 | real(r8) :: precab(pcols) ! rate of precipitation (kg / (m**2 * s)) |
---|
465 | integer i ! work variable |
---|
466 | integer iter ! #iterations for precipitation calculation |
---|
467 | integer k ! work variable |
---|
468 | integer l ! work variable |
---|
469 | |
---|
470 | real(r8) cldm(pcols) ! mean cloud fraction over the time step |
---|
471 | real(r8) cldmax(pcols) ! max cloud fraction above |
---|
472 | real(r8) coef(pcols) ! conversion time scale for condensate to rain |
---|
473 | real(r8) cwm(pcols) ! cwat mixing ratio at midpoint of time step |
---|
474 | real(r8) cwn(pcols) ! cwat mixing ratio at end |
---|
475 | real(r8) denom ! work variable |
---|
476 | real(r8) dqsdt ! change in sat spec. hum. wrt temperature |
---|
477 | real(r8) es(pcols) ! sat. vapor pressure |
---|
478 | real(r8) fracw(pcols,pver) ! relative importance of collection of liquid by rain |
---|
479 | real(r8) fsaci(pcols,pver) ! relative importance of collection of ice by snow |
---|
480 | real(r8) fsacw(pcols,pver) ! relative importance of collection of liquid by snow |
---|
481 | real(r8) fsaut(pcols,pver) ! relative importance of ice auto conversion |
---|
482 | real(r8) fwaut(pcols,pver) ! relative importance of warm cloud autoconversion |
---|
483 | real(r8) gamma(pcols) ! d qs / dT |
---|
484 | real(r8) icwc(pcols) ! in-cloud water content (kg/kg) |
---|
485 | real(r8) mincld ! a small cloud fraction to avoid / zero |
---|
486 | real(r8) omeps ! 1 minus epsilon |
---|
487 | real(r8),parameter ::omsm=0.99999_r8 ! a number just less than unity (for rounding) |
---|
488 | real(r8) prprov(pcols) ! provisional value of precip at btm of layer |
---|
489 | real(r8) prtmp ! work variable |
---|
490 | real(r8) q(pcols,pver) ! mixing ratio before time step ignoring condensate |
---|
491 | real(r8) qs(pcols) ! spec. hum. of water vapor |
---|
492 | real(r8) qsn, esn ! work variable |
---|
493 | real(r8) qsp(pcols,pver) ! sat pt mixing ratio |
---|
494 | real(r8) qtl(pcols) ! tendency which would saturate the grid box in deltat |
---|
495 | real(r8) qtmp, ttmp ! work variable |
---|
496 | real(r8) relhum1(pcols) ! relative humidity |
---|
497 | real(r8) relhum(pcols) ! relative humidity |
---|
498 | !!$ real(r8) tc ! crit temp of transition to ice |
---|
499 | real(r8) t(pcols,pver) ! temp before time step ignoring condensate |
---|
500 | real(r8) tsp(pcols,pver) ! sat pt temperature |
---|
501 | real(r8) pol ! work variable |
---|
502 | real(r8) cdt ! work variable |
---|
503 | real(r8) wtthick ! work variable |
---|
504 | |
---|
505 | ! Extra local work space for cloud scheme modification |
---|
506 | |
---|
507 | real(r8) cpohl !Cp/Hlatv |
---|
508 | real(r8) hlocp !Hlatv/Cp |
---|
509 | real(r8) dto2 !0.5*deltat (delta=2.0*dt) |
---|
510 | real(r8) calpha(pcols) !alpha of new C - E scheme formulation |
---|
511 | real(r8) cbeta (pcols) !beta of new C - E scheme formulation |
---|
512 | real(r8) cbetah(pcols) !beta_hat at saturation portion |
---|
513 | real(r8) cgamma(pcols) !gamma of new C - E scheme formulation |
---|
514 | real(r8) cgamah(pcols) !gamma_hat at saturation portion |
---|
515 | real(r8) rcgama(pcols) !gamma/gamma_hat |
---|
516 | real(r8) csigma(pcols) !sigma of new C - E scheme formulation |
---|
517 | real(r8) cmec1 (pcols) !c1 of new C - E scheme formulation |
---|
518 | real(r8) cmec2 (pcols) !c2 of new C - E scheme formulation |
---|
519 | real(r8) cmec3 (pcols) !c3 of new C - E scheme formulation |
---|
520 | real(r8) cmec4 (pcols) !c4 of new C - E scheme formulation |
---|
521 | real(r8) cmeres(pcols) !residual cond of over-sat after cme and evapprec |
---|
522 | real(r8) ctmp !a scalar representation of cmeres |
---|
523 | real(r8) clrh2o ! Ratio of latvap to water vapor gas const |
---|
524 | real(r8) ice(pcols,pver) ! ice mixing ratio |
---|
525 | real(r8) liq(pcols,pver) ! liquid mixing ratio |
---|
526 | real(r8) rcwn(pcols,2,pver), rliq(pcols,2,pver), rice(pcols,2,pver) |
---|
527 | real(r8) cwnsave(pcols,2,pver), cmesave(pcols,2,pver) |
---|
528 | real(r8) prodprecsave(pcols,2,pver) |
---|
529 | logical error_found |
---|
530 | ! |
---|
531 | !------------------------------------------------------------ |
---|
532 | ! |
---|
533 | clrh2o = hlatv/rh2o ! Ratio of latvap to water vapor gas const |
---|
534 | omeps = 1.0_r8 - epsqs |
---|
535 | #ifdef PERGRO |
---|
536 | mincld = 1.e-4_r8 |
---|
537 | iter = 1 ! number of times to iterate the precipitation calculation |
---|
538 | #else |
---|
539 | mincld = 1.e-4_r8 |
---|
540 | iter = 2 |
---|
541 | #endif |
---|
542 | ! omsm = 0.99999 |
---|
543 | cpohl = cp/hlatv |
---|
544 | hlocp = hlatv/cp |
---|
545 | dto2=0.5_r8*deltat |
---|
546 | ! |
---|
547 | ! Constant for computing rate of evaporation of precipitation: |
---|
548 | ! |
---|
549 | !!$ conke = 1.e-5 |
---|
550 | !!$ conke = 1.e-6 |
---|
551 | ! |
---|
552 | ! initialize a few single level fields |
---|
553 | ! |
---|
554 | do i = 1,ncol |
---|
555 | precip(i) = 0.0_r8 |
---|
556 | precab(i) = 0.0_r8 |
---|
557 | snowab(i) = 0.0_r8 |
---|
558 | cldmax(i) = 0.0_r8 |
---|
559 | end do |
---|
560 | ! |
---|
561 | ! initialize multi-level fields |
---|
562 | ! |
---|
563 | do k = 1,pver |
---|
564 | do i = 1,ncol |
---|
565 | q(i,k) = qn(i,k) |
---|
566 | t(i,k) = tn(i,k) |
---|
567 | ! q(i,k)=qn(i,k)-qtend(i,k)*deltat |
---|
568 | ! t(i,k)=tn(i,k)-ttend(i,k)*deltat |
---|
569 | end do |
---|
570 | end do |
---|
571 | cme (:ncol,:) = 0._r8 |
---|
572 | evapprec(:ncol,:) = 0._r8 |
---|
573 | prodprec(:ncol,:) = 0._r8 |
---|
574 | evapsnow(:ncol,:) = 0._r8 |
---|
575 | prodsnow(:ncol,:) = 0._r8 |
---|
576 | evapheat(:ncol,:) = 0._r8 |
---|
577 | meltheat(:ncol,:) = 0._r8 |
---|
578 | prfzheat(:ncol,:) = 0._r8 |
---|
579 | ice2pr(:ncol,:) = 0._r8 |
---|
580 | liq2pr(:ncol,:) = 0._r8 |
---|
581 | liq2snow(:ncol,:) = 0._r8 |
---|
582 | fwaut(:ncol,:) = 0._r8 |
---|
583 | fsaut(:ncol,:) = 0._r8 |
---|
584 | fracw(:ncol,:) = 0._r8 |
---|
585 | fsacw(:ncol,:) = 0._r8 |
---|
586 | fsaci(:ncol,:) = 0._r8 |
---|
587 | rkflxprc(:ncol,:) = 0._r8 |
---|
588 | rkflxsnw(:ncol,:) = 0._r8 |
---|
589 | |
---|
590 | ! |
---|
591 | ! find the wet bulb temp and saturation value |
---|
592 | ! for the provisional t and q without condensation |
---|
593 | ! |
---|
594 | call findsp (lchnk, ncol, qn, tn, p, tsp, qsp) |
---|
595 | do 800 k = k1mb,pver |
---|
596 | call vqsatd (t(1,k), p(1,k), es, qs, gamma, ncol) |
---|
597 | do i = 1,ncol |
---|
598 | relhum(i) = q(i,k)/qs(i) |
---|
599 | ! |
---|
600 | cldm(i) = max(cldn(i,k),mincld) |
---|
601 | ! |
---|
602 | ! the max cloud fraction above this level |
---|
603 | ! |
---|
604 | cldmax(i) = max(cldmax(i), cldm(i)) |
---|
605 | |
---|
606 | ! define the coefficients for C - E calculation |
---|
607 | |
---|
608 | calpha(i) = 1.0_r8/qs(i) |
---|
609 | cbeta (i) = q(i,k)/qs(i)**2*gamma(i)*cpohl |
---|
610 | cbetah(i) = 1.0_r8/qs(i)*gamma(i)*cpohl |
---|
611 | cgamma(i) = calpha(i)+hlatv*cbeta(i)/cp |
---|
612 | cgamah(i) = calpha(i)+hlatv*cbetah(i)/cp |
---|
613 | rcgama(i) = cgamma(i)/cgamah(i) |
---|
614 | |
---|
615 | if(cldm(i) > mincld) then |
---|
616 | icwc(i) = max(0._r8,cwat(i,k)/cldm(i)) |
---|
617 | else |
---|
618 | icwc(i) = 0.0_r8 |
---|
619 | endif |
---|
620 | !PJR the above logic give zero icwc with nonzero cwat, dont like it! |
---|
621 | !PJR generates problems with csigma |
---|
622 | !PJR set the icwc to a very small number, so we can start from zero cloud cover and make some clouds |
---|
623 | ! icwc(i) = max(1.e-8_r8,cwat(i,k)/cldm(i)) |
---|
624 | |
---|
625 | ! |
---|
626 | ! initial guess of evaporation, will be updated within iteration |
---|
627 | ! |
---|
628 | evapprec(i,k) = conke*(1._r8 - cldm(i))*sqrt(precab(i)) & |
---|
629 | *(1._r8 - min(relhum(i),1._r8)) |
---|
630 | |
---|
631 | ! |
---|
632 | ! zero cmeres before iteration for each level |
---|
633 | ! |
---|
634 | cmeres(i)=0.0_r8 |
---|
635 | |
---|
636 | end do |
---|
637 | do i = 1,ncol |
---|
638 | ! |
---|
639 | ! fractions of ice at this level |
---|
640 | ! |
---|
641 | !!$ tc = t(i,k) - t0 |
---|
642 | !!$ fice(i,k) = max(0._r8,min(-tc*0.05,1.0_r8)) |
---|
643 | ! |
---|
644 | ! calculate the cooling due to a phase change of the rainwater |
---|
645 | ! from above |
---|
646 | ! |
---|
647 | if (t(i,k) >= t0) then |
---|
648 | meltheat(i,k) = -hlatf * snowab(i) * gravit/pdel(i,k) |
---|
649 | snowab(i) = 0._r8 |
---|
650 | else |
---|
651 | meltheat(i,k) = 0._r8 |
---|
652 | endif |
---|
653 | end do |
---|
654 | |
---|
655 | ! |
---|
656 | ! calculate cme and formation of precip. |
---|
657 | ! |
---|
658 | ! The cloud microphysics is highly nonlinear and coupled with cme |
---|
659 | ! Both rain processes and cme are calculated iteratively. |
---|
660 | ! |
---|
661 | do 100 l = 1,iter |
---|
662 | |
---|
663 | do i = 1,ncol |
---|
664 | |
---|
665 | ! |
---|
666 | ! calculation of cme has 4 scenarios |
---|
667 | ! ================================== |
---|
668 | ! |
---|
669 | if(relhum(i) > rhu00(i,k)) then |
---|
670 | |
---|
671 | ! 1. whole grid saturation |
---|
672 | ! ======================== |
---|
673 | if(relhum(i) >= 0.999_r8 .or. cldm(i) >= 0.999_r8 ) then |
---|
674 | cme(i,k)=(calpha(i)*qtend(i,k)-cbetah(i)*ttend(i,k))/cgamah(i) |
---|
675 | |
---|
676 | ! 2. fractional saturation |
---|
677 | ! ======================== |
---|
678 | else |
---|
679 | if (rhdfda(i,k) .eq. 0. .and. icwc(i).eq.0.) then |
---|
680 | write (iulog,*) ' cldwat.F90: empty rh cloud ', i, k, lchnk |
---|
681 | write (iulog,*) ' relhum, iter ', relhum(i), l, rhu00(i,k), cldm(i), cldn(i,k) |
---|
682 | call endrun () |
---|
683 | endif |
---|
684 | csigma(i) = 1.0_r8/(rhdfda(i,k)+cgamma(i)*icwc(i)) |
---|
685 | cmec1(i) = (1.0_r8-cldm(i))*csigma(i)*rhdfda(i,k) |
---|
686 | cmec2(i) = cldm(i)*calpha(i)/cgamah(i)+(1.0_r8-rcgama(i)*cldm(i))* & |
---|
687 | csigma(i)*calpha(i)*icwc(i) |
---|
688 | cmec3(i) = cldm(i)*cbetah(i)/cgamah(i) + & |
---|
689 | (cbeta(i)-rcgama(i)*cldm(i)*cbetah(i))*csigma(i)*icwc(i) |
---|
690 | cmec4(i) = csigma(i)*cgamma(i)*icwc(i) |
---|
691 | |
---|
692 | ! Q=C-E=-C1*Al + C2*Aq - C3* At + C4*Er |
---|
693 | |
---|
694 | cme(i,k) = -cmec1(i)*lctend(i,k) + cmec2(i)*qtend(i,k) & |
---|
695 | -cmec3(i)*ttend(i,k) + cmec4(i)*evapprec(i,k) |
---|
696 | endif |
---|
697 | |
---|
698 | ! 3. when rh < rhu00, evaporate existing cloud water |
---|
699 | ! ================================================== |
---|
700 | else if(cwat(i,k) > 0.0_r8)then |
---|
701 | ! liquid water should be evaporated but not to exceed |
---|
702 | ! saturation point. if qn > qsp, not to evaporate cwat |
---|
703 | cme(i,k)=-min(max(0._r8,qsp(i,k)-qn(i,k)),cwat(i,k))/deltat |
---|
704 | |
---|
705 | ! 4. no condensation nor evaporation |
---|
706 | ! ================================== |
---|
707 | else |
---|
708 | cme(i,k)=0.0_r8 |
---|
709 | endif |
---|
710 | |
---|
711 | |
---|
712 | end do !end loop for cme update |
---|
713 | |
---|
714 | ! Because of the finite time step, |
---|
715 | ! place a bound here not to exceed wet bulb point |
---|
716 | ! and not to evaporate more than available water |
---|
717 | ! |
---|
718 | do i = 1, ncol |
---|
719 | qtmp = qn(i,k) - cme(i,k)*deltat |
---|
720 | |
---|
721 | ! possibilities to have qtmp > qsp |
---|
722 | ! |
---|
723 | ! 1. if qn > qs(tn), it condenses; |
---|
724 | ! if after applying cme, qtmp > qsp, more condensation is applied. |
---|
725 | ! |
---|
726 | ! 2. if qn < qs, evaporation should not exceed qsp, |
---|
727 | |
---|
728 | if(qtmp > qsp(i,k)) then |
---|
729 | cme(i,k) = cme(i,k) + (qtmp-qsp(i,k))/deltat |
---|
730 | endif |
---|
731 | |
---|
732 | ! |
---|
733 | ! if net evaporation, it should not exceed available cwat |
---|
734 | ! |
---|
735 | if(cme(i,k) < -cwat(i,k)/deltat) & |
---|
736 | cme(i,k) = -cwat(i,k)/deltat |
---|
737 | ! |
---|
738 | ! addition of residual condensation from previous step of iteration |
---|
739 | ! |
---|
740 | cme(i,k) = cme(i,k) + cmeres(i) |
---|
741 | |
---|
742 | end do |
---|
743 | |
---|
744 | ! limit cme for roundoff errors |
---|
745 | do i = 1, ncol |
---|
746 | cme(i,k) = cme(i,k)*omsm |
---|
747 | end do |
---|
748 | |
---|
749 | do i = 1,ncol |
---|
750 | ! |
---|
751 | ! as a safe limit, condensation should not reduce grid mean rh below rhu00 |
---|
752 | ! |
---|
753 | if(cme(i,k) > 0.0_r8 .and. relhum(i) > rhu00(i,k) ) & |
---|
754 | cme(i,k) = min(cme(i,k), (qn(i,k)-qs(i)*rhu00(i,k))/deltat) |
---|
755 | ! |
---|
756 | ! initial guess for cwm (mean cloud water over time step) if 1st iteration |
---|
757 | ! |
---|
758 | if(l < 2) then |
---|
759 | cwm(i) = max(cwat(i,k)+cme(i,k)*dto2, 0._r8) |
---|
760 | endif |
---|
761 | |
---|
762 | enddo |
---|
763 | |
---|
764 | ! provisional precipitation falling through model layer |
---|
765 | do i = 1,ncol |
---|
766 | !!$ prprov(i) = precab(i) + prodprec(i,k)*pdel(i,k)/gravit |
---|
767 | ! rain produced in this layer not too effective in collection process |
---|
768 | wtthick = max(0._r8,min(0.5_r8,((zi(i,k)-zi(i,k+1))/1000._r8)**2)) |
---|
769 | prprov(i) = precab(i) + wtthick*prodprec(i,k)*pdel(i,k)/gravit |
---|
770 | end do |
---|
771 | |
---|
772 | ! calculate conversion of condensate to precipitation by cloud microphysics |
---|
773 | call findmcnew (lchnk ,ncol , & |
---|
774 | k ,prprov ,snowab, t ,p , & |
---|
775 | cwm ,cldm ,cldmax ,fice(1,k),coef , & |
---|
776 | fwaut(1,k),fsaut(1,k),fracw(1,k),fsacw(1,k),fsaci(1,k), & |
---|
777 | seaicef, snowh) |
---|
778 | ! |
---|
779 | ! calculate the precip rate |
---|
780 | ! |
---|
781 | error_found = .false. |
---|
782 | do i = 1,ncol |
---|
783 | if (cldm(i) > 0) then |
---|
784 | ! |
---|
785 | ! first predict the cloud water |
---|
786 | ! |
---|
787 | cdt = coef(i)*deltat |
---|
788 | if (cdt > 0.01_r8) then |
---|
789 | pol = cme(i,k)/coef(i) ! production over loss |
---|
790 | cwn(i) = max(0._r8,(cwat(i,k)-pol)*exp(-cdt)+ pol) |
---|
791 | else |
---|
792 | cwn(i) = max(0._r8,(cwat(i,k) + cme(i,k)*deltat)/(1+cdt)) |
---|
793 | endif |
---|
794 | ! |
---|
795 | ! now back out the tendency of net rain production |
---|
796 | ! |
---|
797 | prodprec(i,k) = max(0._r8,cme(i,k)-(cwn(i)-cwat(i,k))/deltat) |
---|
798 | else |
---|
799 | prodprec(i,k) = 0.0_r8 |
---|
800 | cwn(i) = 0._r8 |
---|
801 | endif |
---|
802 | |
---|
803 | ! provisional calculation of conversion terms |
---|
804 | ice2pr(i,k) = prodprec(i,k)*(fsaut(i,k)+fsaci(i,k)) |
---|
805 | liq2pr(i,k) = prodprec(i,k)*(fwaut(i,k)+fsacw(i,k)+fracw(i,k)) |
---|
806 | !old liq2snow(i,k) = prodprec(i,k)*fsacw(i,k) |
---|
807 | |
---|
808 | ! revision suggested by Jim McCaa |
---|
809 | ! it controls the amount of snow hitting the sfc |
---|
810 | ! by forcing a lot of conversion of cloud liquid to snow phase |
---|
811 | ! it might be better done later by an explicit representation of |
---|
812 | ! rain accreting ice (and freezing), or by an explicit freezing of raindrops |
---|
813 | liq2snow(i,k) = max(prodprec(i,k)*fsacw(i,k), fsnow(i,k)*liq2pr(i,k)) |
---|
814 | |
---|
815 | ! bounds |
---|
816 | nice2pr = min(ice2pr(i,k),(cwat(i,k)+cme(i,k)*deltat)*fice(i,k)/deltat) |
---|
817 | nliq2pr = min(liq2pr(i,k),(cwat(i,k)+cme(i,k)*deltat)*(1._r8-fice(i,k))/deltat) |
---|
818 | ! write(iulog,*) ' prodprec ', i, k, prodprec(i,k) |
---|
819 | ! write(iulog,*) ' nliq2pr, nice2pr ', nliq2pr, nice2pr |
---|
820 | if (liq2pr(i,k).ne.0._r8) then |
---|
821 | nliq2snow = liq2snow(i,k)*nliq2pr/liq2pr(i,k) ! correction |
---|
822 | else |
---|
823 | nliq2snow = liq2snow(i,k) |
---|
824 | endif |
---|
825 | |
---|
826 | ! avoid roundoff problems generating negatives |
---|
827 | nliq2snow = nliq2snow*omsm |
---|
828 | nliq2pr = nliq2pr*omsm |
---|
829 | nice2pr = nice2pr*omsm |
---|
830 | |
---|
831 | ! final estimates of conversion to precip and snow |
---|
832 | prodprec(i,k) = (nliq2pr + nice2pr) |
---|
833 | prodsnow(i,k) = (nice2pr + nliq2snow) |
---|
834 | |
---|
835 | rcwn(i,l,k) = cwat(i,k) + (cme(i,k)- prodprec(i,k))*deltat |
---|
836 | rliq(i,l,k) = (cwat(i,k) + cme(i,k)*deltat)*(1._r8-fice(i,k)) - nliq2pr * deltat |
---|
837 | rice(i,l,k) = (cwat(i,k) + cme(i,k)*deltat)* fice(i,k) - nice2pr *deltat |
---|
838 | |
---|
839 | ! Save for sanity check later... |
---|
840 | ! Putting sanity checks inside loops 100 and 800 screws up the |
---|
841 | ! IBM compiler for reasons as yet unknown. TBH |
---|
842 | cwnsave(i,l,k) = cwn(i) |
---|
843 | cmesave(i,l,k) = cme(i,k) |
---|
844 | prodprecsave(i,l,k) = prodprec(i,k) |
---|
845 | ! End of save for sanity check later... |
---|
846 | |
---|
847 | ! final version of condensate to precip terms |
---|
848 | liq2pr(i,k) = nliq2pr |
---|
849 | liq2snow(i,k) = nliq2snow |
---|
850 | ice2pr(i,k) = nice2pr |
---|
851 | |
---|
852 | cwn(i) = rcwn(i,l,k) |
---|
853 | ! |
---|
854 | ! update any remaining provisional values |
---|
855 | ! |
---|
856 | cwm(i) = (cwn(i) + cwat(i,k))*0.5_r8 |
---|
857 | ! |
---|
858 | ! update in cloud water |
---|
859 | ! |
---|
860 | if(cldm(i) > mincld) then |
---|
861 | icwc(i) = cwm(i)/cldm(i) |
---|
862 | else |
---|
863 | icwc(i) = 0.0_r8 |
---|
864 | endif |
---|
865 | !PJR the above logic give zero icwc with nonzero cwat, dont like it! |
---|
866 | !PJR generates problems with csigma |
---|
867 | !PJR set the icwc to a very small number, so we can start from zero cloud cover and make some clouds |
---|
868 | ! icwc(i) = max(1.e-8_r8,cwm(i)/cldm(i)) |
---|
869 | |
---|
870 | end do ! end of do i = 1,ncol |
---|
871 | |
---|
872 | ! |
---|
873 | ! calculate provisional value of cloud water for |
---|
874 | ! evaporation of precipitate (evapprec) calculation |
---|
875 | ! |
---|
876 | do i = 1,ncol |
---|
877 | qtmp = qn(i,k) - cme(i,k)*deltat |
---|
878 | ttmp = tn(i,k) + deltat/cp * ( meltheat(i,k) & |
---|
879 | + (hlatv + hlatf*fice(i,k)) * cme(i,k) ) |
---|
880 | esn = estblf(ttmp) |
---|
881 | qsn = min(epsqs*esn/(p(i,k) - omeps*esn),1._r8) |
---|
882 | qtl(i) = max((qsn - qtmp)/deltat,0._r8) |
---|
883 | relhum1(i) = qtmp/qsn |
---|
884 | end do |
---|
885 | ! |
---|
886 | do i = 1,ncol |
---|
887 | #ifdef PERGRO |
---|
888 | evapprec(i,k) = conke*(1._r8 - max(cldm(i),mincld))* & |
---|
889 | sqrt(precab(i))*(1._r8 - min(relhum1(i),1._r8)) |
---|
890 | #else |
---|
891 | evapprec(i,k) = conke*(1._r8 - cldm(i))*sqrt(precab(i)) & |
---|
892 | *(1._r8 - min(relhum1(i),1._r8)) |
---|
893 | #endif |
---|
894 | ! |
---|
895 | ! limit the evaporation to the amount which is entering the box |
---|
896 | ! or saturates the box |
---|
897 | ! |
---|
898 | prtmp = precab(i)*gravit/pdel(i,k) |
---|
899 | evapprec(i,k) = min(evapprec(i,k), prtmp, qtl(i))*omsm |
---|
900 | #ifdef PERGRO |
---|
901 | ! zeroing needed for pert growth |
---|
902 | evapprec(i,k) = 0._r8 |
---|
903 | #endif |
---|
904 | ! |
---|
905 | ! Partition evaporation of precipitate between rain and snow using |
---|
906 | ! the fraction of snow falling into the box. Determine the heating |
---|
907 | ! due to evaporation. Note that evaporation is positive (loss of precip, |
---|
908 | ! gain of vapor) and that heating is negative. |
---|
909 | if (evapprec(i,k) > 0._r8) then |
---|
910 | evapsnow(i,k) = evapprec(i,k) * snowab(i) / precab(i) |
---|
911 | evapheat(i,k) = -hlatv * evapprec(i,k) - hlatf * evapsnow(i,k) |
---|
912 | else |
---|
913 | evapsnow(i,k) = 0._r8 |
---|
914 | evapheat(i,k) = 0._r8 |
---|
915 | end if |
---|
916 | ! Account for the latent heat of fusion for liquid drops collected by falling snow |
---|
917 | prfzheat(i,k) = hlatf * liq2snow(i,k) |
---|
918 | end do |
---|
919 | |
---|
920 | ! now remove the residual of any over-saturation. Normally, |
---|
921 | ! the oversaturated water vapor should have been removed by |
---|
922 | ! cme formulation plus constraints by wet bulb tsp/qsp |
---|
923 | ! as computed above. However, because of non-linearity, |
---|
924 | ! addition of (cme-evapprec) to update t and q may still cause |
---|
925 | ! a very small amount of over saturation. It is called a |
---|
926 | ! residual of over-saturation because theoretically, cme |
---|
927 | ! should have taken care of all of large scale condensation. |
---|
928 | ! |
---|
929 | |
---|
930 | do i = 1,ncol |
---|
931 | qtmp = qn(i,k)-(cme(i,k)-evapprec(i,k))*deltat |
---|
932 | ttmp = tn(i,k) + deltat/cp * ( meltheat(i,k) + evapheat(i,k) + prfzheat(i,k) & |
---|
933 | + (hlatv + hlatf*fice(i,k)) * cme(i,k) ) |
---|
934 | esn = estblf(ttmp) |
---|
935 | qsn = min(epsqs*esn/(p(i,k) - omeps*esn),1._r8) |
---|
936 | ! |
---|
937 | !Upper stratosphere and mesosphere, qsn calculated |
---|
938 | !above may be negative. Here just to skip it instead |
---|
939 | !of resetting it to 1 as in aqsat |
---|
940 | ! |
---|
941 | if(qtmp > qsn .and. qsn > 0) then |
---|
942 | !calculate dqsdt, a more precise calculation |
---|
943 | !which taking into account different range of T |
---|
944 | !can be found in aqsatd.F. Here follows |
---|
945 | !cond.F to calculate it. |
---|
946 | ! |
---|
947 | denom = (p(i,k)-omeps*esn)*ttmp*ttmp |
---|
948 | dqsdt = clrh2o*qsn*p(i,k)/denom |
---|
949 | ! |
---|
950 | !now extra condensation to bring air to just saturation |
---|
951 | ! |
---|
952 | ctmp = (qtmp-qsn)/(1._r8+hlocp*dqsdt)/deltat |
---|
953 | cme(i,k) = cme(i,k)+ctmp |
---|
954 | ! |
---|
955 | ! save residual on cmeres to addtion to cme on entering next iteration |
---|
956 | ! cme exit here contain the residual but overrided if back to iteration |
---|
957 | ! |
---|
958 | cmeres(i) = ctmp |
---|
959 | else |
---|
960 | cmeres(i) = 0.0_r8 |
---|
961 | endif |
---|
962 | end do |
---|
963 | |
---|
964 | 100 continue ! end of do l = 1,iter |
---|
965 | |
---|
966 | ! |
---|
967 | ! precipitation |
---|
968 | ! |
---|
969 | do i = 1,ncol |
---|
970 | precip(i) = precip(i) + pdel(i,k)/gravit * (prodprec(i,k) - evapprec(i,k)) |
---|
971 | precab(i) = precab(i) + pdel(i,k)/gravit * (prodprec(i,k) - evapprec(i,k)) |
---|
972 | if(precab(i).lt.0._r8) precab(i)=0._r8 |
---|
973 | ! snowab(i) = snowab(i) + pdel(i,k)/gravit * (prodprec(i,k)*fice(i,k) - evapsnow(i,k)) |
---|
974 | snowab(i) = snowab(i) + pdel(i,k)/gravit * (prodsnow(i,k) - evapsnow(i,k)) |
---|
975 | |
---|
976 | ! If temperature above freezing, all precip is rain flux. if temperature below freezing, all precip is snow flux. |
---|
977 | rkflxprc(i,k+1) = precab(i) - snowab(i) |
---|
978 | rkflxsnw(i,k+1) = snowab(i) |
---|
979 | |
---|
980 | !!$ if ((precab(i)) < 1.e-10) then |
---|
981 | !!$ precab(i) = 0. |
---|
982 | !!$ snowab(i) = 0. |
---|
983 | !!$ endif |
---|
984 | end do |
---|
985 | 800 continue ! level loop (k=1,pver) |
---|
986 | |
---|
987 | ! begin sanity checks |
---|
988 | error_found = .false. |
---|
989 | do k = k1mb,pver |
---|
990 | do l = 1,iter |
---|
991 | do i = 1,ncol |
---|
992 | if (abs(rcwn(i,l,k)).lt.1.e-300_r8) rcwn(i,l,k) = 0._r8 |
---|
993 | if (abs(rliq(i,l,k)).lt.1.e-300_r8) rliq(i,l,k) = 0._r8 |
---|
994 | if (abs(rice(i,l,k)).lt.1.e-300_r8) rice(i,l,k) = 0._r8 |
---|
995 | if (rcwn(i,l,k).lt.0._r8) error_found = .true. |
---|
996 | if (rliq(i,l,k).lt.0._r8) error_found = .true. |
---|
997 | if (rice(i,l,k).lt.0._r8) error_found = .true. |
---|
998 | enddo |
---|
999 | enddo |
---|
1000 | enddo |
---|
1001 | if (error_found) then |
---|
1002 | do k = k1mb,pver |
---|
1003 | do l = 1,iter |
---|
1004 | do i = 1,ncol |
---|
1005 | if (rcwn(i,l,k).lt.0._r8) then |
---|
1006 | write(iulog,*) ' prob with neg rcwn1 ', rcwn(i,l,k), & |
---|
1007 | cwnsave(i,l,k) |
---|
1008 | write(iulog,*) ' cwat, cme*deltat, prodprec*deltat ', & |
---|
1009 | cwat(i,k), cmesave(i,l,k)*deltat, & |
---|
1010 | prodprecsave(i,l,k)*deltat, & |
---|
1011 | (cmesave(i,l,k)-prodprecsave(i,l,k))*deltat |
---|
1012 | call endrun('PCOND') |
---|
1013 | endif |
---|
1014 | if (rliq(i,l,k).lt.0._r8) then |
---|
1015 | write(iulog,*) ' prob with neg rliq1 ', rliq(i,l,k) |
---|
1016 | call endrun('PCOND') |
---|
1017 | endif |
---|
1018 | if (rice(i,l,k).lt.0._r8) then |
---|
1019 | write(iulog,*) ' prob with neg rice ', rice(i,l,k) |
---|
1020 | call endrun('PCOND') |
---|
1021 | endif |
---|
1022 | enddo |
---|
1023 | enddo |
---|
1024 | enddo |
---|
1025 | end if |
---|
1026 | ! end sanity checks |
---|
1027 | |
---|
1028 | return |
---|
1029 | end subroutine pcond |
---|
1030 | |
---|
1031 | !############################################################################## |
---|
1032 | |
---|
1033 | subroutine findmcnew (lchnk ,ncol , & |
---|
1034 | k ,precab ,snowab, t ,p , & |
---|
1035 | cwm ,cldm ,cldmax ,fice ,coef , & |
---|
1036 | fwaut ,fsaut ,fracw ,fsacw ,fsaci , & |
---|
1037 | seaicef ,snowh ) |
---|
1038 | !----------------------------------------------------------------------- |
---|
1039 | ! |
---|
1040 | ! Purpose: |
---|
1041 | ! calculate the conversion of condensate to precipitate |
---|
1042 | ! |
---|
1043 | ! Method: |
---|
1044 | ! See: Rasch, P. J, and J. E. Kristjansson, A Comparison of the CCM3 |
---|
1045 | ! model climate using diagnosed and |
---|
1046 | ! predicted condensate parameterizations, 1998, J. Clim., 11, |
---|
1047 | ! pp1587---1614. |
---|
1048 | ! |
---|
1049 | ! Author: P. Rasch |
---|
1050 | ! |
---|
1051 | !----------------------------------------------------------------------- |
---|
1052 | use phys_grid, only: get_rlat_all_p |
---|
1053 | use comsrf, only: landm |
---|
1054 | ! |
---|
1055 | ! input args |
---|
1056 | ! |
---|
1057 | integer, intent(in) :: lchnk ! chunk identifier |
---|
1058 | integer, intent(in) :: ncol ! number of atmospheric columns |
---|
1059 | integer, intent(in) :: k ! level index |
---|
1060 | |
---|
1061 | real(r8), intent(in) :: precab(pcols) ! rate of precipitation from above (kg / (m**2 * s)) |
---|
1062 | real(r8), intent(in) :: t(pcols,pver) ! temperature (K) |
---|
1063 | real(r8), intent(in) :: p(pcols,pver) ! pressure (Pa) |
---|
1064 | real(r8), intent(in) :: cldm(pcols) ! cloud fraction |
---|
1065 | real(r8), intent(in) :: cldmax(pcols) ! max cloud fraction above this level |
---|
1066 | real(r8), intent(in) :: cwm(pcols) ! condensate mixing ratio (kg/kg) |
---|
1067 | real(r8), intent(in) :: fice(pcols) ! fraction of cwat that is ice |
---|
1068 | real(r8), intent(in) :: seaicef(pcols) ! sea ice fraction |
---|
1069 | real(r8), intent(in) :: snowab(pcols) ! rate of snow from above (kg / (m**2 * s)) |
---|
1070 | real(r8), intent(in) :: snowh(pcols) ! Snow depth over land, water equivalent (m) |
---|
1071 | |
---|
1072 | ! output arguments |
---|
1073 | real(r8), intent(out) :: coef(pcols) ! conversion rate (1/s) |
---|
1074 | real(r8), intent(out) :: fwaut(pcols) ! relative importance of liquid autoconversion (a diagnostic) |
---|
1075 | real(r8), intent(out) :: fsaut(pcols) ! relative importance of ice autoconversion (a diagnostic) |
---|
1076 | real(r8), intent(out) :: fracw(pcols) ! relative importance of rain accreting liquid (a diagnostic) |
---|
1077 | real(r8), intent(out) :: fsacw(pcols) ! relative importance of snow accreting liquid (a diagnostic) |
---|
1078 | real(r8), intent(out) :: fsaci(pcols) ! relative importance of snow accreting ice (a diagnostic) |
---|
1079 | |
---|
1080 | ! work variables |
---|
1081 | |
---|
1082 | integer i |
---|
1083 | integer ii |
---|
1084 | integer ind(pcols) |
---|
1085 | integer ncols |
---|
1086 | |
---|
1087 | real(r8), parameter :: degrad = 57.296_r8 ! divide by this to convert degrees to radians |
---|
1088 | real(r8) capn ! local cloud particles / cm3 |
---|
1089 | real(r8) capnoice ! local cloud particles when not over sea ice / cm3 |
---|
1090 | real(r8) ciaut ! coefficient of autoconversion of ice (1/s) |
---|
1091 | real(r8) cldloc(pcols) ! non-zero amount of cloud |
---|
1092 | real(r8) cldpr(pcols) ! assumed cloudy volume occupied by rain and cloud |
---|
1093 | real(r8) con1 ! work constant |
---|
1094 | real(r8) con2 ! work constant |
---|
1095 | real(r8) csacx ! constant used for snow accreting liquid or ice |
---|
1096 | !!$ real(r8) dtice ! interval for transition from liquid to ice |
---|
1097 | real(r8) icemr(pcols) ! in-cloud ice mixing ratio |
---|
1098 | real(r8) icrit ! threshold for autoconversion of ice |
---|
1099 | real(r8) liqmr(pcols) ! in-cloud liquid water mixing ratio |
---|
1100 | real(r8) pracw ! rate of rain accreting water |
---|
1101 | real(r8) prlloc(pcols) ! local rain flux in mm/day |
---|
1102 | real(r8) prscgs(pcols) ! local snow amount in cgs units |
---|
1103 | real(r8) psaci ! rate of collection of ice by snow (lin et al 1983) |
---|
1104 | real(r8) psacw ! rate of collection of liquid by snow (lin et al 1983) |
---|
1105 | real(r8) psaut ! rate of autoconversion of ice condensate |
---|
1106 | real(r8) ptot ! total rate of conversion |
---|
1107 | real(r8) pwaut ! rate of autoconversion of liquid condensate |
---|
1108 | real(r8) r3l ! volume radius |
---|
1109 | real(r8) rainmr(pcols) ! in-cloud rain mixing ratio |
---|
1110 | real(r8) rat1 ! work constant |
---|
1111 | real(r8) rat2 ! work constant |
---|
1112 | !!$ real(r8) rdtice ! recipricol of dtice |
---|
1113 | real(r8) rho(pcols) ! density (mks units) |
---|
1114 | real(r8) rhocgs ! density (cgs units) |
---|
1115 | real(r8) rlat(pcols) ! latitude (radians) |
---|
1116 | real(r8) snowfr ! fraction of precipate existing as snow |
---|
1117 | real(r8) totmr(pcols) ! in-cloud total condensate mixing ratio |
---|
1118 | real(r8) vfallw ! fall speed of precipitate as liquid |
---|
1119 | real(r8) wp ! weight factor used in calculating pressure dep of autoconversion |
---|
1120 | real(r8) wsi ! weight factor for sea ice |
---|
1121 | real(r8) wt ! fraction of ice |
---|
1122 | real(r8) wland ! fraction of land |
---|
1123 | |
---|
1124 | ! real(r8) csaci |
---|
1125 | ! real(r8) csacw |
---|
1126 | ! real(r8) cwaut |
---|
1127 | ! real(r8) efact |
---|
1128 | ! real(r8) lamdas |
---|
1129 | ! real(r8) lcrit |
---|
1130 | ! real(r8) rcwm |
---|
1131 | ! real(r8) r3lc2 |
---|
1132 | ! real(r8) snowmr(pcols) |
---|
1133 | ! real(r8) vfalls |
---|
1134 | |
---|
1135 | real(8) ftot |
---|
1136 | |
---|
1137 | ! inline statement functions |
---|
1138 | real(r8) heavy, heavym, a1, a2, heavyp, heavymp |
---|
1139 | heavy(a1,a2) = max(0._r8,sign(1._r8,a1-a2)) ! heavyside function |
---|
1140 | heavym(a1,a2) = max(0.01_r8,sign(1._r8,a1-a2)) ! modified heavyside function |
---|
1141 | ! |
---|
1142 | ! New heavyside functions to perhaps address error growth problems |
---|
1143 | ! |
---|
1144 | heavyp(a1,a2) = a1/(a2+a1+1.e-36_r8) |
---|
1145 | heavymp(a1,a2) = (a1+0.01_r8*a2)/(a2+a1+1.e-36_r8) |
---|
1146 | |
---|
1147 | ! |
---|
1148 | ! find all the points where we need to do the microphysics |
---|
1149 | ! and set the output variables to zero |
---|
1150 | ! |
---|
1151 | ncols = 0 |
---|
1152 | do i = 1,ncol |
---|
1153 | coef(i) = 0._r8 |
---|
1154 | fwaut(i) = 0._r8 |
---|
1155 | fsaut(i) = 0._r8 |
---|
1156 | fracw(i) = 0._r8 |
---|
1157 | fsacw(i) = 0._r8 |
---|
1158 | fsaci(i) = 0._r8 |
---|
1159 | liqmr(i) = 0._r8 |
---|
1160 | rainmr(i) = 0._r8 |
---|
1161 | if (cwm(i) > 1.e-20_r8) then |
---|
1162 | ncols = ncols + 1 |
---|
1163 | ind(ncols) = i |
---|
1164 | endif |
---|
1165 | end do |
---|
1166 | |
---|
1167 | !cdir nodep |
---|
1168 | !DIR$ CONCURRENT |
---|
1169 | do ii = 1,ncols |
---|
1170 | i = ind(ii) |
---|
1171 | ! |
---|
1172 | ! the local cloudiness at this level |
---|
1173 | ! |
---|
1174 | cldloc(i) = max(cldmin,cldm(i)) |
---|
1175 | ! |
---|
1176 | ! a weighted mean between max cloudiness above, and this layer |
---|
1177 | ! |
---|
1178 | cldpr(i) = max(cldmin,(cldmax(i)+cldm(i))*0.5_r8) |
---|
1179 | ! |
---|
1180 | ! decompose the suspended condensate into |
---|
1181 | ! an incloud liquid and ice phase component |
---|
1182 | ! |
---|
1183 | totmr(i) = cwm(i)/cldloc(i) |
---|
1184 | icemr(i) = totmr(i)*fice(i) |
---|
1185 | liqmr(i) = totmr(i)*(1._r8-fice(i)) |
---|
1186 | ! |
---|
1187 | ! density |
---|
1188 | ! |
---|
1189 | rho(i) = p(i,k)/(287._r8*t(i,k)) |
---|
1190 | rhocgs = rho(i)*1.e-3_r8 ! density in cgs units |
---|
1191 | ! |
---|
1192 | ! decompose the precipitate into a liquid and ice phase |
---|
1193 | ! |
---|
1194 | if (t(i,k) > t0) then |
---|
1195 | vfallw = convfw/sqrt(rho(i)) |
---|
1196 | rainmr(i) = precab(i)/(rho(i)*vfallw*cldpr(i)) |
---|
1197 | snowfr = 0 |
---|
1198 | ! snowmr(i) |
---|
1199 | else |
---|
1200 | snowfr = 1 |
---|
1201 | rainmr(i) = 0._r8 |
---|
1202 | endif |
---|
1203 | ! rainmr(i) = (precab(i)-snowab(i))/(rho(i)*vfallw*cldpr(i)) |
---|
1204 | ! |
---|
1205 | ! local snow amount in cgs units |
---|
1206 | ! |
---|
1207 | prscgs(i) = precab(i)/cldpr(i)*0.1_r8*snowfr |
---|
1208 | ! prscgs(i) = snowab(i)/cldpr(i)*0.1 |
---|
1209 | ! |
---|
1210 | ! local rain amount in mm/day |
---|
1211 | ! |
---|
1212 | prlloc(i) = precab(i)*86400._r8/cldpr(i) |
---|
1213 | end do |
---|
1214 | |
---|
1215 | con1 = 1._r8/(1.333_r8*pi)**0.333_r8 * 0.01_r8 ! meters |
---|
1216 | ! |
---|
1217 | ! calculate the conversion terms |
---|
1218 | ! |
---|
1219 | call get_rlat_all_p(lchnk, ncol, rlat) |
---|
1220 | |
---|
1221 | !cdir nodep |
---|
1222 | !DIR$ CONCURRENT |
---|
1223 | do ii = 1,ncols |
---|
1224 | i = ind(ii) |
---|
1225 | rhocgs = rho(i)*1.e-3_r8 ! density in cgs units |
---|
1226 | ! |
---|
1227 | ! exponential temperature factor |
---|
1228 | ! |
---|
1229 | ! efact = exp(0.025*(t(i,k)-t0)) |
---|
1230 | ! |
---|
1231 | ! some temperature dependent constants |
---|
1232 | ! |
---|
1233 | !!$ wt = min(1._r8,max(0._r8,(t0-t(i,k))*rdtice)) |
---|
1234 | wt = fice(i) |
---|
1235 | icrit = icritc*wt + icritw*(1-wt) |
---|
1236 | ! |
---|
1237 | ! jrm Reworked droplet number concentration algorithm |
---|
1238 | ! Start with pressure-dependent value appropriate for continental air |
---|
1239 | ! Note: reltab has a temperature dependence here |
---|
1240 | capn = capnw + (capnc-capnw) * min(1._r8,max(0._r8,1.0_r8-(p(i,k)-0.8_r8*p(i,pver))/(0.2_r8*p(i,pver)))) |
---|
1241 | ! Modify for snow depth over land |
---|
1242 | capn = capn + (capnc-capn) * min(1.0_r8,max(0.0_r8,snowh(i)*10._r8)) |
---|
1243 | ! Ramp between polluted value over land to clean value over ocean. |
---|
1244 | capn = capn + (capnc-capn) * min(1.0_r8,max(0.0_r8,1.0_r8-landm(i,lchnk))) |
---|
1245 | ! Ramp between the resultant value and a sea ice value in the presence of ice. |
---|
1246 | capn = capn + (capnsi-capn) * min(1.0_r8,max(0.0_r8,seaicef(i))) |
---|
1247 | ! end jrm |
---|
1248 | ! |
---|
1249 | #ifdef DEBUG2 |
---|
1250 | if ( (lat(i) == latlook(1)) .or. (lat(i) == latlook(2)) ) then |
---|
1251 | if (i == ilook(1)) then |
---|
1252 | write(iulog,*) ' findmcnew: lat, k, seaicef, landm, wp, capnoice, capn ', & |
---|
1253 | lat(i), k, seaicef(i), landm(i,lat(i)), wp, capnoice, capn |
---|
1254 | endif |
---|
1255 | endif |
---|
1256 | #endif |
---|
1257 | |
---|
1258 | ! |
---|
1259 | ! useful terms in following calculations |
---|
1260 | ! |
---|
1261 | rat1 = rhocgs/rhow |
---|
1262 | rat2 = liqmr(i)/capn |
---|
1263 | con2 = (rat1*rat2)**0.333_r8 |
---|
1264 | ! |
---|
1265 | ! volume radius |
---|
1266 | ! |
---|
1267 | ! r3l = (rhocgs*liqmr(i)/(1.333*pi*capn*rhow))**0.333 * 0.01 ! meters |
---|
1268 | r3l = con1*con2 |
---|
1269 | ! |
---|
1270 | ! critical threshold for autoconversion if modified for mixed phase |
---|
1271 | ! clouds to mimic a bergeron findeisen process |
---|
1272 | ! r3lc2 = r3lcrit*(1.-0.5*fice(i)*(1-fice(i))) |
---|
1273 | ! |
---|
1274 | ! autoconversion of liquid |
---|
1275 | ! |
---|
1276 | ! cwaut = 2.e-4 |
---|
1277 | ! cwaut = 1.e-3 |
---|
1278 | ! lcrit = 2.e-4 |
---|
1279 | ! lcrit = 5.e-4 |
---|
1280 | ! pwaut = max(0._r8,liqmr(i)-lcrit)*cwaut |
---|
1281 | ! |
---|
1282 | ! pwaut is following tripoli and cotton (and many others) |
---|
1283 | ! we reduce the autoconversion below critpr, because these are regions where |
---|
1284 | ! the drop size distribution is likely to imply much smaller collector drops than |
---|
1285 | ! those relevant for a cloud distribution corresponding to the value of effc = 0.55 |
---|
1286 | ! suggested by cotton (see austin 1995 JAS, baker 1993) |
---|
1287 | |
---|
1288 | ! easy to follow form |
---|
1289 | ! pwaut = capc*liqmr(i)**2*rhocgs/rhow |
---|
1290 | ! $ *(liqmr(i)*rhocgs/(rhow*capn))**(.333) |
---|
1291 | ! $ *heavy(r3l,r3lcrit) |
---|
1292 | ! $ *max(0.10_r8,min(1._r8,prlloc(i)/critpr)) |
---|
1293 | ! somewhat faster form |
---|
1294 | #define HEAVYNEW |
---|
1295 | #ifdef HEAVYNEW |
---|
1296 | !#ifdef PERGRO |
---|
1297 | pwaut = capc*liqmr(i)**2*rat1*con2*heavymp(r3l,r3lcrit) * & |
---|
1298 | max(0.10_r8,min(1._r8,prlloc(i)/critpr)) |
---|
1299 | #else |
---|
1300 | pwaut = capc*liqmr(i)**2*rat1*con2*heavym(r3l,r3lcrit)* & |
---|
1301 | max(0.10_r8,min(1._r8,prlloc(i)/critpr)) |
---|
1302 | #endif |
---|
1303 | ! |
---|
1304 | ! autoconversion of ice |
---|
1305 | ! |
---|
1306 | ! ciaut = ciautb*efact |
---|
1307 | ciaut = ciautb |
---|
1308 | ! psaut = capc*totmr(i)**2*rhocgs/rhoi |
---|
1309 | ! $ *(totmr(i)*rhocgs/(rhoi*capn))**(.333) |
---|
1310 | ! |
---|
1311 | ! autoconversion of ice condensate |
---|
1312 | ! |
---|
1313 | #ifdef PERGRO |
---|
1314 | psaut = heavyp(icemr(i),icrit)*icemr(i)*ciaut |
---|
1315 | #else |
---|
1316 | psaut = max(0._r8,icemr(i)-icrit)*ciaut |
---|
1317 | #endif |
---|
1318 | ! |
---|
1319 | ! collection of liquid by rain |
---|
1320 | ! |
---|
1321 | ! pracw = cracw*rho(i)*liqmr(i)*rainmr(i) !(beheng 1994) |
---|
1322 | pracw = cracw*rho(i)*sqrt(rho(i))*liqmr(i)*rainmr(i) !(tripoli and cotton) |
---|
1323 | !! pracw = 0. |
---|
1324 | ! |
---|
1325 | ! the following lines calculate the slope parameter and snow mixing ratio |
---|
1326 | ! from the precip rate using the equations found in lin et al 83 |
---|
1327 | ! in the most natural form, but it is expensive, so after some tedious |
---|
1328 | ! algebraic manipulation you can use the cheaper form found below |
---|
1329 | ! vfalls = c*gam4pd/(6*lamdas**d)*sqrt(rhonot/rhocgs) |
---|
1330 | ! $ *0.01 ! convert from cm/s to m/s |
---|
1331 | ! snowmr(i) = snowfr*precab(i)/(rho(i)*vfalls*cldpr(i)) |
---|
1332 | ! snowmr(i) = ( prscgs(i)*mcon02 * (rhocgs**mcon03) )**mcon04 |
---|
1333 | ! lamdas = (prhonos/max(rhocgs*snowmr(i),small))**0.25 |
---|
1334 | ! csacw = mcon01*sqrt(rhonot/rhocgs)/(lamdas**thrpd) |
---|
1335 | ! |
---|
1336 | ! coefficient for collection by snow independent of phase |
---|
1337 | ! |
---|
1338 | csacx = mcon07*rhocgs**mcon08*prscgs(i)**mcon05 |
---|
1339 | |
---|
1340 | ! |
---|
1341 | ! collection of liquid by snow (lin et al 1983) |
---|
1342 | ! |
---|
1343 | psacw = csacx*liqmr(i)*esw |
---|
1344 | #ifdef PERGRO |
---|
1345 | ! this is necessary for pergro |
---|
1346 | psacw = 0._r8 |
---|
1347 | #endif |
---|
1348 | ! |
---|
1349 | ! collection of ice by snow (lin et al 1983) |
---|
1350 | ! |
---|
1351 | psaci = csacx*icemr(i)*esi |
---|
1352 | ! |
---|
1353 | ! total conversion of condensate to precipitate |
---|
1354 | ! |
---|
1355 | ptot = pwaut + psaut + pracw + psacw + psaci |
---|
1356 | ! |
---|
1357 | ! the recipricol of cloud water amnt (or zero if no cloud water) |
---|
1358 | ! |
---|
1359 | ! rcwm = totmr(i)/(max(totmr(i),small)**2) |
---|
1360 | ! |
---|
1361 | ! turn the tendency back into a loss rate (1/seconds) |
---|
1362 | ! |
---|
1363 | if (totmr(i) > 0._r8) then |
---|
1364 | coef(i) = ptot/totmr(i) |
---|
1365 | else |
---|
1366 | coef(i) = 0._r8 |
---|
1367 | endif |
---|
1368 | |
---|
1369 | if (ptot.gt.0._r8) then |
---|
1370 | fwaut(i) = pwaut/ptot |
---|
1371 | fsaut(i) = psaut/ptot |
---|
1372 | fracw(i) = pracw/ptot |
---|
1373 | fsacw(i) = psacw/ptot |
---|
1374 | fsaci(i) = psaci/ptot |
---|
1375 | else |
---|
1376 | fwaut(i) = 0._r8 |
---|
1377 | fsaut(i) = 0._r8 |
---|
1378 | fracw(i) = 0._r8 |
---|
1379 | fsacw(i) = 0._r8 |
---|
1380 | fsaci(i) = 0._r8 |
---|
1381 | endif |
---|
1382 | |
---|
1383 | ftot = fwaut(i)+fsaut(i)+fracw(i)+fsacw(i)+fsaci(i) |
---|
1384 | ! if (abs(ftot-1._r8).gt.1.e-14_r8.and.ftot.ne.0._r8) then |
---|
1385 | ! write(iulog,*) ' something is wrong in findmcnew ', ftot, & |
---|
1386 | ! fwaut(i),fsaut(i),fracw(i),fsacw(i),fsaci(i) |
---|
1387 | ! write(iulog,*) ' unscaled ', ptot, & |
---|
1388 | ! pwaut,psaut,pracw,psacw,psaci |
---|
1389 | ! write(iulog,*) ' totmr, liqmr, icemr ', totmr(i), liqmr(i), icemr(i) |
---|
1390 | ! call endrun() |
---|
1391 | ! endif |
---|
1392 | end do |
---|
1393 | #ifdef DEBUG |
---|
1394 | i = icollook(nlook) |
---|
1395 | if (lchnk == lchnklook(nlook) ) then |
---|
1396 | write(iulog,*) |
---|
1397 | write(iulog,*) '------', k, i, lchnk |
---|
1398 | write(iulog,*) ' liqmr, rainmr,precab ', liqmr(i), rainmr(i), precab(i)*8.64e4_r8 |
---|
1399 | write(iulog,*) ' frac: waut,saut,racw,sacw,saci ', & |
---|
1400 | fwaut(i), fsaut(i), fracw(i), fsacw(i), fsaci(i) |
---|
1401 | endif |
---|
1402 | #endif |
---|
1403 | |
---|
1404 | return |
---|
1405 | end subroutine findmcnew |
---|
1406 | |
---|
1407 | !############################################################################## |
---|
1408 | |
---|
1409 | subroutine findsp (lchnk, ncol, q, t, p, tsp, qsp) |
---|
1410 | !----------------------------------------------------------------------- |
---|
1411 | ! |
---|
1412 | ! Purpose: |
---|
1413 | ! find the wet bulb temperature for a given t and q |
---|
1414 | ! in a longitude height section |
---|
1415 | ! wet bulb temp is the temperature and spec humidity that is |
---|
1416 | ! just saturated and has the same enthalpy |
---|
1417 | ! if q > qs(t) then tsp > t and qsp = qs(tsp) < q |
---|
1418 | ! if q < qs(t) then tsp < t and qsp = qs(tsp) > q |
---|
1419 | ! |
---|
1420 | ! Method: |
---|
1421 | ! a Newton method is used |
---|
1422 | ! first guess uses an algorithm provided by John Petch from the UKMO |
---|
1423 | ! we exclude points where the physical situation is unrealistic |
---|
1424 | ! e.g. where the temperature is outside the range of validity for the |
---|
1425 | ! saturation vapor pressure, or where the water vapor pressure |
---|
1426 | ! exceeds the ambient pressure, or the saturation specific humidity is |
---|
1427 | ! unrealistic |
---|
1428 | ! |
---|
1429 | ! Author: P. Rasch |
---|
1430 | ! |
---|
1431 | !----------------------------------------------------------------------- |
---|
1432 | ! |
---|
1433 | ! input arguments |
---|
1434 | ! |
---|
1435 | integer, intent(in) :: lchnk ! chunk identifier |
---|
1436 | integer, intent(in) :: ncol ! number of atmospheric columns |
---|
1437 | |
---|
1438 | real(r8), intent(in) :: q(pcols,pver) ! water vapor (kg/kg) |
---|
1439 | real(r8), intent(in) :: t(pcols,pver) ! temperature (K) |
---|
1440 | real(r8), intent(in) :: p(pcols,pver) ! pressure (Pa) |
---|
1441 | ! |
---|
1442 | ! output arguments |
---|
1443 | ! |
---|
1444 | real(r8), intent(out) :: tsp(pcols,pver) ! saturation temp (K) |
---|
1445 | real(r8), intent(out) :: qsp(pcols,pver) ! saturation mixing ratio (kg/kg) |
---|
1446 | ! |
---|
1447 | ! local variables |
---|
1448 | ! |
---|
1449 | integer i ! work variable |
---|
1450 | integer k ! work variable |
---|
1451 | logical lflg ! work variable |
---|
1452 | integer iter ! work variable |
---|
1453 | integer l ! work variable |
---|
1454 | logical :: error_found |
---|
1455 | |
---|
1456 | real(r8) omeps ! 1 minus epsilon |
---|
1457 | real(r8) trinv ! work variable |
---|
1458 | real(r8) es ! sat. vapor pressure |
---|
1459 | real(r8) desdt ! change in sat vap pressure wrt temperature |
---|
1460 | ! real(r8) desdp ! change in sat vap pressure wrt pressure |
---|
1461 | real(r8) dqsdt ! change in sat spec. hum. wrt temperature |
---|
1462 | real(r8) dgdt ! work variable |
---|
1463 | real(r8) g ! work variable |
---|
1464 | real(r8) weight(pcols) ! work variable |
---|
1465 | real(r8) hlatsb ! (sublimation) |
---|
1466 | real(r8) hlatvp ! (vaporization) |
---|
1467 | real(r8) hltalt(pcols,pver) ! lat. heat. of vap. |
---|
1468 | real(r8) tterm ! work var. |
---|
1469 | real(r8) qs ! spec. hum. of water vapor |
---|
1470 | real(r8) tc ! crit temp of transition to ice |
---|
1471 | |
---|
1472 | ! work variables |
---|
1473 | real(r8) t1, q1, dt, dq |
---|
1474 | real(r8) dtm, dqm |
---|
1475 | real(r8) qvd, a1, tmp |
---|
1476 | real(r8) rair |
---|
1477 | real(r8) r1b, c1, c2, c3 |
---|
1478 | real(r8) denom |
---|
1479 | real(r8) dttol |
---|
1480 | real(r8) dqtol |
---|
1481 | integer doit(pcols) |
---|
1482 | real(r8) enin(pcols), enout(pcols) |
---|
1483 | real(r8) tlim(pcols) |
---|
1484 | |
---|
1485 | omeps = 1.0_r8 - epsqs |
---|
1486 | trinv = 1.0_r8/ttrice |
---|
1487 | a1 = 7.5_r8*log(10._r8) |
---|
1488 | rair = 287.04_r8 |
---|
1489 | c3 = rair*a1/cp |
---|
1490 | dtm = 0._r8 ! needed for iter=0 blowup with f90 -ei |
---|
1491 | dqm = 0._r8 ! needed for iter=0 blowup with f90 -ei |
---|
1492 | dttol = 1.e-4_r8 ! the relative temp error tolerance required to quit the iteration |
---|
1493 | dqtol = 1.e-4_r8 ! the relative moisture error tolerance required to quit the iteration |
---|
1494 | ! tmin = 173.16 ! the coldest temperature we can deal with |
---|
1495 | ! |
---|
1496 | ! max number of times to iterate the calculation |
---|
1497 | iter = 8 |
---|
1498 | ! |
---|
1499 | do k = k1mb,pver |
---|
1500 | |
---|
1501 | ! |
---|
1502 | ! first guess on the wet bulb temperature |
---|
1503 | ! |
---|
1504 | do i = 1,ncol |
---|
1505 | |
---|
1506 | #ifdef DEBUG |
---|
1507 | if ( (lchnk == lchnklook(nlook) ) .and. (i == icollook(nlook) ) ) then |
---|
1508 | write(iulog,*) ' ' |
---|
1509 | write(iulog,*) ' level, t, q, p', k, t(i,k), q(i,k), p(i,k) |
---|
1510 | endif |
---|
1511 | #endif |
---|
1512 | ! limit the temperature range to that relevant to the sat vap pres tables |
---|
1513 | #if ( ! defined WACCM_PHYS ) |
---|
1514 | tlim(i) = min(max(t(i,k),173._r8),373._r8) |
---|
1515 | #else |
---|
1516 | tlim(i) = min(max(t(i,k),128._r8),373._r8) |
---|
1517 | #endif |
---|
1518 | es = estblf(tlim(i)) |
---|
1519 | denom = p(i,k) - omeps*es |
---|
1520 | qs = epsqs*es/denom |
---|
1521 | doit(i) = 0 |
---|
1522 | enout(i) = 1._r8 |
---|
1523 | ! make sure a meaningful calculation is possible |
---|
1524 | if (p(i,k) > 5._r8*es .and. qs > 0._r8 .and. qs < 0.5_r8) then |
---|
1525 | ! |
---|
1526 | ! Saturation specific humidity |
---|
1527 | ! |
---|
1528 | qs = min(epsqs*es/denom,1._r8) |
---|
1529 | ! |
---|
1530 | ! "generalized" analytic expression for t derivative of es |
---|
1531 | ! accurate to within 1 percent for 173.16 < t < 373.16 |
---|
1532 | ! |
---|
1533 | ! Weighting of hlat accounts for transition from water to ice |
---|
1534 | ! polynomial expression approximates difference between es over |
---|
1535 | ! water and es over ice from 0 to -ttrice (C) (min of ttrice is |
---|
1536 | ! -40): required for accurate estimate of es derivative in transition |
---|
1537 | ! range from ice to water also accounting for change of hlatv with t |
---|
1538 | ! above freezing where const slope is given by -2369 j/(kg c) = cpv - cw |
---|
1539 | ! |
---|
1540 | tc = tlim(i) - t0 |
---|
1541 | lflg = (tc >= -ttrice .and. tc < 0.0_r8) |
---|
1542 | weight(i) = min(-tc*trinv,1.0_r8) |
---|
1543 | hlatsb = hlatv + weight(i)*hlatf |
---|
1544 | hlatvp = hlatv - 2369.0_r8*tc |
---|
1545 | if (tlim(i) < t0) then |
---|
1546 | hltalt(i,k) = hlatsb |
---|
1547 | else |
---|
1548 | hltalt(i,k) = hlatvp |
---|
1549 | end if |
---|
1550 | enin(i) = cp*tlim(i) + hltalt(i,k)*q(i,k) |
---|
1551 | |
---|
1552 | ! make a guess at the wet bulb temp using a UKMO algorithm (from J. Petch) |
---|
1553 | tmp = q(i,k) - qs |
---|
1554 | c1 = hltalt(i,k)*c3 |
---|
1555 | c2 = (tlim(i) + 36._r8)**2 |
---|
1556 | r1b = c2/(c2 + c1*qs) |
---|
1557 | qvd = r1b*tmp |
---|
1558 | tsp(i,k) = tlim(i) + ((hltalt(i,k)/cp)*qvd) |
---|
1559 | #ifdef DEBUG |
---|
1560 | if ( (lchnk == lchnklook(nlook) ) .and. (i == icollook(nlook) ) ) then |
---|
1561 | write(iulog,*) ' relative humidity ', q(i,k)/qs |
---|
1562 | write(iulog,*) ' first guess ', tsp(i,k) |
---|
1563 | endif |
---|
1564 | #endif |
---|
1565 | es = estblf(tsp(i,k)) |
---|
1566 | qsp(i,k) = min(epsqs*es/(p(i,k) - omeps*es),1._r8) |
---|
1567 | else |
---|
1568 | doit(i) = 1 |
---|
1569 | tsp(i,k) = tlim(i) |
---|
1570 | qsp(i,k) = q(i,k) |
---|
1571 | enin(i) = 1._r8 |
---|
1572 | endif |
---|
1573 | end do ! end do i |
---|
1574 | ! |
---|
1575 | ! now iterate on first guess |
---|
1576 | ! |
---|
1577 | do l = 1, iter |
---|
1578 | dtm = 0 |
---|
1579 | dqm = 0 |
---|
1580 | do i = 1,ncol |
---|
1581 | if (doit(i) == 0) then |
---|
1582 | es = estblf(tsp(i,k)) |
---|
1583 | ! |
---|
1584 | ! Saturation specific humidity |
---|
1585 | ! |
---|
1586 | qs = min(epsqs*es/(p(i,k) - omeps*es),1._r8) |
---|
1587 | ! |
---|
1588 | ! "generalized" analytic expression for t derivative of es |
---|
1589 | ! accurate to within 1 percent for 173.16 < t < 373.16 |
---|
1590 | ! |
---|
1591 | ! Weighting of hlat accounts for transition from water to ice |
---|
1592 | ! polynomial expression approximates difference between es over |
---|
1593 | ! water and es over ice from 0 to -ttrice (C) (min of ttrice is |
---|
1594 | ! -40): required for accurate estimate of es derivative in transition |
---|
1595 | ! range from ice to water also accounting for change of hlatv with t |
---|
1596 | ! above freezing where const slope is given by -2369 j/(kg c) = cpv - cw |
---|
1597 | ! |
---|
1598 | tc = tsp(i,k) - t0 |
---|
1599 | lflg = (tc >= -ttrice .and. tc < 0.0_r8) |
---|
1600 | weight(i) = min(-tc*trinv,1.0_r8) |
---|
1601 | hlatsb = hlatv + weight(i)*hlatf |
---|
1602 | hlatvp = hlatv - 2369.0_r8*tc |
---|
1603 | if (tsp(i,k) < t0) then |
---|
1604 | hltalt(i,k) = hlatsb |
---|
1605 | else |
---|
1606 | hltalt(i,k) = hlatvp |
---|
1607 | end if |
---|
1608 | if (lflg) then |
---|
1609 | tterm = pcf(1) + tc*(pcf(2) + tc*(pcf(3)+tc*(pcf(4) + tc*pcf(5)))) |
---|
1610 | else |
---|
1611 | tterm = 0.0_r8 |
---|
1612 | end if |
---|
1613 | desdt = hltalt(i,k)*es/(rgasv*tsp(i,k)*tsp(i,k)) + tterm*trinv |
---|
1614 | dqsdt = (epsqs + omeps*qs)/(p(i,k) - omeps*es)*desdt |
---|
1615 | ! g = cp*(tlim(i)-tsp(i,k)) + hltalt(i,k)*q(i,k)- hltalt(i,k)*qsp(i,k) |
---|
1616 | g = enin(i) - (cp*tsp(i,k) + hltalt(i,k)*qsp(i,k)) |
---|
1617 | dgdt = -(cp + hltalt(i,k)*dqsdt) |
---|
1618 | t1 = tsp(i,k) - g/dgdt |
---|
1619 | dt = abs(t1 - tsp(i,k))/t1 |
---|
1620 | tsp(i,k) = max(t1,tmin) |
---|
1621 | es = estblf(tsp(i,k)) |
---|
1622 | q1 = min(epsqs*es/(p(i,k) - omeps*es),1._r8) |
---|
1623 | dq = abs(q1 - qsp(i,k))/max(q1,1.e-12_r8) |
---|
1624 | qsp(i,k) = q1 |
---|
1625 | #ifdef DEBUG |
---|
1626 | if ( (lchnk == lchnklook(nlook) ) .and. (i == icollook(nlook) ) ) then |
---|
1627 | write(iulog,*) ' rel chg lev, iter, t, q ', k, l, dt, dq, g |
---|
1628 | endif |
---|
1629 | #endif |
---|
1630 | dtm = max(dtm,dt) |
---|
1631 | dqm = max(dqm,dq) |
---|
1632 | ! if converged at this point, exclude it from more iterations |
---|
1633 | if (dt < dttol .and. dq < dqtol) then |
---|
1634 | doit(i) = 2 |
---|
1635 | endif |
---|
1636 | enout(i) = cp*tsp(i,k) + hltalt(i,k)*qsp(i,k) |
---|
1637 | ! bail out if we are too near the end of temp range |
---|
1638 | #if ( ! defined WACCM_PHYS ) |
---|
1639 | if (tsp(i,k) < 174.16_r8) then |
---|
1640 | #else |
---|
1641 | if (tsp(i,k) < 130.16_r8) then |
---|
1642 | #endif |
---|
1643 | doit(i) = 4 |
---|
1644 | endif |
---|
1645 | else |
---|
1646 | endif |
---|
1647 | end do ! do i = 1,ncol |
---|
1648 | |
---|
1649 | if (dtm < dttol .and. dqm < dqtol) then |
---|
1650 | go to 10 |
---|
1651 | endif |
---|
1652 | |
---|
1653 | end do ! do l = 1,iter |
---|
1654 | 10 continue |
---|
1655 | |
---|
1656 | error_found = .false. |
---|
1657 | if (dtm > dttol .or. dqm > dqtol) then |
---|
1658 | do i = 1,ncol |
---|
1659 | if (doit(i) == 0) error_found = .true. |
---|
1660 | end do |
---|
1661 | if (error_found) then |
---|
1662 | do i = 1,ncol |
---|
1663 | if (doit(i) == 0) then |
---|
1664 | write(iulog,*) ' findsp not converging at point i, k ', i, k |
---|
1665 | write(iulog,*) ' t, q, p, enin ', t(i,k), q(i,k), p(i,k), enin(i) |
---|
1666 | write(iulog,*) ' tsp, qsp, enout ', tsp(i,k), qsp(i,k), enout(i) |
---|
1667 | call endrun ('FINDSP') |
---|
1668 | endif |
---|
1669 | end do |
---|
1670 | endif |
---|
1671 | endif |
---|
1672 | do i = 1,ncol |
---|
1673 | if (doit(i) == 2 .and. abs((enin(i)-enout(i))/(enin(i)+enout(i))) > 1.e-4_r8) then |
---|
1674 | error_found = .true. |
---|
1675 | endif |
---|
1676 | end do |
---|
1677 | if (error_found) then |
---|
1678 | do i = 1,ncol |
---|
1679 | if (doit(i) == 2 .and. abs((enin(i)-enout(i))/(enin(i)+enout(i))) > 1.e-4_r8) then |
---|
1680 | write(iulog,*) ' the enthalpy is not conserved for point ', & |
---|
1681 | i, k, enin(i), enout(i) |
---|
1682 | write(iulog,*) ' t, q, p, enin ', t(i,k), q(i,k), p(i,k), enin(i) |
---|
1683 | write(iulog,*) ' tsp, qsp, enout ', tsp(i,k), qsp(i,k), enout(i) |
---|
1684 | call endrun ('FINDSP') |
---|
1685 | endif |
---|
1686 | end do |
---|
1687 | endif |
---|
1688 | |
---|
1689 | end do ! level loop (k=1,pver) |
---|
1690 | |
---|
1691 | return |
---|
1692 | end subroutine findsp |
---|
1693 | |
---|
1694 | #endif |
---|
1695 | |
---|
1696 | end module cldwat |
---|