Radiation Scheme in CCPP
module_bfmicrophysics.f
Go to the documentation of this file.
2 !
3  USE machine , ONLY : kind_phys
4  USE funcphys
5  USE physcons, cp => con_cp, rd => con_rd, rv => con_rv &
6  &, t0c => con_t0c, hvap => con_hvap, hfus => con_hfus &
7  &, eps => con_eps, epsm1 => con_epsm1 &
8  &, eps1 => con_fvirt, pi => con_pi, grav => con_g
9  implicit none
10 !
11 !--- Common block of constants used in column microphysics
12 !
13  real,private :: abfr, cbfr, ciacw, ciacr, c_n0r0, &
14  &CN0r0, CN0r_DMRmin, CN0r_DMRmax, CRACW, CRAUT, ESW0, &
15  &QAUTx, RFmax, RQR_DR1, RQR_DR2, RQR_DR3, RQR_DRmin, &
16  &RQR_DRmax, RR_DRmin, RR_DR1, RR_DR2, RR_DR3, RR_DRmax
17 !
18  integer, private :: mic_step
19 !
20 !--- Common block for lookup table used in calculating growth rates of
21 ! nucleated ice crystals growing in water saturated conditions
22 !--- Discretized growth rates of small ice crystals after their nucleation
23 ! at 1 C intervals from -1 C to -35 C, based on calculations by Miller
24 ! and Young (1979, JAS) after 600 s of growth. Resultant growth rates
25 ! are multiplied by physics time step in GSMCONST.
26 !
27  INTEGER, PRIVATE,PARAMETER :: my_t1=1, my_t2=35
28  REAL,PRIVATE,DIMENSION(MY_T1:MY_T2) :: my_growth
29 !
30 !--- Parameters for ice lookup tables, which establish the range of mean ice
31 ! particle diameters; from a minimum mean diameter of 0.05 mm (DMImin) to a
32 ! maximum mean diameter of 1.00 mm (DMImax). The tables store solutions
33 ! at 1 micron intervals (DelDMI) of mean ice particle diameter.
34 !
35  REAL, PRIVATE,PARAMETER :: dmimin=.05e-3, dmimax=1.e-3, &
36  & xmimin=1.e6*dmimin, xmimax=1.e6*dmimax,&
37  & deldmi=1.e-6
38  INTEGER, PRIVATE,PARAMETER :: mdimin=xmimin, mdimax=xmimax
39 !
40 !--- Various ice lookup tables
41 !
42  REAL, PRIVATE,DIMENSION(MDImin:MDImax) :: &
43  & ACCRI,MASSI,SDENS,VSNOWI,VENTI1,VENTI2
44 !
45 !--- Mean rain drop diameters varying from 50 microns (0.05 mm) to 450 microns
46 ! (0.45 mm), assuming an exponential size distribution.
47 !
48  REAL, PRIVATE,PARAMETER :: dmrmin=.05e-3, dmrmax=.45e-3, &
49  & xmrmin=1.e6*dmrmin, xmrmax=1.e6*dmrmax,&
50  & deldmr=1.e-6, nlimin=100.
51 ! &, NLImin=100., NLImax=20.E3
52  INTEGER, PRIVATE,PARAMETER :: mdrmin=xmrmin, mdrmax=xmrmax
53 !
54 !--- Factor of 1.5 for RECImin, RESNOWmin, & RERAINmin accounts for
55 ! integrating exponential distributions for effective radius
56 ! (i.e., the r**3/r**2 moments).
57 !
58 ! INTEGER, PRIVATE, PARAMETER :: INDEXSmin=300
59 !! INTEGER, PRIVATE, PARAMETER :: INDEXSmin=200
60  INTEGER, PRIVATE, PARAMETER :: indexsmin=100
61  REAL, PRIVATE, PARAMETER :: rerainmin=1.5*xmrmin &
62 ! &, RECImin=1.5*XMImin, RESNOWmin=1.5*INDEXSmin, RECWmin=8.0
63 ! &, RECImin=1.5*XMImin, RESNOWmin=1.5*INDEXSmin, RECWmin=7.5
64  &, recimin=1.5*xmimin, resnowmin=1.5*indexsmin, recwmin=10.
65 ! &, RECImin=1.5*XMImin, RESNOWmin=1.5*INDEXSmin, RECWmin=15.
66 ! &, RECImin=1.5*XMImin, RESNOWmin=1.5*INDEXSmin, RECWmin=5.
67 
68 !
69 !--- Various rain lookup tables
70 !--- Rain lookup tables for mean rain drop diameters from DMRmin to DMRmax,
71 ! assuming exponential size distributions for the rain drops
72 !
73  REAL, PRIVATE,DIMENSION(MDRmin:MDRmax):: &
74  & ACCRR,MASSR,RRATE,VRAIN,VENTR1,VENTR2
75 !
76 !--- Common block for riming tables
77 !--- VEL_RF - velocity increase of rimed particles as functions of crude
78 ! particle size categories (at 0.1 mm intervals of mean ice particle
79 ! sizes) and rime factor (different values of Rime Factor of 1.1**N,
80 ! where N=0 to Nrime).
81 !
82  INTEGER, PRIVATE,PARAMETER :: nrime=40
83  REAL, DIMENSION(2:9,0:Nrime),PRIVATE :: vel_rf
84 !
85 !--- The following variables are for microphysical statistics
86 !
87  INTEGER, PARAMETER :: itlo=-60, ithi=40
88  INTEGER nstats(itlo:ithi,4)
89  REAL qmax(itlo:ithi,5), qtot(itlo:ithi,22)
90 !
91  REAL, PRIVATE, PARAMETER :: &
92 ! & T_ICE=-10., T_ICE_init=-5. !- Ver1
93 !!! &, T_ICE=-20. !- Ver2
94  & t_ice=-40., t_ice_init=-15. !- Ver2
95 ! & T_ICE=-30., T_ICE_init=-5. !- Ver2
96 !
97 ! Some other miscellaneous parameters
98 !
99  REAL, PRIVATE, PARAMETER :: thom=t_ice, tnw=50., toler=1.0e-20 &
100 ! REAL, PRIVATE, PARAMETER :: Thom=T_ICE, TNW=50., TOLER=5.E-7
101 ! REAL, PRIVATE, PARAMETER :: Thom=-35., TNW=50., TOLER=5.E-7
102 
103 ! &, emisCU=.75/1.66 ! Used for convective cloud l/w emissivity
104 
105 ! Assume fixed cloud ice effective radius
106  &, recice=recimin &
107  &, epsq=1.0e-20 &
108 ! &, EPSQ=1.E-12 &
109  &, flg0p1=0.1, flg0p2=0.2, flg1p0=1.0
110 !
111 !
112  CONTAINS
113 !
114 !#######################################################################
115 !------- Initialize constants & lookup tables for microphysics ---------
116 !#######################################################################
117 !
118  SUBROUTINE gsmconst (DTPG,mype,first)
119 !
120  implicit none
121 !-------------------------------------------------------------------------------
122 !--- SUBPROGRAM DOCUMENTATION BLOCK
123 ! PRGRMMR: Ferrier ORG: W/NP22 DATE: February 2001
124 !-------------------------------------------------------------------------------
125 ! ABSTRACT:
126 ! * Reads various microphysical lookup tables used in COLUMN_MICRO
127 ! * Lookup tables were created "offline" and are read in during execution
128 ! * Creates lookup tables for saturation vapor pressure w/r/t water & ice
129 !-------------------------------------------------------------------------------
130 !
131 ! USAGE: CALL GSMCONST FROM SUBROUTINE GSMDRIVE AT MODEL START TIME
132 !
133 ! INPUT ARGUMENT LIST:
134 ! DTPH - physics time step (s)
135 !
136 ! OUTPUT ARGUMENT LIST:
137 ! NONE
138 !
139 ! OUTPUT FILES:
140 ! NONE
141 !
142 !
143 ! SUBROUTINES:
144 ! MY_GROWTH_RATES - lookup table for growth of nucleated ice
145 !
146 ! UNIQUE: NONE
147 !
148 ! LIBRARY: NONE
149 !
150 !
151 ! ATTRIBUTES:
152 ! LANGUAGE: FORTRAN 90
153 ! MACHINE : IBM SP
154 !
155  integer mype
156  real dtpg
157  logical first
158 !
159 !--- Parameters & data statement for local calculations
160 !
161  REAL, PARAMETER :: C1=1./3., dmr1=.1e-3, dmr2=.2e-3, dmr3=.32e-3, &
162  & n0r0=8.e6, n0s0=4.e6, rhol=1000., rhos=100., &
163  & xmr1=1.e6*dmr1, xmr2=1.e6*dmr2, xmr3=1.e6*dmr3
164  INTEGER, PARAMETER :: MDR1=xmr1, mdr2=xmr2, mdr3=xmr3
165 !
166  real dtph, bbfr
167  integer i
168 !
169 !--- Added on 5/16/01 for Moorthi
170 !
171  logical, parameter :: read_lookup=.false., write_lookup=.false.
172 !
173 !------------------------------------------------------------------------
174 ! ************* Parameters used in ETA model -- Not used in Global Model *****
175 !
176 !--- DPHD, DLMD are delta latitude and longitude at the model (NOT geodetic) equator
177 ! => "DX" is the hypotenuse of the model zonal & meridional grid increments.
178 !
179 ! DX=111.*(DPHD**2+DLMD**2)**.5 ! Resolution at MODEL equator (km)
180 ! DX=MIN(100., MAX(5., DX) )
181 !
182 !--- Assume the following functional relationship for key constants that
183 ! depend on grid resolution from DXmin (5 km) to DXmax (100 km) resolution:
184 !
185 ! DXmin=5.
186 ! DXmax=100.
187 ! DX=MIN(DXmax, MAX(DXmin, DX) )
188 !
189 !--- EXtune determines the degree to which the coefficients change with resolution.
190 ! The larger EXtune is, the more sensitive the parameter.
191 !
192 ! EXtune=1.
193 
194 !
195 !--- FXtune ==> F(DX) is the grid-resolution tuning parameter (from 0 to 1)
196 !
197 ! FXtune=((DXmax-DX)/(DXmax-DXmin))**EXtune
198 ! FXtune=MAX(0., MIN(1., FXtune))
199 !
200 !--- Calculate grid-averaged RH for the onset of condensation (RHgrd) based on
201 ! simple ***ASSUMED*** (user-specified) values at DXmax and at DXmin.
202 !
203 ! RH_DXmax=.90 !-- 90% RH at DXmax=100 km
204 ! RH_DXmin=.98 !-- 98% RH at DXmin=5 km
205 !
206 !--- Note that RHgrd is right now fixed throughout the domain!!
207 !
208 ! RHgrd=RH_DXmax+(RH_DXmin-RH_DXmax)*FXtune
209 ! ********************************************************************************
210 !
211 !
212  if (first) then
213 !
214 !--- Read in various lookup tables
215 !
216  if ( read_lookup ) then
217  OPEN (unit=1,file="eta_micro_lookup.dat",form="UNFORMATTED")
218  READ(1) ventr1
219  READ(1) ventr2
220  READ(1) accrr
221  READ(1) massr
222  READ(1) vrain
223  READ(1) rrate
224  READ(1) venti1
225  READ(1) venti2
226  READ(1) accri
227  READ(1) massi
228  READ(1) vsnowi
229  READ(1) vel_rf
230 ! read(1) my_growth ! Applicable only for DTPH=180 s for offline testing
231  CLOSE (1)
232  else
233  CALL ice_lookup ! Lookup tables for ice
234  CALL rain_lookup ! Lookup tables for rain
235  if (write_lookup) then
236  open(unit=1,file='micro_lookup.dat',form='unformatted')
237  write(1) ventr1
238  write(1) ventr2
239  write(1) accrr
240  write(1) massr
241  write(1) vrain
242  write(1) rrate
243  write(1) venti1
244  write(1) venti2
245  write(1) accri
246  write(1) massi
247  write(1) vsnowi
248  write(1) vel_rf
249 ! write(1) my_growth ! Applicable only for DTPH=180 s ????
250  CLOSE (1)
251  endif
252  endif
253 !!
254 !--- Constants associated with Biggs (1953) freezing of rain, as parameterized
255 ! following Lin et al. (JCAM, 1983) & Reisner et al. (1998, QJRMS).
256 !
257  abfr=-0.66
258  bbfr=100.
259  cbfr=20.*pi*pi*bbfr*rhol*1.e-21
260 !
261 !--- QAUT0 is the threshold cloud content for autoconversion to rain
262 ! needed for droplets to reach a diameter of 20 microns (following
263 ! Manton and Cotton, 1977; Banta and Hanson, 1987, JCAM). It is
264 ! **STRONGLY** affected by the assumed droplet number concentrations
265 ! XNCW! For example, QAUT0=1.2567, 0.8378, or 0.4189 g/m**3 for
266 ! droplet number concentrations of 300, 200, and 100 cm**-3, respectively.
267 !
268 !--- Calculate grid-averaged XNCW based on simple ***ASSUMED*** (user-specified)
269 ! values at DXmax and at DXmin.
270 !
271 ! XNCW_DXmax=50.E6 !-- 50 /cm**3 at DXmax=100 km
272 ! XNCW_DXmin=200.E6 !-- 200 /cm**3 at DXmin=5 km
273 !
274 !--- Note that XNCW is right now fixed throughout the domain!!
275 !
276 ! XNCW=XNCW_DXmax+(XNCW_DXmin-XNCW_DXmax)*FXtune
277 !
278 ! QAUT0=PI*RHOL*XNCW*(20.E-6)**3/6.
279  qautx=pi*rhol*1.0e6*(20.e-6)**3/6.
280 !
281 !--- Based on rain lookup tables for mean diameters from 0.05 to 0.45 mm
282 ! * Four different functional relationships of mean drop diameter as
283 ! a function of rain rate (RR), derived based on simple fits to
284 ! mass-weighted fall speeds of rain as functions of mean diameter
285 ! from the lookup tables.
286 !
287  rr_drmin=n0r0*rrate(mdrmin) ! RR for mean drop diameter of .05 mm
288  rr_dr1=n0r0*rrate(mdr1) ! RR for mean drop diameter of .10 mm
289  rr_dr2=n0r0*rrate(mdr2) ! RR for mean drop diameter of .20 mm
290  rr_dr3=n0r0*rrate(mdr3) ! RR for mean drop diameter of .32 mm
291  rr_drmax=n0r0*rrate(mdrmax) ! RR for mean drop diameter of .45 mm
292 !
293  rqr_drmin=n0r0*massr(mdrmin) ! Rain content for mean drop diameter of .05 mm
294  rqr_dr1=n0r0*massr(mdr1) ! Rain content for mean drop diameter of .10 mm
295  rqr_dr2=n0r0*massr(mdr2) ! Rain content for mean drop diameter of .20 mm
296  rqr_dr3=n0r0*massr(mdr3) ! Rain content for mean drop diameter of .32 mm
297  rqr_drmax=n0r0*massr(mdrmax) ! Rain content for mean drop diameter of .45 mm
298  c_n0r0=pi*rhol*n0r0
299  cn0r0=1.e6/c_n0r0**.25
300  cn0r_dmrmin=1./(pi*rhol*dmrmin**4)
301  cn0r_dmrmax=1./(pi*rhol*dmrmax**4)
302 !
303  endif ! If (first) then loop ends here
304 !
305 ! Find out what microphysics time step should be
306 !
307  mic_step = max(1, int(dtpg/600.0+0.5))
308 ! mic_step = max(1, int(dtpg/300.0+0.5))
309  dtph = dtpg / mic_step
310  if (mype == 0) print *,' DTPG=',dtpg,' mic_step=',mic_step &
311  &, ' dtph=',dtph
312 !
313 !--- Calculates coefficients for growth rates of ice nucleated in water
314 ! saturated conditions, scaled by physics time step (lookup table)
315 !
316  CALL my_growth_rates (dtph)
317 !
318 !--- CIACW is used in calculating riming rates
319 ! The assumed effective collection efficiency of cloud water rimed onto
320 ! ice is =0.5 below:
321 !
322 !Moor CIACW=DTPH*0.25*PI*0.5*(1.E5)**C1 ! commented on 20050422
323 ! ice is =0.1 below:
324  ciacw=dtph*0.25*pi*0.1*(1.e5)**c1
325 ! CIACW = 0.0 ! Brad's suggestion 20040614
326 !
327 !--- CIACR is used in calculating freezing of rain colliding with large ice
328 ! The assumed collection efficiency is 1.0
329 !
330  ciacr=pi*dtph
331 !
332 !--- CRACW is used in calculating collection of cloud water by rain (an
333 ! assumed collection efficiency of 1.0)
334 !
335 !Moor CRACW=DTPH*0.25*PI*1.0 ! commented on 20050422
336 !
337 ! assumed collection efficiency of 0.1)
338  cracw=dtph*0.25*pi*0.1
339 ! CRACW = 0.0 ! Brad's suggestion 20040614
340 !
341  esw0=fpvsl(t0c) ! Saturation vapor pressure at 0C
342  rfmax=1.1**nrime ! Maximum rime factor allowed
343 !
344 !------------------------------------------------------------------------
345 !--------------- Constants passed through argument list -----------------
346 !------------------------------------------------------------------------
347 !
348 !--- Important parameters for self collection (autoconversion) of
349 ! cloud water to rain.
350 !
351 !--- CRAUT is proportional to the rate that cloud water is converted by
352 ! self collection to rain (autoconversion rate)
353 !
354  craut=1.-exp(-1.e-3*dtph)
355 !
356 ! IF (MYPE == 0)
357 ! & WRITE(6,"(A, A,F6.2,A, A,F5.4, A,F7.3,A, A,F6.2,A, A,F6.3,A)")
358 ! & 'KEY MICROPHYSICAL PARAMETERS FOR '
359 ! & ,'DX=',DX,' KM:'
360 ! & ,' FXtune=',FXtune
361 ! & ,' RHgrd=',100.*RHgrd,' %'
362 ! & ,' NCW=',1.E-6*XNCW,' /cm**3'
363 ! & ,' QAUT0=',1.E3*QAUT0,' g/kg'
364 !
365 !--- For calculating snow optical depths by considering bulk density of
366 ! snow based on emails from Q. Fu (6/27-28/01), where optical
367 ! depth (T) = 1.5*SWP/(Reff*DENS), SWP is snow water path, Reff
368 ! is effective radius, and DENS is the bulk density of snow.
369 !
370 ! SWP (kg/m**2)=(1.E-3 kg/g)*SWPrad, SWPrad in g/m**2 used in radiation
371 ! T = 1.5*1.E3*SWPrad/(Reff*DENS)
372 !
373 ! See derivation for MASSI(INDEXS), note equal to RHO*QSNOW/NSNOW
374 !
375 ! SDENS=1.5e3/DENS, DENS=MASSI(INDEXS)/[PI*(1.E-6*INDEXS)**3]
376 !
377  DO i=mdimin,mdimax
378 !MoorthiSDENS(I)=PI*1.5E-15*FLOAT(I*I*I)/MASSI(I)
379  sdens(i)=pi*1.0e-15*float(i*i*i)/massi(i)
380  ENDDO
381 !
382 !-----------------------------------------------------------------------
383 !
384  END subroutine gsmconst
385 
386 !
387 !#######################################################################
388 !--- Sets up lookup table for calculating initial ice crystal growth ---
389 !#######################################################################
390 !
391  SUBROUTINE my_growth_rates (DTPH)
392 !
393  implicit none
394 !
395 !--- Below are tabulated values for the predicted mass of ice crystals
396 ! after 600 s of growth in water saturated conditions, based on
397 ! calculations from Miller and Young (JAS, 1979). These values are
398 ! crudely estimated from tabulated curves at 600 s from Fig. 6.9 of
399 ! Young (1993). Values at temperatures colder than -27C were
400 ! assumed to be invariant with temperature.
401 !
402 !--- Used to normalize Miller & Young (1979) calculations of ice growth
403 ! over large time steps using their tabulated values at 600 s.
404 ! Assumes 3D growth with time**1.5 following eq. (6.3) in Young (1993).
405 !
406  real dtph, dt_ice
407  REAL MY_600(my_t1:my_t2)
408 !
409 !-- 20090714: These values are in g and need to be converted to kg below
410  DATA my_600 / &
411  & 5.5e-8, 1.4e-7, 2.8e-7, 6.e-7, 3.3e-6, & ! -1 to -5 deg C
412  & 2.e-6, 9.e-7, 8.8e-7, 8.2e-7, 9.4e-7, & ! -6 to -10 deg C
413  & 1.2e-6, 1.85e-6, 5.5e-6, 1.5e-5, 1.7e-5, & ! -11 to -15 deg C
414  & 1.5e-5, 1.e-5, 3.4e-6, 1.85e-6, 1.35e-6, & ! -16 to -20 deg C
415  & 1.05e-6, 1.e-6, 9.5e-7, 9.0e-7 , 9.5e-7, & ! -21 to -25 deg C
416  & 9.5e-7, 9.e-7, 9.e-7, 9.e-7, 9.e-7, & ! -26 to -30 deg C
417  & 9.e-7, 9.e-7, 9.e-7, 9.e-7, 9.e-7 / ! -31 to -35 deg C
418 !
419 !-----------------------------------------------------------------------
420 !
421  dt_ice=(dtph/600.)**1.5
422 ! MY_GROWTH=DT_ICE*MY_600 ! original version
423  my_growth=dt_ice*my_600*1.e-3 !-- 20090714: Convert from g to kg
424 !
425 !-----------------------------------------------------------------------
426 !
427  END subroutine my_growth_rates
428 !
429 !#######################################################################
430 !--------------- Creates lookup tables for ice processes ---------------
431 !#######################################################################
432 !
433  subroutine ice_lookup
434 !
435  implicit none
436 !-----------------------------------------------------------------------------------
437 !
438 !---- Key diameter values in mm
439 !
440 !-----------------------------------------------------------------------------------
441 !
442 !---- Key concepts:
443 ! - Actual physical diameter of particles (D)
444 ! - Ratio of actual particle diameters to mean diameter (x=D/MD)
445 ! - Mean diameter of exponentially distributed particles, which is the
446 ! same as 1./LAMDA of the distribution (MD)
447 ! - All quantitative relationships relating ice particle characteristics as
448 ! functions of their diameter (e.g., ventilation coefficients, normalized
449 ! accretion rates, ice content, and mass-weighted fall speeds) are a result
450 ! of using composite relationships for ice crystals smaller than 1.5 mm
451 ! diameter merged with analogous relationships for larger sized aggregates.
452 ! Relationships are derived as functions of mean ice particle sizes assuming
453 ! exponential size spectra and assuming the properties of ice crystals at
454 ! sizes smaller than 1.5 mm and aggregates at larger sizes.
455 !
456 !-----------------------------------------------------------------------------------
457 !
458 !---- Actual ice particle diameters for which integrated distributions are derived
459 ! - DminI - minimum diameter for integration (.02 mm, 20 microns)
460 ! - DmaxI - maximum diameter for integration (2 cm)
461 ! - DdelI - interval for integration (1 micron)
462 !
463  real, parameter :: DminI=.02e-3, dmaxi=20.e-3, ddeli=1.e-6, &
464  & ximin=1.e6*dmini, ximax=1.e6*dmaxi
465  integer, parameter :: IDImin=ximin, idimax=ximax
466 !
467 !---- Meaning of the following arrays:
468 ! - diam - ice particle diameter (m)
469 ! - mass - ice particle mass (kg)
470 ! - vel - ice particle fall speeds (m/s)
471 ! - vent1 - 1st term in ice particle ventilation factor
472 ! - vent2 - 2nd term in ice particle ventilation factor
473 !
474  real diam(idimin:idimax),mass(idimin:idimax),vel(idimin:idimax), &
475  & vent1(IDImin:IDImax),vent2(IDImin:IDImax)
476 !
477 !-----------------------------------------------------------------------------------
478 !
479 !---- Found from trial & error that the m(D) & V(D) mass & velocity relationships
480 ! between the ice crystals and aggregates overlapped & merged near a particle
481 ! diameter sizes of 1.5 mm. Thus, ice crystal relationships are used for
482 ! sizes smaller than 1.5 mm and aggregate relationships for larger sizes.
483 !
484  real, parameter :: d_crystal_max=1.5
485 !
486 !---- The quantity xmax represents the maximum value of "x" in which the
487 ! integrated values are calculated. For xmax=20., this means that
488 ! integrated ventilation, accretion, mass, and precipitation rates are
489 ! calculated for ice particle sizes less than 20.*mdiam, the mean particle diameter.
490 !
491  real, parameter :: xmax=20.
492 !
493 !-----------------------------------------------------------------------------------
494 !
495 !---- Meaning of the following arrays:
496 ! - mdiam - mean diameter (m)
497 ! - VENTI1 - integrated quantity associated w/ ventilation effects
498 ! (capacitance only) for calculating vapor deposition onto ice
499 ! - VENTI2 - integrated quantity associated w/ ventilation effects
500 ! (with fall speed) for calculating vapor deposition onto ice
501 ! - ACCRI - integrated quantity associated w/ cloud water collection by ice
502 ! - MASSI - integrated quantity associated w/ ice mass
503 ! - VSNOWI - mass-weighted fall speed of snow, used to calculate precip rates
504 !
505 !--- Mean ice-particle diameters varying from 50 microns to 1000 microns (1 mm),
506 ! assuming an exponential size distribution.
507 !
508  real mdiam
509 !
510 !-----------------------------------------------------------------------------------
511 !------------- Constants & parameters for ventilation factors of ice ---------------
512 !-----------------------------------------------------------------------------------
513 !
514 !---- These parameters are used for calculating the ventilation factors for ice
515 ! crystals between 0.2 and 1.5 mm diameter (Hall and Pruppacher, JAS, 1976).
516 ! From trial & error calculations, it was determined that the ventilation
517 ! factors of smaller ice crystals could be approximated by a simple linear
518 ! increase in the ventilation coefficient from 1.0 at 50 microns (.05 mm) to
519 ! 1.1 at 200 microns (0.2 mm), rather than using the more complex function of
520 ! 1.0 + .14*(Sc**.33*Re**.5)**2 recommended by Hall & Pruppacher.
521 !
522  real, parameter :: cvent1i=.86, cvent2i=.28
523 !
524 !---- These parameters are used for calculating the ventilation factors for larger
525 ! aggregates, where D>=1.5 mm (see Rutledge and Hobbs, JAS, 1983;
526 ! Thorpe and Mason, 1966).
527 !
528  real, parameter :: cvent1a=.65, cvent2a=.44
529 !
530  real m_agg,m_bullet,m_column,m_ice,m_plate
531 !
532 !---- Various constants
533 !
534  real, parameter :: c1=2./3., cexp=1./3.
535 !
536  logical :: iprint
537  logical, parameter :: print_diag=.false.
538 !
539 !-----------------------------------------------------------------------------------
540 !- Constants & parameters for calculating the increase in fall speed of rimed ice --
541 !-----------------------------------------------------------------------------------
542 !
543 !---- Constants & arrays for estimating increasing fall speeds of rimed ice.
544 ! Based largely on theory and results from Bohm (JAS, 1989, 2419-2427).
545 !
546 !-------------------- Standard atmosphere conditions at 1000 mb --------------------
547 !
548  real, parameter :: t_std=288., dens_std=1000.e2/(287.04*288.)
549 !
550 !---- These "bulk densities" are the actual ice densities in the ice portion of the
551 ! lattice. They are based on text associated w/ (12) on p. 2425 of Bohm (JAS,
552 ! 1989). Columns, plates, & bullets are assumed to have an average bulk density
553 ! of 850 kg/m**3. Aggregates were assumed to have a slightly larger bulk density
554 ! of 600 kg/m**3 compared with dendrites (i.e., the least dense, most "lacy" &
555 ! tenous ice crystal, which was assumed to be ~500 kg/m**3 in Bohm).
556 !
557  real, parameter :: dens_crystal=850., dens_agg=600.
558 !
559 !--- A value of Nrime=40 for a logarithmic ratio of 1.1 yields a maximum rime factor
560 ! of 1.1**40 = 45.26 that is resolved in these tables. This allows the largest
561 ! ice particles with a mean diameter of MDImax=1000 microns to achieve bulk
562 ! densities of 900 kg/m**3 for rimed ice.
563 !
564 ! integer, parameter :: Nrime=40
565  real m_rime, &
566  & rime_factor(0:Nrime), rime_vel(0:Nrime), &
567  & vel_rime(IDImin:IDImax,Nrime), ivel_rime(MDImin:MDImax,Nrime)
568 !
569  integer i, j, jj, k, icount
570  real c2, cbulk, cbulk_ice, px, dynvis_std, crime1 &
571  &, crime2, crime3, crime4, crime5, d, c_avg, c_agg &
572  &, c_bullet, c_column, c_plate, cl_agg, cl_bullet &
573  &, cl_column, cl_plate, v_agg, v_bullet, v_column &
574  &, v_plate, wd, ecc_column &
575  &, cvent1, cvent2, crime_best, rime_m1, rime_m2 &
576  &, x_rime, re_rime, smom3, pratei, expf &
577  &, bulk_dens, xmass, xmdiam, ecc_plate, dx
578 !
579 !-----------------------------------------------------------------------------------
580 !----------------------------- BEGIN EXECUTION -------------------------------------
581 !-----------------------------------------------------------------------------------
582 !
583 !
584  c2=1./sqrt(3.)
585 ! pi=acos(-1.)
586  cbulk=6./pi
587  cbulk_ice=900.*pi/6. ! Maximum bulk ice density allowed of 900 kg/m**3
588  px=.4**cexp ! Convert fall speeds from 400 mb (Starr & Cox) to 1000 mb
589 !
590 !--------------------- Dynamic viscosity (1000 mb, 288 K) --------------------------
591 !
592  dynvis_std=1.496e-6*t_std**1.5/(t_std+120.)
593  crime1=pi/24.
594  crime2=8.*9.81*dens_std/(pi*dynvis_std**2)
595  crime3=crime1*dens_crystal
596  crime4=crime1*dens_agg
597  crime5=dynvis_std/dens_std
598  do i=0,nrime
599  rime_factor(i)=1.1**i
600  enddo
601 !
602 !#######################################################################
603 ! Characteristics as functions of actual ice particle diameter
604 !#######################################################################
605 !
606 !---- M(D) & V(D) for 3 categories of ice crystals described by Starr
607 !---- & Cox (1985).
608 !
609 !---- Capacitance & characteristic lengths for Reynolds Number calculations
610 !---- are based on Young (1993; p. 144 & p. 150). c-axis & a-axis
611 !---- relationships are from Heymsfield (JAS, 1972; Table 1, p. 1351).
612 !
613  icount=60
614 !
615  if (print_diag) &
616  & write(7,"(2a)") '---- Increase in fall speeds of rimed ice', &
617  & ' particles as function of ice particle diameter ----'
618  do i=idimin,idimax
619  if (icount == 60 .and. print_diag) then
620  write(6,"(/2a/3a)") 'Particle masses (mg), fall speeds ', &
621  & '(m/s), and ventilation factors', &
622  & ' D(mm) CR_mass Mass_bull Mass_col Mass_plat ', &
623  & ' Mass_agg CR_vel V_bul CR_col CR_pla Aggreg', &
624  & ' Vent1 Vent2 '
625  write(7,"(3a)") ' <----------------------------------',&
626  & '--------------- Rime Factor --------------------------', &
627  & '--------------------------->'
628  write(7,"(a,23f5.2)") ' D(mm)',(rime_factor(k), k=1,5), &
629  & (rime_factor(k), k=6,40,2)
630  icount=0
631  endif
632  d=(float(i)+.5)*1.e-3 ! in mm
633  c_avg=0.
634  c_agg=0.
635  c_bullet=0.
636  c_column=0.
637  c_plate=0.
638  cl_agg=0.
639  cl_bullet=0.
640  cl_column=0.
641  cl_plate=0.
642  m_agg=0.
643  m_bullet=0.
644  m_column=0.
645  m_plate=0.
646  v_agg=0.
647  v_bullet=0.
648  v_column=0.
649  v_plate=0.
650  if (d < d_crystal_max) then
651 !
652 !---- This block of code calculates bulk characteristics based on average
653 ! characteristics of bullets, plates, & column ice crystals <1.5 mm size
654 !
655 !---- Mass-diameter relationships from Heymsfield (1972) & used
656 ! in Starr & Cox (1985), units in mg
657 !---- "d" is maximum dimension size of crystal in mm,
658 !
659 ! Mass of pure ice for spherical particles, used as an upper limit for the
660 ! mass of small columns (<~ 80 microns) & plates (<~ 35 microns)
661 !
662  m_ice=.48*d**3 ! Mass of pure ice for spherical particle
663 !
664  m_bullet=min(.044*d**3, m_ice)
665  m_column=min(.017*d**1.7, m_ice)
666  m_plate=min(.026*d**2.5, m_ice)
667 !
668  mass(i)=m_bullet+m_column+m_plate
669 !
670 !---- These relationships are from Starr & Cox (1985), applicable at 400 mb
671 !---- "d" is maximum dimension size of crystal in mm, dx in microns
672 !
673  dx=1000.*d ! Convert from mm to microns
674  if (dx <= 200.) then
675  v_column=8.114e-5*dx**1.585
676  v_bullet=5.666e-5*dx**1.663
677  v_plate=1.e-3*dx
678  else if (dx <= 400.) then
679  v_column=4.995e-3*dx**.807
680  v_bullet=3.197e-3*dx**.902
681  v_plate=1.48e-3*dx**.926
682  else if (dx <= 600.) then
683  v_column=2.223e-2*dx**.558
684  v_bullet=2.977e-2*dx**.529
685  v_plate=9.5e-4*dx
686  else if (dx <= 800.) then
687  v_column=4.352e-2*dx**.453
688  v_bullet=2.144e-2*dx**.581
689  v_plate=3.161e-3*dx**.812
690  else
691  v_column=3.833e-2*dx**.472
692  v_bullet=3.948e-2*dx**.489
693  v_plate=7.109e-3*dx**.691
694  endif
695 !
696 !---- Reduce fall speeds from 400 mb to 1000 mb
697 !
698  v_column=px*v_column
699  v_bullet=px*v_bullet
700  v_plate=px*v_plate
701 !
702 !---- DIFFERENT VERSION! CALCULATES MASS-WEIGHTED CRYSTAL FALL SPEEDS
703 !
704  vel(i)=(m_bullet*v_bullet+m_column*v_column+m_plate*v_plate)/ &
705  & mass(i)
706  mass(i)=mass(i)/3.
707 !
708 !---- Shape factor and characteristic length of various ice habits,
709 ! capacitance is equal to 4*PI*(Shape factor)
710 ! See Young (1993, pp. 143-152 for guidance)
711 !
712 !---- Bullets:
713 !
714 !---- Shape factor for bullets (Heymsfield, 1975)
715  c_bullet=.5*d
716 !---- Length-width functions for bullets from Heymsfield (JAS, 1972)
717  if (d > 0.3) then
718  wd=.25*d**.7856 ! Width (mm); a-axis
719  else
720  wd=.185*d**.552
721  endif
722 !---- Characteristic length for bullets (see first multiplicative term on right
723 ! side of eq. 7 multiplied by crystal width on p. 821 of Heymsfield, 1975)
724  cl_bullet=.5*pi*wd*(.25*wd+d)/(d+wd)
725 !
726 !---- Plates:
727 !
728 !---- Length-width function for plates from Heymsfield (JAS, 1972)
729  wd=.0449*d**.449 ! Width or thickness (mm); c-axis
730 !---- Eccentricity & shape factor for thick plates following Young (1993, p. 144)
731  ecc_plate=sqrt(1.-wd*wd/(d*d)) ! Eccentricity
732  c_plate=d*ecc_plate/asin(ecc_plate) ! Shape factor
733 !---- Characteristic length for plates following Young (1993, p. 150, eq. 6.6)
734  cl_plate=d+2.*wd ! Characteristic lengths for plates
735 !
736 !---- Columns:
737 !
738 !---- Length-width function for columns from Heymsfield (JAS, 1972)
739  if (d > 0.2) then
740  wd=.1973*d**.414 ! Width (mm); a-axis
741  else
742  wd=.5*d ! Width (mm); a-axis
743  endif
744 !---- Eccentricity & shape factor for columns following Young (1993, p. 144)
745  ecc_column=sqrt(1.-wd*wd/(d*d)) ! Eccentricity
746  c_column=ecc_column*d/log((1.+ecc_column)*d/wd) ! Shape factor
747 !---- Characteristic length for columns following Young (1993, p. 150, eq. 6.7)
748  cl_column=(wd+2.*d)/(c1+c2*d/wd) ! Characteristic lengths for columns
749 !
750 !---- Convert shape factor & characteristic lengths from mm to m for
751 ! ventilation calculations
752 !
753  c_bullet=.001*c_bullet
754  c_plate=.001*c_plate
755  c_column=.001*c_column
756  cl_bullet=.001*cl_bullet
757  cl_plate=.001*cl_plate
758  cl_column=.001*cl_column
759 !
760 !---- Make a smooth transition between a ventilation coefficient of 1.0 at 50 microns
761 ! to 1.1 at 200 microns
762 !
763  if (d > 0.2) then
764  cvent1=cvent1i
765  cvent2=cvent2i/3.
766  else
767  cvent1=1.0+.1*max(0., d-.05)/.15
768  cvent2=0.
769  endif
770 !
771 !---- Ventilation factors for ice crystals:
772 !
773  vent1(i)=cvent1*(c_bullet+c_plate+c_column)/3.
774  vent2(i)=cvent2*(c_bullet*sqrt(v_bullet*cl_bullet) &
775  & +c_plate*sqrt(v_plate*cl_plate) &
776  & +c_column*sqrt(v_column*cl_column) )
777  crime_best=crime3 ! For calculating Best No. of rimed ice crystals
778  else
779 !
780 !---- This block of code calculates bulk characteristics based on average
781 ! characteristics of unrimed aggregates >= 1.5 mm using Locatelli &
782 ! Hobbs (JGR, 1974, 2185-2197) data.
783 !
784 !----- This category is a composite of aggregates of unrimed radiating
785 !----- assemblages of dendrites or dendrites; aggregates of unrimed
786 !----- radiating assemblages of plates, side planes, bullets, & columns;
787 !----- aggregates of unrimed side planes (mass in mg, velocity in m/s)
788 !
789  m_agg=(.073*d**1.4+.037*d**1.9+.04*d**1.4)/3.
790  v_agg=(.8*d**.16+.69*d**.41+.82*d**.12)/3.
791  mass(i)=m_agg
792  vel(i)=v_agg
793 !
794 !---- Assume spherical aggregates
795 !
796 !---- Shape factor is the same as for bullets, = D/2
797  c_agg=.001*.5*d ! Units of m
798 !---- Characteristic length is surface area divided by perimeter
799 ! (.25*PI*D**2)/(PI*D**2) = D/4
800  cl_agg=.5*c_agg ! Units of m
801 !
802 !---- Ventilation factors for aggregates:
803 !
804  vent1(i)=cvent1a*c_agg
805  vent2(i)=cvent2a*c_agg*sqrt(v_agg*cl_agg)
806  crime_best=crime4 ! For calculating Best No. of rimed aggregates
807  endif
808 !
809 !---- Convert from shape factor to capacitance for ventilation factors
810 !
811  vent1(i)=4.*pi*vent1(i)
812  vent2(i)=4.*pi*vent2(i)
813  diam(i)=1.e-3*d ! Convert from mm to m
814  mass(i)=1.e-6*mass(i) ! Convert from mg to kg
815 !
816 !---- Calculate increase in fall speeds of individual rimed ice particles
817 !
818  do k=0,nrime
819 !---- Mass of rimed ice particle associated with rime_factor(k)
820  rime_m1=rime_factor(k)*mass(i)
821  rime_m2=cbulk_ice*diam(i)**3
822  m_rime=min(rime_m1, rime_m2)
823 !---- Best Number (X) of rimed ice particle combining eqs. (8) & (12) in Bohm
824  x_rime=crime2*m_rime*(crime_best/m_rime)**.25
825 !---- Reynolds Number for rimed ice particle using eq. (11) in Bohm
826  re_rime=8.5*(sqrt(1.+.1519*sqrt(x_rime))-1.)**2
827  rime_vel(k)=crime5*re_rime/diam(i)
828  enddo
829  do k=1,nrime
830  vel_rime(i,k)=rime_vel(k)/rime_vel(0)
831  enddo
832  if (print_diag) then
833  !
834  !---- Determine if statistics should be printed out.
835  !
836  iprint=.false.
837  if (d <= 1.) then
838  if (mod(i,10) == 0) iprint=.true.
839  else
840  if (mod(i,100) == 0) iprint=.true.
841  endif
842  if (iprint) then
843  write(6,"(f7.4,5e11.4,1x,5f7.4,1x,2e11.4)") &
844  & d,1.e6*mass(i),m_bullet,m_column,m_plate,m_agg, &
845  & vel(i),v_bullet,v_column,v_plate,v_agg, &
846  & vent1(i),vent2(i)
847  write(7,"(f7.4,23f5.2)") d,(vel_rime(i,k), k=1,5), &
848  & (vel_rime(i,k), k=6,40,2)
849  icount=icount+1
850  endif
851  endif
852  enddo
853 !
854 !#######################################################################
855 ! Characteristics as functions of mean particle diameter
856 !#######################################################################
857 !
858  venti1=0.
859  venti2=0.
860  accri=0.
861  massi=0.
862  vsnowi=0.
863  vel_rf=0.
864  ivel_rime=0.
865  icount=0
866  if (print_diag) then
867  icount=60
868  write(6,"(/2a)") '------------- Statistics as functions of ', &
869  & 'mean particle diameter -------------'
870  write(7,"(/2a)") '------ Increase in fall speeds of rimed ice', &
871  & ' particles as functions of mean particle diameter -----'
872  endif
873  do j=mdimin,mdimax
874  if (icount == 60 .AND. print_diag) then
875  write(6,"(/2a)") 'D(mm) Vent1 Vent2 ', &
876  & 'Accrete Mass Vel Dens '
877  write(7,"(3a)") ' <----------------------------------', &
878  & '--------------- Rime Factor --------------------------', &
879  & '--------------------------->'
880  write(7,"(a,23f5.2)") 'D(mm)',(rime_factor(k), k=1,5), &
881  & (rime_factor(k), k=6,40,2)
882  icount=0
883  endif
884  mdiam=deldmi*float(j) ! in m
885  smom3=0.
886  pratei=0.
887  rime_vel=0. ! Note that this array is being reused!
888  do i=idimin,idimax
889  dx=diam(i)/mdiam
890  if (dx <= xmax) then ! To prevent arithmetic underflows
891  expf=exp(-dx)*ddeli
892  venti1(j)=venti1(j)+vent1(i)*expf
893  venti2(j)=venti2(j)+vent2(i)*expf
894  accri(j)=accri(j)+diam(i)*diam(i)*vel(i)*expf
895  xmass=mass(i)*expf
896  do k=1,nrime
897  rime_vel(k)=rime_vel(k)+xmass*vel_rime(i,k)
898  enddo
899  massi(j)=massi(j)+xmass
900  pratei=pratei+xmass*vel(i)
901  smom3=smom3+diam(i)**3*expf
902  else
903  exit
904  endif
905  enddo
906  !
907  !--- Increased fall velocities functions of mean diameter (j),
908  ! normalized by ice content, and rime factor (k)
909  !
910  do k=1,nrime
911  ivel_rime(j,k)=rime_vel(k)/massi(j)
912  enddo
913  !
914  !--- Increased fall velocities functions of ice content at 0.1 mm
915  ! intervals (j_100) and rime factor (k); accumulations here
916  !
917  jj=j/100
918  if (jj >= 2 .AND. jj <= 9) then
919  do k=1,nrime
920  vel_rf(jj,k)=vel_rf(jj,k)+ivel_rime(j,k)
921  enddo
922  endif
923  bulk_dens=cbulk*massi(j)/smom3
924  venti1(j)=venti1(j)/mdiam
925  venti2(j)=venti2(j)/mdiam
926  accri(j)=accri(j)/mdiam
927  vsnowi(j)=pratei/massi(j)
928  massi(j)=massi(j)/mdiam
929  if (mod(j,10) == 0 .AND. print_diag) then
930  xmdiam=1.e3*mdiam
931  write(6,"(f6.3,4e11.4,f6.3,f8.3)") xmdiam,venti1(j),venti2(j),&
932  & accri(j),massi(j),vsnowi(j),bulk_dens
933  write(7,"(f6.3,23f5.2)") xmdiam,(ivel_rime(j,k), k=1,5), &
934  & (ivel_rime(j,k), k=6,40,2)
935  icount=icount+1
936  endif
937  enddo
938 !
939 !--- Average increase in fall velocities rimed ice as functions of mean
940 ! particle diameter (j, only need 0.1 mm intervals) and rime factor (k)
941 !
942  if (print_diag) then
943  write(7,"(/2a)") ' ------- Increase in fall speeds of rimed ', &
944  & 'ice particles at reduced, 0.1-mm intervals --------'
945  write(7,"(3a)") ' <----------------------------------', &
946  & '--------------- Rime Factor --------------------------', &
947  & '--------------------------->'
948  write(7,"(a,23f5.2)") 'D(mm)',(rime_factor(k), k=1,5), &
949  & (rime_factor(k), k=6,40,2)
950  endif
951  do j=2,9
952  vel_rf(j,0)=1.
953  do k=1,nrime
954  vel_rf(j,k)=.01*vel_rf(j,k)
955  enddo
956  if (print_diag) write(7,"(f4.1,2x,23f5.2)") 0.1*j, &
957  & (vel_rf(j,k), k=1,5),(vel_rf(j,k), k=6,40,2)
958  enddo
959 !
960 !-----------------------------------------------------------------------------------
961 !
962  end subroutine ice_lookup
963 !
964 !#######################################################################
965 !-------------- Creates lookup tables for rain processes ---------------
966 !#######################################################################
967 !
968  subroutine rain_lookup
969  implicit none
970 !
971 !--- Parameters & arrays for fall speeds of rain as a function of rain drop
972 ! diameter. These quantities are integrated over exponential size
973 ! distributions of rain drops at 1 micron intervals (DdelR) from minimum
974 ! drop sizes of .05 mm (50 microns, DminR) to maximum drop sizes of 10 mm
975 ! (DmaxR).
976 !
977  real, parameter :: DminR=.05e-3, dmaxr=10.e-3, ddelr=1.e-6, &
978  & xrmin=1.e6*dminr, xrmax=1.e6*dmaxr
979  integer, parameter :: IDRmin=xrmin, idrmax=xrmax
980  real diam(idrmin:idrmax), vel(idrmin:idrmax)
981 !
982 !--- Parameters rain lookup tables, which establish the range of mean drop
983 ! diameters; from a minimum mean diameter of 0.05 mm (DMRmin) to a
984 ! maximum mean diameter of 0.45 mm (DMRmax). The tables store solutions
985 ! at 1 micron intervals (DelDMR) of mean drop diameter.
986 !
987  real mdiam, mass
988 !
989  logical, parameter :: print_diag=.false.
990 !
991  real d, cmass, pi2, expf
992  integer i, j, i1, i2
993 !
994 !-----------------------------------------------------------------------
995 !------- Fall speeds of rain as function of rain drop diameter ---------
996 !-----------------------------------------------------------------------
997 !
998  do i=idrmin,idrmax
999  diam(i)=float(i)*ddelr
1000  d=100.*diam(i) ! Diameter in cm
1001  if (d <= .42) then
1002  !
1003  !--- Rutledge & Hobbs (1983); vel (m/s), d (cm)
1004  !
1005  vel(i)=max(0., -.267+51.5*d-102.25*d*d+75.7*d**3)
1006  else if (d > 0.42 .and. d <= .58) then
1007  !
1008  !--- Linear interpolation of Gunn & Kinzer (1949) data
1009  !
1010  vel(i)=8.92+.25/(.58-.42)*(d-.42)
1011  else
1012  vel(i)=9.17
1013  endif
1014  enddo
1015  do i=1,100
1016  i1=(i-1)*100+idrmin
1017  i2=i1+90
1018  !
1019  !--- Print out rain fall speeds only for D<=5.8 mm (.58 cm)
1020  !
1021  if (diam(i1) > .58e-2) exit
1022  if (print_diag) then
1023  write(6,"(1x)")
1024  write(6,"('D(mm)-> ',10f7.3)") (1000.*diam(j), j=i1,i2,10)
1025  write(6,"('V(m/s)-> ',10f7.3)") (vel(j), j=i1,i2,10)
1026  endif
1027  enddo
1028 !
1029 !-----------------------------------------------------------------------
1030 !------------------- Lookup tables for rain processes ------------------
1031 !-----------------------------------------------------------------------
1032 !
1033 ! pi=acos(-1.)
1034  pi2=2.*pi
1035  cmass=1000.*pi/6.
1036  if (print_diag) then
1037  write(6,"(/'Diam - Mean diameter (mm)' &
1038  & /'VENTR1 - 1st ventilation coefficient (m**2)' &
1039  & /'VENTR2 - 2nd ventilation coefficient (m**3/s**.5)' &
1040  & /'ACCRR - accretion moment (m**4/s)' &
1041  & /'RHO*QR - mass content (kg/m**3) for N0r=8e6' &
1042  & /'RRATE - rain rate moment (m**5/s)' &
1043  & /'VR - mass-weighted rain fall speed (m/s)' &
1044  & /' Diam VENTR1 VENTR2 ACCRR ', &
1045  & 'RHO*QR RRATE VR')")
1046  endif
1047  do j=mdrmin,mdrmax
1048  mdiam=float(j)*deldmr
1049  ventr2(j)=0.
1050  accrr(j)=0.
1051  massr(j)=0.
1052  rrate(j)=0.
1053  do i=idrmin,idrmax
1054  expf=exp(-diam(i)/mdiam)*ddelr
1055  ventr2(j)=ventr2(j)+diam(i)**1.5*vel(i)**.5*expf
1056  accrr(j)=accrr(j)+diam(i)*diam(i)*vel(i)*expf
1057  massr(j)=massr(j)+diam(i)**3*expf
1058  rrate(j)=rrate(j)+diam(i)**3*vel(i)*expf
1059  enddo
1060  !
1061  !--- Derived based on ventilation, F(D)=0.78+.31*Schmidt**(1/3)*Reynold**.5,
1062  ! where Reynold=(V*D*rho/dyn_vis), V is velocity, D is particle diameter,
1063  ! rho is air density, & dyn_vis is dynamic viscosity. Only terms
1064  ! containing velocity & diameter are retained in these tables.
1065  !
1066  ventr1(j)=.78*pi2*mdiam**2
1067  ventr2(j)=.31*pi2*ventr2(j)
1068  !
1069  massr(j)=cmass*massr(j)
1070  rrate(j)=cmass*rrate(j)
1071  vrain(j)=rrate(j)/massr(j)
1072  if (print_diag) write(6,"(f6.3,5g12.5,f6.3)") 1000.*mdiam, &
1073  & ventr1(j),ventr2(j),accrr(j),8.e6*massr(j),rrate(j),vrain(j)
1074  enddo
1075 !
1076 !-----------------------------------------------------------------------
1077 !
1078  end subroutine rain_lookup
1079 !
1080 !###############################################################################
1081 ! ***** VERSION OF MICROPHYSICS DESIGNED FOR HIGHER RESOLUTION MESO ETA MODEL
1082 ! (1) Represents sedimentation by preserving a portion of the precipitation
1083 ! through top-down integration from cloud-top. Modified procedure to
1084 ! Zhao and Carr (1997).
1085 ! (2) Microphysical equations are modified to be less sensitive to time
1086 ! steps by use of Clausius-Clapeyron equation to account for changes in
1087 ! saturation mixing ratios in response to latent heating/cooling.
1088 ! (3) Prevent spurious temperature oscillations across 0C due to
1089 ! microphysics.
1090 ! (4) Uses lookup tables for: calculating two different ventilation
1091 ! coefficients in condensation and deposition processes; accretion of
1092 ! cloud water by precipitation; precipitation mass; precipitation rate
1093 ! (and mass-weighted precipitation fall speeds).
1094 ! (5) Assumes temperature-dependent variation in mean diameter of large ice
1095 ! (Houze et al., 1979; Ryan et al., 1996).
1096 ! -> 8/22/01: This relationship has been extended to colder temperatures
1097 ! to parameterize smaller large-ice particles down to mean sizes of MDImin,
1098 ! which is 50 microns reached at -55.9C.
1099 ! (6) Attempts to differentiate growth of large and small ice, mainly for
1100 ! improved transition from thin cirrus to thick, precipitating ice
1101 ! anvils.
1102 ! -> 8/22/01: This feature has been diminished by effectively adjusting to
1103 ! ice saturation during depositional growth at temperatures colder than
1104 ! -10C. Ice sublimation is calculated more explicitly. The logic is
1105 ! that sources of are either poorly understood (e.g., nucleation for NWP)
1106 ! or are not represented in the Eta model (e.g., detrainment of ice from
1107 ! convection). Otherwise the model is too wet compared to the radiosonde
1108 ! observations based on 1 Feb - 18 March 2001 retrospective runs.
1109 ! (7) Top-down integration also attempts to treat mixed-phase processes,
1110 ! allowing a mixture of ice and water. Based on numerous observational
1111 ! studies, ice growth is based on nucleation at cloud top &
1112 ! subsequent growth by vapor deposition and riming as the ice particles
1113 ! fall through the cloud. Effective nucleation rates are a function
1114 ! of ice supersaturation following Meyers et al. (JAM, 1992).
1115 ! -> 8/22/01: The simulated relative humidities were far too moist compared
1116 ! to the rawinsonde observations. This feature has been substantially
1117 ! diminished, limited to a much narrower temperature range of 0 to -10C.
1118 ! (8) Depositional growth of newly nucleated ice is calculated for large time
1119 ! steps using Fig. 8 of Miller and Young (JAS, 1979), at 1 deg intervals
1120 ! using their ice crystal masses calculated after 600 s of growth in water
1121 ! saturated conditions. The growth rates are normalized by time step
1122 ! assuming 3D growth with time**1.5 following eq. (6.3) in Young (1993).
1123 ! -> 8/22/01: This feature has been effectively limited to 0 to -10C.
1124 ! (9) Ice precipitation rates can increase due to increase in response to
1125 ! cloud water riming due to (a) increased density & mass of the rimed
1126 ! ice, and (b) increased fall speeds of rimed ice.
1127 ! -> 8/22/01: This feature has been effectively limited to 0 to -10C.
1128 !###############################################################################
1129 !###############################################################################
1130 !
1131  SUBROUTINE gsmcolumn ( ARAING, ASNOWG, DTPG, I_index, J_index, &
1132  & lsfc, p_col, qi_col, qr_col, qv_col, qw_col, rimef_col, t_col, &
1133  & thick_col, wc_col, lm, rhc_col, xncw, flgmin, print_diag, psfc)
1134 !
1135  implicit none
1136 !
1137 !###############################################################################
1138 !###############################################################################
1139 !
1140 !-------------------------------------------------------------------------------
1141 !----- NOTE: In this version of the Code threading should be done outside!
1142 !-------------------------------------------------------------------------------
1143 !$$$ SUBPROGRAM DOCUMENTATION BLOCK
1144 ! . . .
1145 ! SUBPROGRAM: Grid-scale microphysical processes - condensation & precipitation
1146 ! PRGRMMR: Ferrier ORG: W/NP22 DATE: 08-2001
1147 ! Updated: Moorthi for GFS application
1148 !-------------------------------------------------------------------------------
1149 ! ABSTRACT:
1150 ! * Merges original GSCOND & PRECPD subroutines.
1151 ! * Code has been substantially streamlined and restructured.
1152 ! * Exchange between water vapor & small cloud condensate is calculated using
1153 ! the original Asai (1965, J. Japan) algorithm. See also references to
1154 ! Yau and Austin (1979, JAS), Rutledge and Hobbs (1983, JAS), and Tao et al.
1155 ! (1989, MWR). This algorithm replaces the Sundqvist et al. (1989, MWR)
1156 ! parameterization.
1157 !-------------------------------------------------------------------------------
1158 !
1159 ! USAGE:
1160 ! * CALL GSMCOLUMN FROM SUBROUTINE GSMDRIVE
1161 ! * SUBROUTINE GSMDRIVE CALLED FROM MAIN PROGRAM EBU
1162 !
1163 ! INPUT ARGUMENT LIST:
1164 ! DTPH - physics time step (s)
1165 ! I_index - I index
1166 ! J_index - J index
1167 ! LSFC - Eta level of level above surface, ground
1168 ! P_col - vertical column of model pressure (Pa)
1169 ! QI_col - vertical column of model ice mixing ratio (kg/kg)
1170 ! QR_col - vertical column of model rain ratio (kg/kg)
1171 ! QV_col - vertical column of model water vapor specific humidity (kg/kg)
1172 ! QW_col - vertical column of model cloud water mixing ratio (kg/kg)
1173 ! RimeF_col - vertical column of rime factor for ice in model (ratio, defined below)
1174 ! T_col - vertical column of model temperature (deg K)
1175 ! THICK_col - vertical column of model mass thickness (density*height increment)
1176 ! WC_col - vertical column of model mixing ratio of total condensate (kg/kg)
1177 !
1178 !
1179 ! OUTPUT ARGUMENT LIST:
1180 ! ARAIN - accumulated rainfall at the surface (kg)
1181 ! ASNOW - accumulated snowfall at the surface (kg)
1182 ! QV_col - vertical column of model water vapor specific humidity (kg/kg)
1183 ! WC_col - vertical column of model mixing ratio of total condensate (kg/kg)
1184 ! QW_col - vertical column of model cloud water mixing ratio (kg/kg)
1185 ! QI_col - vertical column of model ice mixing ratio (kg/kg)
1186 ! QR_col - vertical column of model rain ratio (kg/kg)
1187 ! RimeF_col - vertical column of rime factor for ice in model (ratio, defined below)
1188 ! T_col - vertical column of model temperature (deg K)
1189 !
1190 ! OUTPUT FILES:
1191 ! NONE
1192 !
1193 ! Subprograms & Functions called:
1194 ! * Real Function CONDENSE - cloud water condensation
1195 ! * Real Function DEPOSIT - ice deposition (not sublimation)
1196 !
1197 ! UNIQUE: NONE
1198 !
1199 ! LIBRARY: NONE
1200 !
1201 ! ATTRIBUTES:
1202 ! LANGUAGE: FORTRAN 90
1203 ! MACHINE : IBM SP
1204 !
1205 !-------------------------------------------------------------------------
1206 !--------------- Arrays & constants in argument list ---------------------
1207 !-------------------------------------------------------------------------
1208 !
1209  integer lm
1210  REAL ARAING, ASNOWG, P_col(lm), QI_col(lm), QR_col(lm), QV_col(lm)&
1211  &, QW_col(LM), RimeF_col(LM), T_col(LM), THICK_col(LM), &
1212  & WC_col(LM), RHC_col(LM), XNCW(LM), ARAIN, ASNOW, dtpg, psfc
1213  real flgmin
1214 !
1215  INTEGER I_index, J_index, LSFC
1216 !
1217 !
1218 !-------------------------------------------------------------------------
1219 !
1220 !--- Mean ice-particle diameters varying from 50 microns to 1000 microns
1221 ! (1 mm), assuming an exponential size distribution.
1222 !
1223 !---- Meaning of the following arrays:
1224 ! - mdiam - mean diameter (m)
1225 ! - VENTI1 - integrated quantity associated w/ ventilation effects
1226 ! (capacitance only) for calculating vapor deposition onto ice
1227 ! - VENTI2 - integrated quantity associated w/ ventilation effects
1228 ! (with fall speed) for calculating vapor deposition onto ice
1229 ! - ACCRI - integrated quantity associated w/ cloud water collection by ice
1230 ! - MASSI - integrated quantity associated w/ ice mass
1231 ! - VSNOWI - mass-weighted fall speed of snow (large ice), used to calculate
1232 ! precipitation rates
1233 !
1234  REAL, PARAMETER :: DMImin=.05e-3, dmimax=1.e-3, deldmi=1.e-6, &
1235  & xmimin=1.e6*dmimin, xmimax=1.e6*dmimax
1236  INTEGER, PARAMETER :: MDImin=xmimin, mdimax=xmimax
1237 !
1238 !-------------------------------------------------------------------------
1239 !------- Key parameters, local variables, & important comments ---------
1240 !-----------------------------------------------------------------------
1241 !
1242 !--- KEY Parameters:
1243 !
1244 !---- Comments on 14 March 2002
1245 ! * Set EPSQ to the universal value of 1.e-12 throughout the code
1246 ! condensate. The value of EPSQ will need to be changed in the other
1247 ! subroutines in order to make it consistent throughout the Eta code.
1248 ! * Set CLIMIT=10.*EPSQ as the lower limit for the total mass of
1249 ! condensate in the current layer and the input flux of condensate
1250 ! from above (TOT_ICE, TOT_ICEnew, TOT_RAIN, and TOT_RAINnew).
1251 !
1252 !-- NLImax - maximum number concentration of large ice crystals (20,000 /m**3, 20 per liter)
1253 !-- NLImin - minimum number concentration of large ice crystals (100 /m**3, 0.1 per liter)
1254 !
1255  REAL, PARAMETER :: RHOL=1000., xls=hvap+hfus &
1256 
1257 ! &, T_ICE=-10. !- Ver1
1258 ! &, T_ICE_init=-5. !- Ver1
1259 !!! &, T_ICE=-20. !- Ver2
1260 ! &, T_ICE=-40. !- Ver2
1261 ! &, T_ICE_init=-15., !- Ver2
1262 !
1263 ! & CLIMIT=10.*EPSQ, EPS1=RV/RD-1., RCP=1./CP,
1264 
1265  &,climit=10.*epsq, rcp=1./cp, &
1266  & rcprv=rcp/rv, rrhol=1./rhol, xls1=xls*rcp, xls2=xls*xls*rcprv, &
1267  & xls3=xls*xls/rv, &
1268  & c1=1./3., c2=1./6., c3=3.31/6., &
1269  & dmr1=.1e-3, dmr2=.2e-3, dmr3=.32e-3, n0r0=8.e6, n0rmin=1.e4, &
1270  & n0s0=4.e6, rho0=1.194, xmr1=1.e6*dmr1, xmr2=1.e6*dmr2, &
1271  & xmr3=1.e6*dmr3, xratio=.025
1272  INTEGER, PARAMETER :: MDR1=xmr1, mdr2=xmr2, mdr3=xmr3
1273 !
1274 !--- If BLEND=1:
1275 ! precipitation (large) ice amounts are estimated at each level as a
1276 ! blend of ice falling from the grid point above and the precip ice
1277 ! present at the start of the time step (see TOT_ICE below).
1278 !--- If BLEND=0:
1279 ! precipitation (large) ice amounts are estimated to be the precip
1280 ! ice present at the start of the time step.
1281 !
1282 !--- Extended to include sedimentation of rain on 2/5/01
1283 !
1284  REAL, PARAMETER :: BLEND=1.
1285 !
1286 !--- This variable is for debugging purposes (if .true.)
1287 !
1288  LOGICAL PRINT_diag
1289 !
1290 !--- Local variables
1291 !
1292  REAL EMAIRI, N0r, NLICE, NSmICE, NLImax, pfac
1293  LOGICAL CLEAR, ICE_logical, DBG_logical, RAIN_logical
1294 
1295  integer lbef, ipass, ixrf, ixs, itdx, idr &
1296  &, index_my, indexr, indexr1, indexs &
1297  &, i, j, k, l, ntimes, item
1298 ! &, i, j, k, my_600, i1, i2, l, ntimes
1299 
1300  real flimass, xlimass, vsnow, qi_min, dum, piloss &
1301  &, tot_ice, xsimass, vel_inc, vrimef, rimef1, dum1 &
1302  &, dum2, fws, denomi, dwv &
1303  &, xrf, qw0, dli, xli, fsmall &
1304  &, prevp, tk2, dtph &
1305  &, pievp, picnd, piacr, pracw &
1306  &, praut, pimlt, qtice, qlice &
1307  &, gammar, flarge, wvqw, dynvis &
1308  &, tfactor, denom, gammas, diffus, therm_cond &
1309  &, wvnew, delv, tnew, tot_icenew, rimef &
1310  &, deli, fwr, crevp, ventr, delt &
1311  &, delw, fir, delr, qsinew, qswnew &
1312  &, budget, wsnew, vrain2, tot_rainnew &
1313  &, qtnew, qt, wcnew, abw &
1314  &, aievp, tcc, denomf, abi &
1315  &, sfactor, pidep_max, didep, ventis, ventil &
1316  &, dievp, rqr, rfactor, dwvr, rr, tot_rain &
1317  &, dwv0, qsw0, prloss, qtrain, vrain1 &
1318  &, qsw, ws, esi, esw, wv, wc, rhgrd, rho &
1319  &, rrho, dtrho, wsgrd, qsi, qswgrd, qsigrd &
1320  &, tk, tc, pp, bldtrh &
1321  &, xlv, xlv1, xlf, xlf1, xlv2, denomw, denomwi &
1322  &, qwnew, pcond, pidep, qrnew, qi, qr, qw &
1323  &, piacw, piacwi, piacwr, qv, dwvi &
1324  &, arainnew, thick, asnownew &
1325  &, qinew, qi_min_0c, QSW_l, QSI_l, QSW0_l, SCHMIT_FAC
1326 
1327 !
1328 !
1329 !#######################################################################
1330 !########################## Begin Execution ############################
1331 !#######################################################################
1332 !
1333  dtph = dtpg / mic_step
1334  araing = 0. ! Total Accumulated rainfall at surface (kg/m**2)
1335  asnowg = 0. ! Total Accumulated snowfall at surface (kg/m**2)
1336 !
1337  do ntimes =1,mic_step
1338 !
1339  qi_min_0c = 10.e3*massi(mdimin) !- Ver5
1340  arain = 0. ! Accumulated rainfall at surface for this step (kg/m**2)
1341  asnow = 0. ! Accumulated snowfall at surface for this step (kg/m**2)
1342 
1343  indexr = mdrmin
1344 !
1345 !-----------------------------------------------------------------------
1346 !
1347  DO l=1,lsfc ! Loop from top (L=1) to surface (L=LSFC)
1348 
1349 !--- Skip this level and go to the next lower level if no condensate
1350 ! and very low specific humidities
1351 !
1352  IF (qv_col(l) > epsq .OR. wc_col(l) > epsq) THEN
1353 !
1354 !-----------------------------------------------------------------------
1355 !------------ Proceed with cloud microphysics calculations -------------
1356 !-----------------------------------------------------------------------
1357 !
1358  tk = t_col(l) ! Temperature (deg K)
1359  tc = tk-t0c ! Temperature (deg C)
1360  pp = p_col(l) ! Pressure (Pa)
1361  qv = qv_col(l) ! Specific humidity of water vapor (kg/kg)
1362 ! WV = QV/(1.-QV) ! Water vapor mixing ratio (kg/kg)
1363  wv = qv ! Water vapor specific humidity (kg/kg)
1364  wc = wc_col(l) ! Grid-scale mixing ratio of total condensate
1365  ! (water or ice; kg/kg)
1366 ! WC = WC/(1.-WC)
1367  rhgrd = rhc_col(l)
1368 
1369 !
1370 ! Pressure dependen scaling factor for flgmin (tunable)
1371 !
1372 !!! pfac = max(0.5, (min(1.0, pp*0.00002))**2) ! commented on 02182011
1373 !go pfac = max(0.5, (sqrt(min(1.0, pp*0.00004))))
1374  pfac = 1.0
1375 !
1376  clear = .true.
1377 !
1378 !--- Check grid-scale saturation when no condensate is present
1379 !
1380  esw = min(pp, fpvsl(tk)) ! Saturation vapor pressure w/r/t water
1381 ! QSW = EPS*ESW/(PP-ESW) ! Saturation mixing ratio w/r/t water
1382  qsw = eps*esw/(pp+epsm1*esw) ! Saturation specific humidity w/r/t water
1383  ws = qsw ! General saturation mixing ratio (water/ice)
1384  qsi = qsw
1385  IF (tc < 0.) THEN
1386  esi = min(pp,fpvsi(tk)) ! Saturation vapor pressure w/r/t ice
1387 ! QSI = EPS*ESI/(PP-ESI) ! Saturation mixing ratio w/r/t water
1388  qsi = eps*esi/(pp+epsm1*esi) ! Saturation specific humidity w/r/t water
1389  ws = qsi ! General saturation mixing ratio (water/ice)
1390  if (pp <= esi) ws = wv / rhgrd
1391  ENDIF
1392 !
1393  dum = min(pp, esw0)
1394  qsw0 = eps*dum/(pp+epsm1*dum) ! Saturation specific Humidity at 0C
1395 !
1396 !--- Effective grid-scale Saturation mixing ratios
1397 !
1398  qswgrd = rhgrd*qsw
1399  qsigrd = rhgrd*qsi
1400  wsgrd = rhgrd*ws
1401  qsw_l = qswgrd
1402  qsi_l = qsigrd
1403  qsw0_l = qsw0*rhgrd
1404 !
1405 !--- Check if air is subsaturated and w/o condensate
1406 !
1407  IF (wv > wsgrd .OR. wc > epsq) clear = .false. ! Cloudy case
1408  IF (arain > climit) THEN ! If any rain is falling into layer from above
1409  clear = .false.
1410  ELSE
1411  arain = 0.
1412  ENDIF
1413 !
1414 !--- Check if any ice is falling into layer from above
1415 !
1416 !--- NOTE that "SNOW" in variable names is synonomous with
1417 ! large, precipitation ice particles
1418 !
1419  IF (asnow > climit) THEN
1420  clear = .false.
1421  ELSE
1422  asnow = 0.
1423  ENDIF
1424 !
1425 !-----------------------------------------------------------------------
1426 !-- Loop to the end if in clear, subsaturated air free of condensate ---
1427 !-----------------------------------------------------------------------
1428 !
1429  IF (.not. clear) THEN
1430 !
1431 !-----------------------------------------------------------------------
1432 !--------- Initialize RHO, THICK & microphysical processes -------------
1433 !-----------------------------------------------------------------------
1434 !
1435 !
1436 !--- Virtual temperature, TV=T*(1./EPS-1)*Q, Q is specific humidity;
1437 ! (see pp. 63-65 in Fleagle & Businger, 1963)
1438 !
1439  rho = pp/(rd*tk*(1.+eps1*qv)) ! Air density (kg/m**3)
1440  rrho = 1./rho ! Reciprocal of air density
1441  dtrho = dtph*rho ! Time step * air density
1442  bldtrh = blend*dtrho ! Blend parameter * time step * air density
1443  thick = thick_col(l) ! Layer thickness = RHO*DZ = -DP/G
1444 !
1445  arainnew = 0. ! Updated accumulated rainfall at surface
1446  asnownew = 0. ! Updated accumulated snowfall at surface
1447  qi = qi_col(l) ! Ice mixing ratio
1448  qinew = 0. ! Updated ice mixing ratio
1449  qr = qr_col(l) ! Rain mixing ratio
1450  qrnew = 0. ! Updated rain ratio
1451  qw = qw_col(l) ! Cloud water mixing ratio
1452  qwnew = 0. ! Updated cloud water ratio
1453 !
1454  pcond = 0. ! Condensation (>0) or evaporation (<0)
1455  ! of cloud water (kg/kg)
1456  pidep = 0. ! Deposition (>0) or sublimation (<0)
1457  ! of ice crystals (kg/kg)
1458  piacw = 0. ! Cloud water collection (riming)
1459  ! by precipitation ice (kg/kg; >0)
1460  piacwi = 0. ! Growth of precip ice by riming (kg/kg; >0)
1461  piacwr = 0. ! Shedding of accreted cloud water
1462  ! to form rain (kg/kg; >0)
1463  piacr = 0. ! Freezing of rain onto large ice
1464  ! at supercooled temps (kg/kg; >0)
1465  picnd = 0. ! Condensation (>0) onto wet, melting
1466  ! ice (kg/kg)
1467  pievp = 0. ! Evaporation (<0) from wet, melting
1468  ! ice (kg/kg)
1469  pimlt = 0. ! Melting ice (kg/kg; >0)
1470  praut = 0. ! Cloud water autoconversion to rain (kg/kg; >0)
1471  pracw = 0. ! Cloud water collection (accretion) by rain (kg/kg; >0)
1472  prevp = 0. ! Rain evaporation (kg/kg; <0)
1473 !
1474 !---------------------------------------------------------------------------
1475 !--- Double check input hydrometeor mixing ratios
1476 !
1477 ! DUM = WC - (QI+QW+QR)
1478 ! DUM1 = ABS(DUM)
1479 ! DUM2 = TOLER * MIN(WC, QI+QW+QR)
1480 ! IF (DUM1 > DUM2) THEN
1481 ! WRITE(6,"(/2(a,i4),a,i2)") '{@ i=',I_index,' j=',J_index,
1482 ! & ' L=',L
1483 ! WRITE(6,"(4(a12,g11.4,1x))")
1484 ! & '{@ TCold=',TC,'P=',.01*PP,'DIFF=',DUM,'WCold=',WC,
1485 ! & '{@ QIold=',QI,'QWold=',QW,'QRold=',QR
1486 ! ENDIF
1487 !
1488 !***********************************************************************
1489 !*********** MAIN MICROPHYSICS CALCULATIONS NOW FOLLOW! ****************
1490 !***********************************************************************
1491 !
1492 !--- Calculate a few variables, which are used more than once below
1493 !
1494 !--- Latent heat of vaporization as a function of temperature from
1495 ! Bolton (1980, JAS)
1496 !
1497  xlv = 3.148e6 - 2370*tk ! Latent heat of vaporization (Lv)
1498  xlf = xls-xlv ! Latent heat of fusion (Lf)
1499  xlv1 = xlv*rcp ! Lv/Cp
1500  xlf1 = xlf*rcp ! Lf/Cp
1501  tk2 = 1./(tk*tk) ! 1./TK**2
1502  xlv2 = xlv*xlv*qsw_l*tk2/rv ! Lv**2*Qsw_l/(Rv*TK**2)
1503  denomw = 1. + xlv2*rcp ! Denominator term, Clausius-Clapeyron correction
1504 !
1505 !--- Basic thermodynamic quantities
1506 ! * DYNVIS - dynamic viscosity [ kg/(m*s) ]
1507 ! * THERM_COND - thermal conductivity [ J/(m*s*K) ]
1508 ! * DIFFUS - diffusivity of water vapor [ m**2/s ]
1509 !
1510 ! TFACTOR = TK**1.5/(TK+120.)
1511  tfactor = tk*sqrt(tk)/(tk+120.)
1512  dynvis = 1.496e-6*tfactor
1513  therm_cond = 2.116e-3*tfactor
1514  diffus = 8.794e-5*tk**1.81/pp
1515  schmit_fac = (rho/(diffus*diffus*dynvis))**c2
1516 !
1517 !--- Air resistance term for the fall speed of ice following the
1518 ! basic research by Heymsfield, Kajikawa, others
1519 !
1520  gammas = (1.e5/pp)**c1
1521 !
1522 !--- Air resistance for rain fall speed (Beard, 1985, JAOT, p.470)
1523 !
1524  gammar = (rho0/rho)**0.4
1525 !
1526 !----------------------------------------------------------------------
1527 !------------- IMPORTANT MICROPHYSICS DECISION TREE -----------------
1528 !----------------------------------------------------------------------
1529 !
1530 !--- Determine if conditions supporting ice are present
1531 !
1532  IF (tc < 0. .OR. qi > epsq .OR. asnow > climit) THEN
1533  ice_logical = .true.
1534  ELSE
1535  ice_logical = .false.
1536  qlice = 0.
1537  qtice = 0.
1538  ENDIF
1539 !
1540 !--- Determine if rain is present
1541 !
1542  rain_logical = .false.
1543  IF (arain > climit .OR. qr > epsq) rain_logical = .true.
1544 !
1545  IF (ice_logical) THEN
1546 !
1547 !--- IMPORTANT: Estimate time-averaged properties.
1548 !
1549 !---
1550 ! * FLARGE - ratio of number of large ice to total (large & small) ice
1551 ! * FSMALL - ratio of number of small ice crystals to large ice particles
1552 ! -> Small ice particles are assumed to have a mean diameter of 50 microns.
1553 ! * XSIMASS - used for calculating small ice mixing ratio
1554 !---
1555 ! * TOT_ICE - total mass (small & large) ice before microphysics,
1556 ! which is the sum of the total mass of large ice in the
1557 ! current layer and the input flux of ice from above
1558 ! * PILOSS - greatest loss (<0) of total (small & large) ice by
1559 ! sublimation, removing all of the ice falling from above
1560 ! and the ice within the layer
1561 ! * RimeF1 - Rime Factor, which is the mass ratio of total (unrimed & rimed)
1562 ! ice mass to the unrimed ice mass (>=1)
1563 ! * VrimeF - the velocity increase due to rime factor or melting (ratio, >=1)
1564 ! * VSNOW - Fall speed of rimed snow w/ air resistance correction
1565 ! * EMAIRI - equivalent mass of air associated layer and with fall of snow into layer
1566 ! * XLIMASS - used for calculating large ice mixing ratio
1567 ! * FLIMASS - mass fraction of large ice
1568 ! * QTICE - time-averaged mixing ratio of total ice
1569 ! * QLICE - time-averaged mixing ratio of large ice
1570 ! * NLICE - time-averaged number concentration of large ice
1571 ! * NSmICE - number concentration of small ice crystals at current level
1572 !---
1573 !--- Assumed number fraction of large ice particles to total (large & small)
1574 ! ice particles, which is based on a general impression of the literature.
1575 !
1576  wvqw = wv + qw ! Water vapor + cloud water
1577 !
1578 !--- 6/19/03 - Deleted some code here ....
1579 !
1580 ! *********************************************************
1581 
1582 ! IF (TC >= 0. .OR. WVQW < QSIgrd) THEN
1583 ! !
1584 ! !--- Eliminate small ice particle contributions for melting & sublimation
1585 ! !
1586 ! FLARGE = FLARGE1
1587 ! ELSE
1588 ! !
1589 ! !--- Enhanced number of small ice particles during depositional growth
1590 ! ! (effective only when 0C > T >= T_ice [-10C] )
1591 ! !
1592 ! FLARGE = FLARGE2
1593 ! !
1594 ! !--- Larger number of small ice particles due to rime splintering
1595 ! !
1596 ! IF (TC >= -8. .AND. TC <= -3.) FLARGE=.5*FLARGE
1597 !
1598 ! ENDIF ! End IF (TC >= 0. .OR. WVQW < QSIgrd)
1599 ! FSMALL=(1.-FLARGE)/FLARGE
1600 ! XSIMASS=RRHO*MASSI(MDImin)*FSMALL
1601 ! *********************************************************
1602 !
1603  IF (qi <= epsq .AND. asnow <= climit) THEN
1604  indexs = mdimin
1605  flarge = 0. !--- Begin 6/19/03 changes
1606  fsmall = 1.
1607  xsimass = rrho*massi(mdimin) !--- End 6/19/03 changes
1608  tot_ice = 0.
1609  piloss = 0.
1610  rimef1 = 1.
1611  vrimef = 1.
1612  vel_inc = gammas
1613  vsnow = 0.
1614  emairi = thick
1615  xlimass = rrho*rimef1*massi(indexs)
1616  flimass = xlimass/(xlimass+xsimass)
1617  qlice = 0.
1618  qtice = 0.
1619  nlice = 0.
1620  nsmice = 0.
1621  ELSE
1622  !
1623  !--- For T<0C mean particle size follows Houze et al. (JAS, 1979, p. 160),
1624  ! converted from Fig. 5 plot of LAMDAs. Similar set of relationships
1625  ! also shown in Fig. 8 of Ryan (BAMS, 1996, p. 66).
1626  !
1627  !--- Begin 6/19/03 changes => allow NLImax to increase & FLARGE to
1628  ! decrease at COLDER temperatures; set FLARGE to zero (i.e., only small
1629  ! ice) if the ice mixing ratio is below QI_min
1630 
1631 ! DUM = MAX(0.05, MIN(1., EXP(.0536*TC)) )
1632  dum = max(0.05, min(1., exp(.0564*tc)) )
1633  indexs = min(mdimax, max(mdimin, int(xmimax*dum) ) )
1634 !
1635  dum = max(flgmin*pfac, dum)
1636 
1637  qi_min = qi_min_0c * dum !- Ver5 ----Not used ----
1638 !! QI_min = QI_min_0C !- Ver5
1639 !!! QI_min = QI_min_0C/DUM !- Ver5
1640 
1641  nlimax = 10.e3/sqrt(dum) !- Ver3
1642  IF (tc < 0.) THEN
1643  flarge = dum !- Ver4
1644  ELSE
1645  flarge = 1.
1646  ENDIF
1647  fsmall = (1.-flarge)/flarge
1648  xsimass = rrho*massi(mdimin)*fsmall
1649  tot_ice = thick*qi + blend*asnow
1650  piloss = -tot_ice/thick
1651  lbef = max(1,l-1)
1652  rimef1 = (rimef_col(l)*thick*qi &
1653  & + rimef_col(lbef)*blend*asnow)/tot_ice
1654  rimef1 = min(rimef1, rfmax)
1655  vsnow = 0.0
1656  DO ipass=0,1
1657  IF (rimef1 .LE. 1.) THEN
1658  rimef1 = 1.
1659  vrimef = 1.
1660  ELSE
1661  ixs = max(2, min(indexs/100, 9))
1662  xrf = 10.492*log(rimef1)
1663  ixrf = max(0, min(int(xrf), nrime))
1664  IF (ixrf .GE. nrime) THEN
1665  vrimef = vel_rf(ixs,nrime)
1666  ELSE
1667  vrimef = vel_rf(ixs,ixrf)+(xrf-float(ixrf))* &
1668  & (vel_rf(ixs,ixrf+1)-vel_rf(ixs,ixrf))
1669  ENDIF
1670  ENDIF ! End IF (RimeF1 <= 1.)
1671  vel_inc = gammas*vrimef
1672  vsnow = vel_inc*vsnowi(indexs)
1673  emairi = thick + bldtrh*vsnow
1674  xlimass = rrho*rimef1*massi(indexs)
1675  flimass = xlimass/(xlimass+xsimass)
1676  qtice = tot_ice/emairi
1677  qlice = flimass*qtice
1678  nlice = qlice/xlimass
1679  nsmice = fsmall*nlice
1680  !
1681  IF ( (nlice >= nlimin .AND. nlice <= nlimax) &
1682  & .OR. ipass == 1) THEN
1683  EXIT
1684  ELSE
1685  IF(tc < 0) THEN
1686  xli = rho*(qtice/dum-xsimass)/rimef1
1687  IF (xli <= massi(mdimin) ) THEN
1688  indexs = mdimin
1689  ELSE IF (xli <= massi(450) ) THEN
1690  dli = 9.5885e5*xli**.42066 ! DLI in microns
1691  indexs = min(mdimax, max(mdimin, int(dli) ) )
1692  ELSE IF (xli <= massi(mdimax) ) THEN
1693  dli = 3.9751e6*xli**.49870 ! DLI in microns
1694  indexs = min(mdimax, max(mdimin, int(dli) ) )
1695  ELSE
1696  indexs = mdimax
1697  ENDIF ! End IF (XLI <= MASSI(MDImin) )
1698  ENDIF ! End IF (TC < 0)
1699 !
1700 !--- Reduce excessive accumulation of ice at upper levels
1701 ! associated with strong grid-resolved ascent
1702 !
1703 !--- Force NLICE to be between NLImin and NLImax
1704 !
1705 !--- 8/22/01: Increase density of large ice if maximum limits
1706 ! are reached for number concentration (NLImax) and mean size
1707 ! (MDImax). Done to increase fall out of ice.
1708 !
1709 !
1710 
1711  dum = max(nlimin, min(nlimax, nlice) )
1712  IF (dum >= nlimax .AND. indexs >= mdimax) &
1713  & rimef1 = rho*(qtice/nlimax-xsimass)/massi(indexs)
1714 !
1715 ! WRITE(6,"(4(a12,g11.4,1x))")
1716 ! & '{$ TC=',TC,'P=',.01*PP,'NLICE=',NLICE,'DUM=',DUM,
1717 ! & '{$ XLI=',XLI,'INDEXS=',FLOAT(INDEXS),'RHO=',RHO,'QTICE=',QTICE,
1718 ! & '{$ XSIMASS=',XSIMASS,'RimeF1=',RimeF1
1719 
1720  ENDIF ! End IF ( (NLICE >=NLImin .AND. NLICE >= NLImax)
1721  ENDDO ! End DO IPASS=0,1
1722  ENDIF ! End IF (QI <= EPSQ .AND. ASNOW <= CLIMIT)
1723  ENDIF ! End IF (ICE_logical)
1724 !
1725 !----------------------------------------------------------------------
1726 !--------------- Calculate individual processes -----------------------
1727 !----------------------------------------------------------------------
1728 !
1729 !--- Cloud water autoconversion to rain and collection by rain
1730 !
1731  IF (qw > epsq .AND. tc >= t_ice) THEN
1732  !
1733  !--- QW0 could be modified based on land/sea properties,
1734  ! presence of convection, etc. This is why QAUT0 and CRAUT
1735  ! are passed into the subroutine as externally determined
1736  ! parameters. Can be changed in the future if desired.
1737  !
1738 ! QW0 = QAUT0*RRHO
1739  qw0 = qautx*rrho*xncw(l)
1740  praut = max(0., qw-qw0)*craut
1741  IF (qlice > epsq) THEN
1742  !
1743  !--- Collection of cloud water by large ice particles ("snow")
1744  ! PIACWI=PIACW for riming, PIACWI=0 for shedding
1745  !
1746 !Moor FWS = MIN(1., CIACW*VEL_INC*NLICE*ACCRI(INDEXS)/PP**C1) ! 20050422
1747  fws = min(0.1, ciacw*vel_inc*nlice*accri(indexs) &
1748  & /pp**c1)
1749  piacw = fws*qw
1750  IF (tc < 0.) piacwi = piacw ! Large ice riming
1751 
1752  ENDIF ! End IF (QLICE > EPSQ)
1753  ENDIF ! End IF (QW > EPSQ .AND. TC >= T_ICE)
1754 !
1755 !----------------------------------------------------------------------
1756 !--- Loop around some of the ice-phase processes if no ice should be present
1757 !----------------------------------------------------------------------
1758 !
1759  IF (ice_logical) THEN
1760 !
1761 !--- Now the pretzel logic of calculating ice deposition
1762 !
1763  IF (tc < t_ice .AND. (wv > qsigrd .OR. qw > epsq)) THEN
1764 !
1765 !--- Adjust to ice saturation at T<T_ICE (-10C) if supersaturated.
1766 ! Sources of ice due to nucleation and convective detrainment are
1767 ! either poorly understood, poorly resolved at typical NWP
1768 ! resolutions, or are not represented (e.g., no detrained
1769 ! condensate in BMJ Cu scheme).
1770 !
1771  pcond = -qw
1772  dum1 = tk + xlv1*pcond ! Updated (dummy) temperature (deg K)
1773  dum2 = wv+qw ! Updated (dummy) water vapor mixing ratio
1774  dum = min(pp,fpvsi(dum1)) ! Updated (dummy) saturation vapor pressure w/r/t ice
1775  dum = rhgrd*eps*dum/(pp+epsm1*dum) ! Updated (dummy) saturation specific humidity w/r/t ice
1776 ! DUM = RHgrd*EPS*DUM/(PP-DUM) ! Updated (dummy) saturation mixing ratio w/r/t ice
1777 
1778  IF (dum2 > dum) pidep = deposit(pp, rhgrd, dum1, dum2)
1779 
1780  dwvi = 0. ! Used only for debugging
1781 !
1782  ELSE IF (tc < 0.) THEN
1783 !
1784 !--- These quantities are handy for ice deposition/sublimation
1785 ! PIDEP_max - max deposition or minimum sublimation to ice saturation
1786 !
1787  denomi = 1. + xls2*qsi_l*tk2
1788  dwvi = min(wvqw,qsw_l)-qsi_l
1789  pidep_max = max(piloss, dwvi/denomi)
1790  IF (qtice > 0.) THEN
1791 !
1792 !--- Calculate ice deposition/sublimation
1793 ! * SFACTOR - [VEL_INC**.5]*[Schmidt**(1./3.)]*[(RHO/DYNVIS)**.5],
1794 ! where Schmidt (Schmidt Number) =DYNVIS/(RHO*DIFFUS)
1795 ! * Units: SFACTOR - s**.5/m ; ABI - m**2/s ; NLICE - m**-3 ;
1796 ! VENTIL, VENTIS - m**-2 ; VENTI1 - m ;
1797 ! VENTI2 - m**2/s**.5 ; DIDEP - unitless
1798 !
1799 ! SFACTOR = VEL_INC**.5*(RHO/(DIFFUS*DIFFUS*DYNVIS))**C2
1800  sfactor = sqrt(vel_inc)*schmit_fac
1801  abi = 1./(rho*xls3*qsi*tk2/therm_cond+1./diffus)
1802 !
1803 !--- VENTIL - Number concentration * ventilation factors for large ice
1804 !--- VENTIS - Number concentration * ventilation factors for small ice
1805 !
1806 !--- Variation in the number concentration of ice with time is not
1807 ! accounted for in these calculations (could be in the future).
1808 !
1809  ventil = (venti1(indexs) + sfactor*venti2(indexs)) &
1810  & * nlice
1811  ventis = (venti1(mdimin) + sfactor*venti2(mdimin)) &
1812  & * nsmice
1813  didep = abi*(ventil+ventis)*dtph
1814 !
1815 !--- Account for change in water vapor supply w/ time
1816 !
1817  IF (didep >= xratio) &
1818  & didep = (1.-exp(-didep*denomi))/denomi
1819  IF (dwvi > 0.) THEN
1820  pidep = min(dwvi*didep, pidep_max)
1821  ELSE IF (dwvi < 0.) THEN
1822  pidep = max(dwvi*didep, pidep_max)
1823  ENDIF
1824 !
1825  ELSE IF (wvqw > qsi_l .AND. tc <= t_ice_init) THEN
1826 !
1827 !--- Ice nucleation in near water-saturated conditions. Ice crystal
1828 ! growth during time step calculated using Miller & Young (1979, JAS).
1829 !--- These deposition rates could drive conditions below water saturation,
1830 ! which is the basis of these calculations. Intended to approximate
1831 ! more complex & computationally intensive calculations.
1832 !
1833  index_my = max(my_t1, min( int(.5-tc), my_t2 ) )
1834 !
1835 !--- DUM1 is the supersaturation w/r/t ice at water-saturated conditions
1836 !
1837 !--- DUM2 is the number of ice crystals nucleated at water-saturated
1838 ! conditions based on Meyers et al. (JAM, 1992).
1839 !
1840 !--- Prevent unrealistically large ice initiation (limited by PIDEP_max)
1841 ! if DUM2 values are increased in future experiments
1842 !
1843  dum1 = qsw/qsi - 1.
1844  dum2 = 1.e3*exp(12.96*dum1-0.639)
1845  pidep = min(pidep_max,dum2*my_growth(index_my)*rrho)
1846 !
1847  ENDIF ! End IF (QTICE > 0.)
1848 !
1849  ENDIF ! End IF (TC < T_ICE .AND. (WV > QSIgrd .OR. QW > EPSQ))
1850 !
1851 !------------------------------------------------------------------------
1852 !
1853  ENDIF ! End of IF (ICE_logical)then loop
1854 !
1855 !------------------------------------------------------------------------
1856 !
1857 !--- Cloud water condensation
1858 !
1859  IF (tc >= t_ice .AND. (qw > epsq .OR. wv > qswgrd)) THEN
1860  IF (piacwi == 0. .AND. pidep == 0.) THEN
1861  pcond = condense(pp, qw, rhgrd, tk, wv)
1862  ELSE !-- Modify cloud condensation in response to ice processes
1863  dum = xlv*qswgrd*rcprv*tk2
1864  denomwi = 1. + xls*dum
1865  denomf = xlf*dum
1866  dum = max(0., pidep)
1867  pcond = (wv-qswgrd-denomwi*dum-denomf*piacwi)/denomw
1868  dum1 = -qw
1869  dum2 = pcond - piacw
1870  IF (dum2 < dum1) THEN !--- Limit cloud water sinks
1871  dum = dum1/dum2
1872  pcond = dum*pcond
1873  piacw = dum*piacw
1874  piacwi = dum*piacwi
1875  ENDIF ! EndIF (DUM2 < DUM1)
1876  ENDIF ! EndIF (PIACWI == 0. .AND. PIDEP == 0.)
1877  ENDIF ! EndIF (TC >= T_ICE .AND. (QW > EPSQ .OR. WV > QSWgrd))
1878 !
1879 !--- Limit freezing of accreted rime to prevent temperature oscillations,
1880 ! a crude Schumann-Ludlam limit (p. 209 of Young, 1993).
1881 !
1882  tcc = tc + xlv1*pcond + xls1*pidep + xlf1*piacwi
1883  IF (tcc > 0.) THEN
1884  piacwi = 0.
1885  tcc = tc + xlv1*pcond + xls1*pidep
1886  ENDIF
1887 !
1888  IF (tc > 0. .AND. tcc > 0. .AND. ice_logical) THEN
1889 !
1890 !--- Calculate melting and evaporation/condensation
1891 ! * Units: SFACTOR - s**.5/m ; ABI - m**2/s ; NLICE - m**-3 ;
1892 ! VENTIL - m**-2 ; VENTI1 - m ;
1893 ! VENTI2 - m**2/s**.5 ; CIEVP - /s
1894 !
1895 ! SFACTOR = VEL_INC**.5*(RHO/(DIFFUS*DIFFUS*DYNVIS))**C2
1896  sfactor = sqrt(vel_inc)*schmit_fac
1897  ventil = nlice*(venti1(indexs)+sfactor*venti2(indexs))
1898  aievp = ventil*diffus*dtph
1899  IF (aievp < xratio) THEN
1900  dievp = aievp
1901  ELSE
1902  dievp = 1. - exp(-aievp)
1903  ENDIF
1904 ! QSW0 = EPS*ESW0/(PP-ESW0)
1905 ! QSW0 = EPS*ESW0/(PP+epsm1*ESW0)
1906 !! dum = min(PP, ESW0)
1907 !! QSW0 = EPS*dum/(PP+epsm1*dum)
1908 ! DWV0 = MIN(WV,QSW)-QSW0
1909  dwv0 = min(wv,qsw_l)-qsw0_l
1910  dum = qw + pcond
1911  IF (wv < qsw_l .AND. dum <= epsq) THEN
1912  !
1913  !--- Evaporation from melting snow (sink of snow) or shedding
1914  ! of water condensed onto melting snow (source of rain)
1915  !
1916  dum = dwv0*dievp
1917  pievp = max( min(0., dum), piloss)
1918  picnd = max(0., dum)
1919  ENDIF ! End IF (WV < QSW_l .AND. DUM <= EPSQ)
1920  pimlt = therm_cond*tcc*ventil*rrho*dtph/xlf
1921  !
1922  !--- Limit melting to prevent temperature oscillations across 0C
1923  !
1924  dum1 = max( 0., (tcc+xlv1*pievp)/xlf1 )
1925  pimlt = min(pimlt, dum1)
1926  !
1927  !--- Limit loss of snow by melting (>0) and evaporation
1928  !
1929  dum = pievp - pimlt
1930  IF (dum < piloss) THEN
1931  dum1 = piloss/dum
1932  pimlt = pimlt*dum1
1933  pievp = pievp*dum1
1934  ENDIF ! End IF (DUM > QTICE)
1935  ENDIF ! End IF (TC > 0. .AND. TCC > 0. .AND. ICE_logical)
1936 !
1937 !--- IMPORTANT: Estimate time-averaged properties.
1938 !
1939 ! * TOT_RAIN - total mass of rain before microphysics, which is the sum of
1940 ! the total mass of rain in the current layer and the input
1941 ! flux of rain from above
1942 ! * VRAIN1 - fall speed of rain into grid from above (with air resistance correction)
1943 ! * QTRAIN - time-averaged mixing ratio of rain (kg/kg)
1944 ! * PRLOSS - greatest loss (<0) of rain, removing all rain falling from
1945 ! above and the rain within the layer
1946 ! * RQR - rain content (kg/m**3)
1947 ! * INDEXR - mean size of rain drops to the nearest 1 micron in size
1948 ! * N0r - intercept of rain size distribution (typically 10**6 m**-4)
1949 !
1950  tot_rain = 0.
1951  vrain1 = 0.
1952  qtrain = 0.
1953  prloss = 0.
1954  rqr = 0.
1955  n0r = 0.
1956  indexr = mdrmin
1957  indexr1 = indexr ! For debugging only
1958  IF (rain_logical) THEN
1959  IF (arain <= 0.) THEN
1960  indexr = mdrmin
1961  vrain1 = 0.
1962  ELSE
1963  !
1964  !--- INDEXR (related to mean diameter) & N0r could be modified
1965  ! by land/sea properties, presence of convection, etc.
1966  !
1967  !--- Rain rate normalized to a density of 1.194 kg/m**3
1968  !
1969  rr = arain / (dtph*gammar)
1970  !
1971  IF (rr <= rr_drmin) THEN
1972  !
1973  !--- Assume fixed mean diameter of rain (0.2 mm) for low rain rates,
1974  ! instead vary N0r with rain rate
1975  !
1976  indexr = mdrmin
1977  ELSE IF (rr <= rr_dr1) THEN
1978  !
1979  !--- Best fit to mass-weighted fall speeds (V) from rain lookup tables
1980  ! for mean diameters (Dr) between 0.05 and 0.10 mm:
1981  ! V(Dr)=5.6023e4*Dr**1.136, V in m/s and Dr in m
1982  ! RR = PI*1000.*N0r0*5.6023e4*Dr**(4+1.136) = 1.408e15*Dr**5.136,
1983  ! RR in kg/(m**2*s)
1984  ! Dr (m) = 1.123e-3*RR**.1947 -> Dr (microns) = 1.123e3*RR**.1947
1985  !
1986  indexr = int( 1.123e3*rr**.1947 + .5 )
1987  indexr = max( mdrmin, min(indexr, mdr1) )
1988 
1989  ELSE IF (rr <= rr_dr2) THEN
1990  !
1991  !--- Best fit to mass-weighted fall speeds (V) from rain lookup tables
1992  ! for mean diameters (Dr) between 0.10 and 0.20 mm:
1993  ! V(Dr)=1.0867e4*Dr**.958, V in m/s and Dr in m
1994  ! RR = PI*1000.*N0r0*1.0867e4*Dr**(4+.958) = 2.731e14*Dr**4.958,
1995  ! RR in kg/(m**2*s)
1996  ! Dr (m) = 1.225e-3*RR**.2017 -> Dr (microns) = 1.225e3*RR**.2017
1997  !
1998  indexr = int( 1.225e3*rr**.2017 + .5 )
1999  indexr = max( mdr1, min(indexr, mdr2) )
2000 
2001  ELSE IF (rr <= rr_dr3) THEN
2002  !
2003  !--- Best fit to mass-weighted fall speeds (V) from rain lookup tables
2004  ! for mean diameters (Dr) between 0.20 and 0.32 mm:
2005  ! V(Dr)=2831.*Dr**.80, V in m/s and Dr in m
2006  ! RR = PI*1000.*N0r0*2831.*Dr**(4+.80) = 7.115e13*Dr**4.80,
2007  ! RR in kg/(m**2*s)
2008  ! Dr (m) = 1.3006e-3*RR**.2083 -> Dr (microns) = 1.3006e3*RR**.2083
2009  !
2010  indexr = int( 1.3006e3*rr**.2083 + .5 )
2011  indexr = max( mdr2, min(indexr, mdr3) )
2012 
2013  ELSE IF (rr <= rr_drmax) THEN
2014  !
2015  !--- Best fit to mass-weighted fall speeds (V) from rain lookup tables
2016  ! for mean diameters (Dr) between 0.32 and 0.45 mm:
2017  ! V(Dr)=944.8*Dr**.6636, V in m/s and Dr in m
2018  ! RR = PI*1000.*N0r0*944.8*Dr**(4+.6636) = 2.3745e13*Dr**4.6636,
2019  ! RR in kg/(m**2*s)
2020  ! Dr (m) = 1.355e-3*RR**.2144 -> Dr (microns) = 1.355e3*RR**.2144
2021  !
2022  indexr = int( 1.355e3*rr**.2144 + .5 )
2023  indexr = max( mdr3, min(indexr, mdrmax) )
2024  ELSE
2025  !
2026  !--- Assume fixed mean diameter of rain (0.45 mm) for high rain rates,
2027  ! instead vary N0r with rain rate
2028  !
2029  indexr = mdrmax
2030  ENDIF ! End IF (RR <= RR_DRmin) etc.
2031 !
2032  vrain1 = gammar*vrain(indexr)
2033  ENDIF ! End IF (ARAIN <= 0.)
2034 !
2035  indexr1 = indexr ! For debugging only
2036  tot_rain = thick*qr+blend*arain
2037  qtrain = tot_rain/(thick+bldtrh*vrain1)
2038  prloss = -tot_rain/thick
2039  rqr = rho*qtrain
2040  !
2041  !--- RQR - time-averaged rain content (kg/m**3)
2042  !
2043  IF (rqr <= rqr_drmin) THEN
2044  n0r = max(n0rmin, cn0r_dmrmin*rqr)
2045  indexr = mdrmin
2046  ELSE IF (rqr >= rqr_drmax) THEN
2047  n0r = cn0r_dmrmax*rqr
2048  indexr = mdrmax
2049  ELSE
2050  n0r = n0r0
2051 ! INDEXR = MAX( XMRmin, MIN(CN0r0*RQR**.25, XMRmax) )
2052  item = cn0r0*sqrt(sqrt(rqr)) ! Moorthi 07/31/08
2053  indexr = max( mdrmin, min(item, mdrmax) ) ! Moorthi 07/31/08
2054  ENDIF
2055  !
2056  IF (tc < t_ice) THEN
2057  piacr = -prloss
2058  ELSE
2059  dwvr = wv - pcond - qsw_l
2060  dum = qw + pcond
2061  IF (dwvr < 0. .AND. dum <= epsq) THEN
2062 !
2063 !--- Rain evaporation
2064 !
2065 ! * RFACTOR - [GAMMAR**.5]*[Schmidt**(1./3.)]*[(RHO/DYNVIS)**.5],
2066 ! where Schmidt (Schmidt Number) =DYNVIS/(RHO*DIFFUS)
2067 !
2068 ! * Units: RFACTOR - s**.5/m ; ABW - m**2/s ; VENTR - m**-2 ;
2069 ! N0r - m**-4 ; VENTR1 - m**2 ; VENTR2 - m**3/s**.5 ;
2070 ! CREVP - unitless
2071 !
2072 ! RFACTOR = GAMMAR**.5*(RHO/(DIFFUS*DIFFUS*DYNVIS))**C2
2073  rfactor = sqrt(gammar)*schmit_fac
2074  abw = 1./(rho*xlv2/therm_cond+1./diffus)
2075 !
2076 !--- Note that VENTR1, VENTR2 lookup tables do not include the
2077 ! 1/Davg multiplier as in the ice tables
2078 !
2079  ventr = n0r*(ventr1(indexr)+rfactor*ventr2(indexr))
2080  crevp = abw*ventr*dtph
2081  IF (crevp < xratio) THEN
2082  dum = dwvr*crevp
2083  ELSE
2084  dum = dwvr*(1.-exp(-crevp*denomw))/denomw
2085  ENDIF
2086  prevp = max(dum, prloss)
2087  ELSE IF (qw > epsq) THEN
2088  fwr = cracw*gammar*n0r*accrr(indexr)
2089 !Moor PRACW = MIN(1.,FWR)*QW ! 20050422
2090  pracw = min(0.1,fwr)*qw
2091  ENDIF ! End IF (DWVr < 0. .AND. DUM <= EPSQ)
2092 !
2093  IF (tc < 0. .AND. tcc < 0.) THEN
2094 !
2095 !--- Biggs (1953) heteorogeneous freezing (e.g., Lin et al., 1983)
2096 ! - Rescaled mean drop diameter from microns (INDEXR) to mm (DUM) to prevent underflow
2097 !
2098  dum = .001*float(indexr)
2099  dum1 = dum * dum
2100  dum = (exp(abfr*tc)-1.)*dum1*dum1*dum1*dum
2101  piacr = min(cbfr*n0r*rrho*dum, qtrain)
2102  IF (qlice > epsq) THEN
2103  !
2104  !--- Freezing of rain by collisions w/ large ice
2105  !
2106  dum = gammar*vrain(indexr)
2107  dum1 = dum-vsnow
2108  !
2109  !--- DUM2 - Difference in spectral fall speeds of rain and
2110  ! large ice, parameterized following eq. (48) on p. 112 of
2111  ! Murakami (J. Meteor. Soc. Japan, 1990)
2112  !
2113  dum2 = (dum1*dum1+.04*dum*vsnow)**.5
2114  dum1 = 5.e-12*indexr*indexr+2.e-12*indexr*indexs &
2115  & +.5e-12*indexs*indexs
2116  fir = min(1., ciacr*nlice*dum1*dum2)
2117  !
2118  !--- Future? Should COLLECTION BY SMALL ICE SHOULD BE INCLUDED???
2119  !
2120  piacr = min(piacr+fir*qtrain, qtrain)
2121  ENDIF ! End IF (QLICE > EPSQ)
2122  dum = prevp - piacr
2123  If (dum < prloss) THEN
2124  dum1 = prloss/dum
2125  prevp = dum1*prevp
2126  piacr = dum1*piacr
2127  ENDIF ! End If (DUM < PRLOSS)
2128  ENDIF ! End IF (TC < 0. .AND. TCC < 0.)
2129  ENDIF ! End IF (TC < T_ICE)
2130  ENDIF ! End IF (RAIN_logical)
2131 !
2132 !----------------------------------------------------------------------
2133 !---------------------- Main Budget Equations -------------------------
2134 !----------------------------------------------------------------------
2135 !
2136 !
2137 !-----------------------------------------------------------------------
2138 !--- Update fields, determine characteristics for next lower layer ----
2139 !-----------------------------------------------------------------------
2140 !
2141 !--- Carefully limit sinks of cloud water
2142 !
2143  dum1 = piacw + praut + pracw - min(0.,pcond)
2144  IF (dum1 > qw) THEN
2145  dum = qw/dum1
2146  piacw = dum*piacw
2147  piacwi = dum*piacwi
2148  praut = dum*praut
2149  pracw = dum*pracw
2150  IF (pcond < 0.) pcond=dum*pcond
2151  ENDIF
2152  piacwr = piacw - piacwi ! TC >= 0C
2153 !
2154 !--- QWnew - updated cloud water mixing ratio
2155 !
2156  delw = pcond - piacw - praut - pracw
2157  qwnew = qw+delw
2158  IF (qwnew <= epsq) qwnew = 0.
2159  IF (qw > 0. .AND. qwnew /= 0.) THEN
2160  dum = qwnew/qw
2161  IF (dum < toler) qwnew = 0.
2162  ENDIF
2163 !
2164 !--- Update temperature and water vapor mixing ratios
2165 !
2166  delt = xlv1 * (pcond+pievp+picnd+prevp) &
2167  & + xls1 * pidep + xlf1*(piacwi+piacr-pimlt)
2168  tnew = tk + delt
2169 !
2170  delv = -pcond - pidep - pievp - picnd - prevp
2171  wvnew = wv + delv
2172 !
2173 !--- Update ice mixing ratios
2174 !
2175 !---
2176 ! * TOT_ICEnew - total mass (small & large) ice after microphysics,
2177 ! which is the sum of the total mass of large ice in the
2178 ! current layer and the flux of ice out of the grid box below
2179 ! * RimeF - Rime Factor, which is the mass ratio of total (unrimed &
2180 ! rimed) ice mass to the unrimed ice mass (>=1)
2181 ! * QInew - updated mixing ratio of total (large & small) ice in layer
2182 ! -> TOT_ICEnew=QInew*THICK+BLDTRH*QLICEnew*VSNOW
2183 ! -> But QLICEnew=QInew*FLIMASS, so
2184 ! -> TOT_ICEnew=QInew*(THICK+BLDTRH*FLIMASS*VSNOW)
2185 ! * ASNOWnew - updated accumulation of snow at bottom of grid cell
2186 !---
2187 !
2188  deli = 0.
2189  rimef = 1.
2190  IF (ice_logical) THEN
2191  deli = pidep + pievp + piacwi + piacr - pimlt
2192  tot_icenew = tot_ice + thick*deli
2193  IF (tot_ice > 0. .AND. tot_icenew /= 0.) THEN
2194  dum = tot_icenew/tot_ice
2195  IF (dum < toler) tot_icenew = 0.
2196  ENDIF
2197  IF (tot_icenew <= climit) THEN
2198  tot_icenew = 0.
2199  rimef = 1.
2200  qinew = 0.
2201  asnownew = 0.
2202  ELSE
2203  !
2204  !--- Update rime factor if appropriate
2205  !
2206  dum = piacwi + piacr
2207  IF (dum <= epsq .AND. pidep <= epsq) THEN
2208  rimef = rimef1
2209  ELSE
2210  !
2211  !--- Rime Factor, RimeF = (Total ice mass)/(Total unrimed ice mass)
2212  ! DUM1 - Total ice mass, rimed & unrimed
2213  ! DUM2 - Estimated mass of *unrimed* ice
2214  !
2215  dum1 = tot_ice + thick*(pidep+dum)
2216  dum2 = tot_ice/rimef1 + thick*pidep
2217  IF (dum2 <= 0.) THEN
2218  rimef = rfmax
2219  ELSE
2220  rimef = min(rfmax, max(1., dum1/dum2) )
2221  ENDIF
2222  ENDIF ! End IF (DUM <= EPSQ .AND. PIDEP <= EPSQ)
2223  qinew = tot_icenew/(thick+bldtrh*flimass*vsnow)
2224  IF (qinew <= epsq) qinew = 0.
2225  IF (qi > 0. .AND. qinew /= 0.) THEN
2226  dum = qinew/qi
2227  IF (dum < toler) qinew = 0.
2228  ENDIF
2229  asnownew = bldtrh*flimass*vsnow*qinew
2230  IF (asnow > 0. .AND. asnownew /= 0.) THEN
2231  dum = asnownew/asnow
2232  IF (dum < toler) asnownew = 0.
2233  ENDIF
2234  ENDIF ! End IF (TOT_ICEnew <= CLIMIT)
2235  ENDIF ! End IF (ICE_logical)
2236 !
2237 !--- Update rain mixing ratios
2238 !
2239 !---
2240 ! * TOT_RAINnew - total mass of rain after microphysics
2241 ! current layer and the input flux of ice from above
2242 ! * VRAIN2 - time-averaged fall speed of rain in grid and below
2243 ! (with air resistance correction)
2244 ! * QRnew - updated rain mixing ratio in layer
2245 ! -> TOT_RAINnew=QRnew*(THICK+BLDTRH*VRAIN2)
2246 ! * ARAINnew - updated accumulation of rain at bottom of grid cell
2247 !---
2248 !
2249  delr = praut + pracw + piacwr - piacr + pimlt &
2250  & + prevp + picnd
2251  tot_rainnew = tot_rain+thick*delr
2252  IF (tot_rain > 0. .AND. tot_rainnew /= 0.) THEN
2253  dum = tot_rainnew/tot_rain
2254  IF (dum < toler) tot_rainnew = 0.
2255  ENDIF
2256  IF (tot_rainnew <= climit) THEN
2257  tot_rainnew = 0.
2258  vrain2 = 0.
2259  qrnew = 0.
2260  arainnew = 0.
2261  ELSE
2262  !
2263  !--- 1st guess time-averaged rain rate at bottom of grid box
2264  !
2265  rr = tot_rainnew/(dtph*gammar)
2266  !
2267  !--- Use same algorithm as above for calculating mean drop diameter
2268  ! (IDR, in microns), which is used to estimate the time-averaged
2269  ! fall speed of rain drops at the bottom of the grid layer. This
2270  ! isn't perfect, but the alternative is solving a transcendental
2271  ! equation that is numerically inefficient and nasty to program
2272  ! (coded in earlier versions of GSMCOLUMN prior to 8-22-01).
2273  !
2274  IF (rr <= rr_drmin) THEN
2275  idr = mdrmin
2276  ELSE IF (rr <= rr_dr1) THEN
2277  idr = int( 1.123e3*rr**.1947 + .5 )
2278  idr = max( mdrmin, min(idr, mdr1) )
2279  ELSE IF (rr <= rr_dr2) THEN
2280  idr = int( 1.225e3*rr**.2017 + .5 )
2281  idr = max( mdr1, min(idr, mdr2) )
2282  ELSE IF (rr <= rr_dr3) THEN
2283  idr = int( 1.3006e3*rr**.2083 + .5 )
2284  idr = max( mdr2, min(idr, mdr3) )
2285  ELSE IF (rr <= rr_drmax) THEN
2286  idr = int( 1.355e3*rr**.2144 + .5 )
2287  idr = max( mdr3, min(idr, mdrmax) )
2288  ELSE
2289  idr = mdrmax
2290  ENDIF ! End IF (RR <= RR_DRmin)
2291  vrain2 = gammar*vrain(idr)
2292  qrnew = tot_rainnew / (thick+bldtrh*vrain2)
2293  IF (qrnew <= epsq) qrnew = 0.
2294  IF (qr > 0. .AND. qrnew /= 0.) THEN
2295  dum = qrnew / qr
2296  IF (dum < toler) qrnew = 0.
2297  ENDIF
2298  arainnew = bldtrh*vrain2*qrnew
2299  IF (arain > 0. .AND. arainnew /= 0.) THEN
2300  dum = arainnew/arain
2301  IF (dum < toler) arainnew = 0.
2302  ENDIF
2303  ENDIF ! End IF (TOT_RAINnew < CLIMIT)
2304 !
2305  wcnew = qwnew + qrnew + qinew
2306 !
2307 !----------------------------------------------------------------------
2308 !-------------- Begin debugging & verification ------------------------
2309 !----------------------------------------------------------------------
2310 !
2311 !--- QT, QTnew - total water (vapor & condensate) before & after microphysics, resp.
2312 !
2313 ! QT=THICK*(QV+WC_col(l))+ARAIN+ASNOW
2314 ! QTnew = THICK*(WVnew/(1.+WVnew)+WCnew/(1.+wcnew))
2315 ! & + ARAINnew + ASNOWnew
2316 
2317  qt = thick*(wv+wc) + arain + asnow
2318  qtnew = thick*(wvnew+wcnew) + arainnew + asnownew
2319  budget = qt-qtnew
2320 !
2321 !--- Additional check on budget preservation, accounting for truncation effects
2322 !
2323  dbg_logical=.false.
2324 ! DUM=ABS(BUDGET)
2325 ! IF (DUM > TOLER) THEN
2326 ! DUM=DUM/MIN(QT, QTnew)
2327 ! IF (DUM > TOLER) DBG_logical=.TRUE.
2328 ! ENDIF
2329 !
2330 ! DUM=(RHgrd+.001)*QSInew
2331 ! IF ( (QWnew > EPSQ .OR. QRnew > EPSQ .OR. WVnew > DUM)
2332 ! & .AND. TC < T_ICE ) DBG_logical=.TRUE.
2333 !
2334 ! IF (TC > 5. .AND. QInewr > EPSQ) DBG_logical=.TRUE.
2335 !
2336  IF ((wvnew < epsq .OR. dbg_logical) .AND. print_diag) THEN
2337 !
2338  WRITE(6,"(/2(a,i4),2(a,i2))") '{} i=',i_index,' j=', &
2339  & j_index, ' L=',l,' LSFC=',lsfc
2340 !
2341  esw = min(pp, fpvsl(tnew))
2342 ! QSWnew = EPS*ESW/(PP-ESW)
2343  qswnew = eps*esw/(pp+epsm1*esw)
2344  IF (tc < 0. .OR. tnew < 0.) THEN
2345  esi = min(pp, fpvsi(tnew))
2346 ! QSInew = EPS*ESI/(PP-ESI)
2347  qsinew = eps*esi/(pp+epsm1*esi)
2348  ELSE
2349  qsi = qsw
2350  qsinew = qswnew
2351  ENDIF
2352  wsnew = qsinew
2353  WRITE(6,"(4(a12,g11.4,1x))") &
2354  & '{} TCold=',tc,'TCnew=',tnew-t0c,'P=',.01*pp,'RHO=',rho, &
2355  & '{} THICK=',thick,'RHold=',wv/ws,'RHnew=',wvnew/wsnew, &
2356  & 'RHgrd=',rhgrd, &
2357  & '{} RHWold=',wv/qsw,'RHWnew=',wvnew/qswnew,'RHIold=',wv/qsi, &
2358  & 'RHInew=',wvnew/qsinew, &
2359  & '{} QSWold=',qsw,'QSWnew=',qswnew,'QSIold=',qsi,'QSInew=',qsinew,&
2360  & '{} WSold=',ws,'WSnew=',wsnew,'WVold=',wv,'WVnew=',wvnew, &
2361  & '{} WCold=',wc,'WCnew=',wcnew,'QWold=',qw,'QWnew=',qwnew, &
2362  & '{} QIold=',qi,'QInew=',qinew,'QRold=',qr,'QRnew=',qrnew, &
2363  & '{} ARAINold=',arain,'ARAINnew=',arainnew,'ASNOWold=',asnow, &
2364  & 'ASNOWnew=',asnownew, &
2365  & '{} TOT_RAIN=',tot_rain,'TOT_RAINnew=',tot_rainnew, &
2366  & 'TOT_ICE=',tot_ice,'TOT_ICEnew=',tot_icenew, &
2367  & '{} BUDGET=',budget,'QTold=',qt,'QTnew=',qtnew
2368 !
2369  WRITE(6,"(4(a12,g11.4,1x))") &
2370  & '{} DELT=',delt,'DELV=',delv,'DELW=',delw,'DELI=',deli, &
2371  & '{} DELR=',delr,'PCOND=',pcond,'PIDEP=',pidep,'PIEVP=',pievp, &
2372  & '{} PICND=',picnd,'PREVP=',prevp,'PRAUT=',praut,'PRACW=',pracw, &
2373  & '{} PIACW=',piacw,'PIACWI=',piacwi,'PIACWR=',piacwr,'PIMLT=', &
2374  & pimlt, &
2375  & '{} PIACR=',piacr
2376 !
2377  IF (ice_logical) WRITE(6,"(4(a12,g11.4,1x))") &
2378  & '{} RimeF1=',rimef1,'GAMMAS=',gammas,'VrimeF=',vrimef, &
2379  & 'VSNOW=',vsnow, &
2380  & '{} INDEXS=',float(indexs),'FLARGE=',flarge,'FSMALL=',fsmall, &
2381  & 'FLIMASS=',flimass, &
2382  & '{} XSIMASS=',xsimass,'XLIMASS=',xlimass,'QLICE=',qlice, &
2383  & 'QTICE=',qtice, &
2384  & '{} NLICE=',nlice,'NSmICE=',nsmice,'PILOSS=',piloss, &
2385  & 'EMAIRI=',emairi, &
2386  & '{} RimeF=',rimef
2387 !
2388  IF (tot_rain > 0. .OR. tot_rainnew > 0.) &
2389  & WRITE(6,"(4(a12,g11.4,1x))") &
2390  & '{} INDEXR1=',float(indexr1),'INDEXR=',float(indexr), &
2391  & 'GAMMAR=',gammar,'N0r=',n0r, &
2392  & '{} VRAIN1=',vrain1,'VRAIN2=',vrain2,'QTRAIN=',qtrain,'RQR=',rqr,&
2393  & '{} PRLOSS=',prloss,'VOLR1=',thick+bldtrh*vrain1, &
2394  & 'VOLR2=',thick+bldtrh*vrain2
2395 !
2396  IF (praut > 0.) WRITE(6,"(a12,g11.4,1x)") '{} QW0=',qw0
2397 !
2398  IF (pracw > 0.) WRITE(6,"(a12,g11.4,1x)") '{} FWR=',fwr
2399 !
2400  IF (piacr > 0.) WRITE(6,"(a12,g11.4,1x)") '{} FIR=',fir
2401 !
2402  dum = pimlt + picnd - prevp - pievp
2403  IF (dum > 0. .or. dwvi /= 0.) &
2404  & WRITE(6,"(4(a12,g11.4,1x))") &
2405  & '{} TFACTOR=',tfactor,'DYNVIS=',dynvis, &
2406  & 'THERM_CON=',therm_cond,'DIFFUS=',diffus
2407 !
2408  IF (prevp < 0.) WRITE(6,"(4(a12,g11.4,1x))") &
2409  & '{} RFACTOR=',rfactor,'ABW=',abw,'VENTR=',ventr,'CREVP=',crevp, &
2410  & '{} DWVr=',dwvr,'DENOMW=',denomw
2411 !
2412  IF (pidep /= 0. .AND. dwvi /= 0.) &
2413  & WRITE(6,"(4(a12,g11.4,1x))") &
2414  & '{} DWVi=',dwvi,'DENOMI=',denomi,'PIDEP_max=',pidep_max, &
2415  & 'SFACTOR=',sfactor, &
2416  & '{} ABI=',abi,'VENTIL=',ventil,'VENTIL1=',venti1(indexs), &
2417  & 'VENTIL2=',sfactor*venti2(indexs), &
2418  & '{} VENTIS=',ventis,'DIDEP=',didep
2419 !
2420  IF (pidep > 0. .AND. pcond /= 0.) &
2421  & WRITE(6,"(4(a12,g11.4,1x))") &
2422  & '{} DENOMW=',denomw,'DENOMWI=',denomwi,'DENOMF=',denomf, &
2423  & 'DUM2=',pcond-piacw
2424 !
2425  IF (fws > 0.) WRITE(6,"(4(a12,g11.4,1x))") '{} FWS=',fws
2426 !
2427  dum = pimlt + picnd - pievp
2428  IF (dum > 0.) WRITE(6,"(4(a12,g11.4,1x))") &
2429  & '{} SFACTOR=',sfactor,'VENTIL=',ventil,'VENTIL1=',venti1(indexs),&
2430  & 'VENTIL2=',sfactor*venti2(indexs), &
2431  & '{} AIEVP=',aievp,'DIEVP=',dievp,'QSW0=',qsw0,'DWV0=',dwv0
2432  !
2433  ENDIF
2434 !
2435 !----------------------------------------------------------------------
2436 !-------------- Water budget statistics & maximum values --------------
2437 !----------------------------------------------------------------------
2438 !
2439  IF (print_diag) THEN
2440  itdx = max( itlo, min( int(tnew-t0c), ithi ) )
2441  IF (qinew > epsq) nstats(itdx,1) = nstats(itdx,1)+1
2442  IF (qinew > epsq .AND. qrnew+qwnew > epsq) &
2443  & nstats(itdx,2) = nstats(itdx,2)+1
2444  IF (qwnew > epsq) nstats(itdx,3) = nstats(itdx,3)+1
2445  IF (qrnew > epsq) nstats(itdx,4) = nstats(itdx,4)+1
2446  !
2447  qmax(itdx,1) = max(qmax(itdx,1), qinew)
2448  qmax(itdx,2) = max(qmax(itdx,2), qwnew)
2449  qmax(itdx,3) = max(qmax(itdx,3), qrnew)
2450  qmax(itdx,4) = max(qmax(itdx,4), asnownew)
2451  qmax(itdx,5) = max(qmax(itdx,5), arainnew)
2452  qtot(itdx,1) = qtot(itdx,1)+qinew*thick
2453  qtot(itdx,2) = qtot(itdx,2)+qwnew*thick
2454  qtot(itdx,3) = qtot(itdx,3)+qrnew*thick
2455  !
2456  qtot(itdx,4) = qtot(itdx,4)+pcond*thick
2457  qtot(itdx,5) = qtot(itdx,5)+picnd*thick
2458  qtot(itdx,6) = qtot(itdx,6)+pievp*thick
2459  qtot(itdx,7) = qtot(itdx,7)+pidep*thick
2460  qtot(itdx,8) = qtot(itdx,8)+prevp*thick
2461  qtot(itdx,9) = qtot(itdx,9)+praut*thick
2462  qtot(itdx,10) = qtot(itdx,10)+pracw*thick
2463  qtot(itdx,11) = qtot(itdx,11)+pimlt*thick
2464  qtot(itdx,12) = qtot(itdx,12)+piacw*thick
2465  qtot(itdx,13) = qtot(itdx,13)+piacwi*thick
2466  qtot(itdx,14) = qtot(itdx,14)+piacwr*thick
2467  qtot(itdx,15) = qtot(itdx,15)+piacr*thick
2468  !
2469  qtot(itdx,16) = qtot(itdx,16)+(wvnew-wv)*thick
2470  qtot(itdx,17) = qtot(itdx,17)+(qwnew-qw)*thick
2471  qtot(itdx,18) = qtot(itdx,18)+(qinew-qi)*thick
2472  qtot(itdx,19) = qtot(itdx,19)+(qrnew-qr)*thick
2473  qtot(itdx,20) = qtot(itdx,20)+(arainnew-arain)
2474  qtot(itdx,21) = qtot(itdx,21)+(asnownew-asnow)
2475  IF (qinew > 0.) &
2476  & qtot(itdx,22) = qtot(itdx,22)+qinew*thick/rimef
2477  !
2478  ENDIF
2479 !
2480 !----------------------------------------------------------------------
2481 !------------------------- Update arrays ------------------------------
2482 !----------------------------------------------------------------------
2483 !
2484  t_col(l) = tnew ! temperature
2485 !
2486 ! QV_col(L) = max(EPSQ, WVnew/(1.+WVnew)) ! specific humidity
2487  qv_col(l) = max(0.0, wvnew ) ! specific humidity
2488  wc_col(l) = max(0.0, wcnew) ! total condensate mixing ratio
2489  qi_col(l) = max(0.0, qinew) ! ice mixing ratio
2490  qr_col(l) = max(0.0, qrnew) ! rain mixing ratio
2491  qw_col(l) = max(0.0, qwnew) ! cloud water mixing ratio
2492  rimef_col(l) = rimef ! rime factor
2493  asnow = asnownew ! accumulated snow
2494  arain = arainnew ! accumulated rain
2495 !
2496 !#######################################################################
2497 !
2498  ENDIF ! End of IF (.NOT. CLEAR) THEN
2499  ENDIF ! End of IF (QV_col(L) <= EPSQ .AND. WC_col(L) <= EPSQ) THEN
2500 !
2501  ENDDO ! ##### End "L" loop through model levels #####
2502 !
2503  araing = araing + arain
2504  asnowg = asnowg + asnow
2505  enddo ! do for ntimes=1,mic_step
2506 !
2507 !#######################################################################
2508 !
2509 !-----------------------------------------------------------------------
2510 !--------------------------- Return to GSMDRIVE -----------------------
2511 !-----------------------------------------------------------------------
2512 !
2513  CONTAINS
2514 ! END SUBROUTINE GSMCOLUMN
2515 !
2516 !#######################################################################
2517 !--------- Produces accurate calculation of cloud condensation ---------
2518 !#######################################################################
2519 !
2520  REAL FUNCTION condense (PP, QW, RHgrd, TK, WV)
2522  implicit none
2523 !
2524 !---------------------------------------------------------------------------------
2525 !------ The Asai (1965) algorithm takes into consideration the release of ------
2526 !------ latent heat in increasing the temperature & in increasing the ------
2527 !------ saturation mixing ratio (following the Clausius-Clapeyron eqn.). ------
2528 !---------------------------------------------------------------------------------
2529 !
2530  real pp, qw, rhgrd, tk, wv
2531  INTEGER, PARAMETER :: HIGH_PRES=kind_phys
2532 ! INTEGER, PARAMETER :: HIGH_PRES=Selected_Real_Kind(15)
2533  REAL (KIND=HIGH_PRES), PARAMETER :: &
2534  & RHLIMIT=.001, RHLIMIT1=-RHLIMIT
2535  REAL, PARAMETER :: RCP=1./cp, rcprv=rcp/rv
2536  REAL (KIND=HIGH_PRES) :: COND, SSAT, WCdum, tsq
2537  real wvdum, tdum, xlv, xlv1, xlv2, ws, dwv, esw, rfac
2538 !
2539 !-----------------------------------------------------------------------
2540 !
2541 !--- LV (T) is from Bolton (JAS, 1980)
2542 !
2543 ! XLV=3.148E6-2370.*TK
2544 ! XLV1=XLV*RCP
2545 ! XLV2=XLV*XLV*RCPRV
2546 !
2547  tdum = tk
2548  wvdum = wv
2549  wcdum = qw
2550  esw = min(pp, fpvsl(tdum)) ! Saturation vapor press w/r/t water
2551 ! WS = RHgrd*EPS*ESW/(PP-ESW) ! Saturation mixing ratio w/r/t water
2552  ws = rhgrd*eps*esw/(pp+epsm1*esw) ! Saturation specific hum w/r/t water
2553  dwv = wvdum - ws ! Deficit grid-scale specific humidity
2554  ssat = dwv / ws ! Supersaturation ratio
2555  condense = 0.
2556  rfac = 0.5 ! converges faster with 0.5
2557  DO WHILE ((ssat < rhlimit1 .AND. wcdum > epsq) &
2558  & .OR. ssat > rhlimit)
2559 !
2560  xlv = 3.148e6-2370.*tdum
2561  xlv1 = xlv*rcp
2562  xlv2 = xlv*xlv*rcprv
2563 !
2564 ! COND = DWV/(1.+XLV2*WS/(Tdum*Tdum)) ! Asai (1965, J. Japan)
2565  tsq = tdum*tdum
2566  cond = rfac*dwv*tsq/(tsq+xlv2*ws) ! Asai (1965, J. Japan)
2567 ! COND = DWV*tsq/(tsq+XLV2*WS) ! Asai (1965, J. Japan)
2568  cond = max(cond, -wcdum) ! Limit cloud water evaporation
2569  tdum = tdum+xlv1*cond ! Updated temperature
2570  wvdum = wvdum-cond ! Updated water vapor mixing ratio
2571  wcdum = wcdum+cond ! Updated cloud water mixing ratio
2572  condense = condense + cond ! Total cloud water condensation
2573  esw = min(pp, fpvsl(tdum)) ! Updated saturation vapor press w/r/t water
2574 ! WS = RHgrd*EPS*ESW/(PP-ESW) ! Updated saturation mixing ratio w/r/t water
2575  ws = rhgrd*eps*esw/(pp+epsm1*esw) ! Updated saturation mixing ratio w/r/t water
2576  dwv = wvdum-ws ! Deficit grid-scale water vapor mixing ratio
2577  ssat = dwv / ws ! Grid-scale supersaturation ratio
2578  rfac = 1.0
2579  ENDDO
2580 
2581  END FUNCTION condense
2582 !
2583 !#######################################################################
2584 !---------------- Calculate ice deposition at T<T_ICE ------------------
2585 !#######################################################################
2586 !
2587  REAL FUNCTION deposit (PP, RHgrd, Tdum, WVdum)
2589  implicit none
2590 !
2591 !--- Also uses the Asai (1965) algorithm, but uses a different target
2592 ! vapor pressure for the adjustment
2593 !
2594  REAL PP, RHgrd, Tdum, WVdum
2595  INTEGER, PARAMETER :: HIGH_PRES=kind_phys
2596 ! INTEGER, PARAMETER :: HIGH_PRES=Selected_Real_Kind(15)
2597  REAL (KIND=HIGH_PRES), PARAMETER :: RHLIMIT=.001, &
2598  & RHLIMIT1=-RHLIMIT
2599  real, PARAMETER :: rcp=1./cp, rcprv=rcp/rv, xls=hvap+hfus &
2600  &, xls1=xls*rcp, xls2=xls*xls*rcprv
2601  REAL (KIND=HIGH_PRES) :: DEP, SSAT
2602  real esi, ws, dwv
2603 !
2604 !-----------------------------------------------------------------------
2605 !
2606  esi=min(pp, fpvsi(tdum)) ! Saturation vapor press w/r/t ice
2607 ! WS=RHgrd*EPS*ESI/(PP-ESI) ! Saturation mixing ratio
2608  ws=rhgrd*eps*esi/(pp+epsm1*esi) ! Saturation mixing ratio
2609  dwv=wvdum-ws ! Deficit grid-scale water vapor mixing ratio
2610  ssat=dwv/ws ! Supersaturation ratio
2611  deposit=0.
2612  DO WHILE (ssat > rhlimit .OR. ssat < rhlimit1)
2613  !
2614  !--- Note that XLVS2=LS*LV/(CP*RV)=LV*WS/(RV*T*T)*(LS/CP*DEP1),
2615  ! where WS is the saturation mixing ratio following Clausius-
2616  ! Clapeyron (see Asai,1965; Young,1993,p.405)
2617  !
2618  dep=dwv/(1.+xls2*ws/(tdum*tdum)) ! Asai (1965, J. Japan)
2619  tdum=tdum+xls1*dep ! Updated temperature
2620  wvdum=wvdum-dep ! Updated ice mixing ratio
2621  deposit=deposit+dep ! Total ice deposition
2622  esi=min(pp, fpvsi(tdum)) ! Updated saturation vapor press w/r/t ice
2623 ! WS=RHgrd*EPS*ESI/(PP-ESI) ! Updated saturation mixing ratio w/r/t ice
2624  ws=rhgrd*eps*esi/(pp+epsm1*esi) ! Updated saturation mixing ratio w/r/t ice
2625  dwv=wvdum-ws ! Deficit grid-scale water vapor mixing ratio
2626  ssat=dwv/ws ! Grid-scale supersaturation ratio
2627  ENDDO
2628  END FUNCTION deposit
2629 !
2630  END SUBROUTINE gsmcolumn
2631 
2632 
2633  SUBROUTINE rsipath(im, ix, ix2, levs, prsl, prsi, t, q, clw &
2634  &, f_ice, f_rain, f_rime, flgmin &
2635  &, cwatp, cicep, rainp, snowp &
2636  &, recwat, rerain, resnow, lprnt, ipr)
2637 !
2638  implicit none
2639 !
2640 !--------------------CLOUD----------------------------------------------
2641  integer im, ix, ix2, levs, ipr
2642  real prsl(ix,levs), prsi(ix,levs+1), t(ix,levs), q(ix,levs) &
2643  &, clw(ix2,levs), f_ice(ix2,levs), f_rain(ix2,levs) &
2644  &, f_rime(ix2,levs) &
2645  &, cwatp(ix,levs), rainp(ix,levs), cicep(ix,levs) &
2646  &, snowp(ix,levs), recwat(ix,levs), resnow(ix,levs) &
2647  &, rerain(ix,levs)
2648  real flgmin
2649  real frice, frrain, qcice, qcwat, qrain, qsnow, qtot, sden &
2650  &, cpath, rho, dsnow, flarge, rimef, xsimass, nlice &
2651  &, tc, recw1, drain, xli, dum, NLImax, pfac, pp &
2652  &, snofac, tem
2653 !
2654  real, parameter :: cexp=1./3.
2655  integer i, l, indexs
2656  logical lprnt
2657 !
2658 
2659  recw1 = 620.3505 / tnw**cexp ! cloud droplet effective radius
2660 
2661  do l=1,levs
2662  do i=1,im
2663  !--- HYDROMETEOR'S OPTICAL PATH
2664  cwatp(i,l) = 0.
2665  cicep(i,l) = 0.
2666  rainp(i,l) = 0.
2667  snowp(i,l) = 0.
2668  !--- HYDROMETEOR'S EFFECTIVE RADIUS
2669  recwat(i,l) = recwmin
2670  rerain(i,l) = rerainmin
2671  resnow(i,l) = resnowmin
2672  ENDDO
2673  ENDDO
2674 
2675  do l=1,levs
2676  DO i=1,im
2677 
2678  ! Assume parameterized condensate is
2679  ! all water for T>=-10C,
2680  ! all ice for T<=-30C,
2681  ! and a linear mixture at -10C > T > -30C
2682  !
2683  ! * Determine hydrometeor composition of total condensate (QTOT)
2684  !
2685 ! pp = prsl(i,l) * 1000.0
2686  pp = prsl(i,l) / prsi(i,levs+1)
2687 ! pfac = max(0.25, sqrt(sqrt(min(1.0, pp*0.000025))))
2688 ! pfac = max(0.5, sqrt(sqrt(min(1.0, pp*0.000025))))
2689 ! pfac = max(0.5, sqrt(sqrt(min(1.0, pp*0.00002))))
2690 ! pfac = max(0.25, sqrt(sqrt(min(1.0, pp*0.00001))))
2691 ! pfac = max(0.25, sqrt(sqrt(min(1.0, pp))))
2692 ! pfac = max(0.1, sqrt(min(1.0, pp*0.00001)))
2693 ! pfac = max(0.5, sqrt(sqrt(min(1.0, pp*0.000033))))
2694 ! pfac = max(0.5, sqrt(sqrt(min(1.0, pp*0.00004))))
2695 !go pfac = max(0.5, (sqrt(min(1.0, pp*0.000025))))
2696  pfac = 1.0
2697  tc = t(i,l) - t0c
2698  qtot = clw(i,l)
2699  IF (qtot > epsq) THEN
2700  qcwat=0.
2701  qcice=0.
2702  qrain=0.
2703  qsnow=0.
2704  frice = max(0.0, min(1.0, f_ice(i,l)))
2705  frrain = max(0.0, min(1.0, f_rain(i,l)))
2706  IF(tc <= thom) then
2707  qcice = qtot
2708  ELSE
2709  qcice = frice * qtot
2710  qcwat = qtot - qcice
2711  qrain = frrain * qcwat
2712  qcwat = qcwat - qrain
2713  ENDIF
2714  !
2715  !--- Air density (RHO), model mass thickness (CPATH)
2716  !
2717  rho = prsl(i,l)/(rd*t(i,l)*(1.+eps1*q(i,l)))
2718  cpath = (prsi(i,l+1)-prsi(i,l))*(1000000.0/grav)
2719 
2720  !! CLOUD WATER
2721  !
2722  !--- Effective radius (RECWAT) & total water path (CWATP)
2723  ! Assume monodisperse distribution of droplets (no factor of 1.5)
2724  !
2725  IF(qcwat > 0.) THEN
2726  recwat(i,l) = max(recwmin, recw1*(rho*qcwat)**cexp)
2727  cwatp(i,l) = cpath*qcwat ! cloud water path
2728 ! tem = 5.0*(1+max(0.0,min(1.0,-0.05*tc)))
2729 ! RECWAT(I,L) = max(RECWAT(I,L), tem)
2730  ENDIF
2731 
2732  !! RAIN
2733  !
2734  !--- Effective radius (RERAIN) & total water path (RAINP)
2735  !--- Factor of 1.5 accounts for r**3/r**2 moments for exponentially
2736  ! distributed drops in effective radius calculations
2737  ! (from M.D. Chou's code provided to Y.-T. Hou)
2738  !
2739  IF(qrain > 0.) THEN
2740  drain = cn0r0*sqrt(sqrt((rho*qrain)))
2741  rerain(i,l) = 1.5*max(xmrmin, min(drain, xmrmax))
2742  rainp(i,l) = cpath*qrain ! rain water path
2743  ENDIF
2744 
2745  !! SNOW (large ice) & CLOUD ICE
2746  !
2747  !--- Effective radius (RESNOW) & total ice path (SNOWP)
2748  !--- Total ice path (CICEP) for cloud ice
2749  !--- Factor of 1.5 accounts for r**3/r**2 moments for exponentially
2750  ! distributed ice particles in effective radius calculations
2751  !
2752  !--- Separation of cloud ice & "snow" uses algorithm from
2753  ! subroutine GSMCOLUMN
2754  !
2755  IF(qcice > 0.) THEN
2756  !
2757  !--- Mean particle size following Houze et al. (JAS, 1979, p. 160),
2758  ! converted from Fig. 5 plot of LAMDAs. An analogous set of
2759  ! relationships also shown by Fig. 8 of Ryan (BAMS, 1996, p. 66),
2760  ! but with a variety of different relationships that parallel the
2761  ! Houze curves.
2762  !
2763 ! DUM=MAX(0.05, MIN(1., EXP(.0536*TC)) )
2764  dum=max(0.05, min(1., exp(.0564*tc)) )
2765  indexs=min(mdimax, max(mdimin, int(xmimax*dum) ) )
2766 ! indexs=max(INDEXSmin, indexs)
2767 ! NLImax=5.E3/sqrt(DUM) !- Ver3
2768  dum=max(flgmin*pfac, dum)
2769 ! DUM=MAX(FLGmin, DUM)
2770 ! NLImax=20.E3 !- Ver3
2771 ! NLImax=50.E3 !- Ver3 => comment this line out
2772  nlimax=10.e3/sqrt(dum) !- Ver3
2773 ! NLImax=5.E3/sqrt(DUM) !- Ver3
2774 ! NLImax=6.E3/sqrt(DUM) !- Ver3
2775 ! NLImax=7.5E3/sqrt(DUM) !- Ver3
2776 ! NLImax=20.E3/DUM !- Ver3
2777 ! NLImax=20.E3/max(0.2,DUM) !- Ver3
2778 ! NLImax=2.0E3/max(0.1,DUM) !- Ver3
2779 ! NLImax=2.5E3/max(0.1,DUM) !- Ver3
2780 ! NLImax=10.E3/max(0.2,DUM) !- Ver3
2781 ! NLImax=4.E3/max(0.2,DUM) !- Ver3
2782 !Moorthi DSNOW = XMImax*EXP(.0536*TC)
2783 !Moorthi INDEXS = MAX(INDEXSmin, MIN(MDImax, INT(DSNOW)))
2784  !
2785  !--- Assumed number fraction of large ice to total (large & small)
2786  ! ice particles, which is based on a general impression of the
2787  ! literature.
2788  !
2789  ! Small ice are assumed to have a mean diameter of 50 microns.
2790  !
2791  IF(tc >= 0.) THEN
2792  flarge=flg1p0
2793  ELSE
2794  flarge = dum
2795  ENDIF
2796 !------------------------Commented by Moorthi -----------------------------
2797 ! ELSEIF (TC >= -25.) THEN
2798 !
2799 !--- Note that absence of cloud water (QCWAT) is used as a quick
2800 ! substitute for calculating water subsaturation as in GSMCOLUMN
2801 !
2802 ! IF(QCWAT <= 0. .OR. TC < -8.
2803 ! & .OR. TC > -3.)THEN
2804 ! FLARGE=FLG0P2
2805 ! ELSE
2806 
2807 !--- Parameterize effects of rime splintering by increasing
2808 ! number of small ice particles
2809 !
2810 ! FLARGE=FLG0P1
2811 ! ENDIF
2812 ! ELSEIF (TC <= -50.) THEN
2813 ! FLARGE=.01
2814 ! ELSE
2815 ! FLARGE=.2*EXP(.1198*(TC+25))
2816 ! ENDIF
2817 !____________________________________________________________________________
2818 
2819  rimef=max(1., f_rime(i,l))
2820  xsimass=massi(mdimin)*(1.-flarge)/flarge
2821 ! if (lprnt) print *,' rimef=',rimef,' xsimass=',xsimass
2822 ! &,' indexs=',indexs,' massi=',massi(indexs),' flarge=',flarge
2823  nlice=rho*qcice/(xsimass+rimef*massi(indexs))
2824  !
2825  !--- From subroutine GSMCOLUMN:
2826  !--- Minimum number concentration for large ice of NLImin=10/m**3
2827  ! at T>=0C. Done in order to prevent unrealistically small
2828  ! melting rates and tiny amounts of snow from falling to
2829  ! unrealistically warm temperatures.
2830  !
2831  IF(tc >= 0.) THEN
2832  nlice=max(nlimin, nlice)
2833  ELSEIF (nlice > nlimax) THEN
2834  !
2835  !--- Ferrier 6/13/01: Prevent excess accumulation of ice
2836  !
2837  xli=(rho*qcice/nlimax-xsimass)/rimef
2838 
2839  IF(xli <= massi(450) ) THEN
2840  dsnow=9.5885e5*xli**.42066
2841  ELSE
2842  dsnow=3.9751e6*xli**.49870
2843  ENDIF
2844 
2845  indexs=min(mdimax, max(indexs, int(dsnow)))
2846  nlice=rho*qcice/(xsimass+rimef*massi(indexs))
2847  ENDIF
2848 
2849 ! if (tc > -20.0 .and. indexs >= indexsmin) then
2850 ! snofac = max(0.0, min(1.0, exp(1.0*(tc+20.0))))
2851 ! if (indexs >= indexsmin) then
2852 ! if (tc > -20.0 .or. indexs >= indexsmin) then
2853 ! if (tc > -40.0) then
2854 ! if (tc >= -40.0 .or. prsl(i,l) > 50.0) then
2855 !! if (tc >= -20.0) then
2856 ! if (tc >= -20.0 .or. prsl(i,l) > 50.0) then
2857 ! if ((tc >= -20.0 .or.
2858 ! & prsi(i,levs+1)-prsi(i,l) < 30.0)
2859  if (prsi(i,levs+1)-prsi(i,l) < 40.0 &
2860 ! if (prsi(i,levs+1)-prsi(i,l) < 70.0
2861  & .and. indexs >= indexsmin) then
2862 ! & prsi(i,levs)-prsl(i,l) < 20.0) then
2863 ! & prsi(i,levs)-prsl(i,l) < 30.0) then
2864 ! & prsi(i,levs)-prsl(i,l) < 40.0) then
2865 ! snofac = max(0.0, min(1.0, 0.05*(tc+40.0)))
2866 ! snofac = max(0.0, min(1.0, 0.1*(tc+25.0)))
2867 ! snofac = max(0.0, min(1.0, 0.0667*(tc+25.0)))
2868 ! if (indexs > indexsmin) then
2869  qsnow = min(qcice, nlice*rimef*massi(indexs)/rho)
2870 ! & * snofac
2871  endif
2872 ! qsnow = qcice
2873  qcice = max(0., qcice-qsnow)
2874 ! qsnow = 0.0
2875  cicep(i,l) = cpath*qcice ! cloud ice path
2876  resnow(i,l) = 1.5*float(indexs)
2877  sden = sdens(indexs)/rimef ! 1/snow density
2878  snowp(i,l) = cpath*qsnow*sden ! snow path / snow density
2879 ! SNOWP (I,L) = CPATH*QSNOW ! snow path / snow density
2880 ! if (lprnt .and. i == ipr) then
2881 ! print *,' L=',L,' snowp=',snowp(i,l),' cpath=',cpath
2882 ! &,' qsnow=',qsnow,' sden=',sden,' rimef=',rimef,' indexs=',indexs
2883 ! &,' sdens=',sdens(indexs),' resnow=',resnow(i,l)
2884 ! &,' qcice=',qcice,' cicep=',cicep(i,l)
2885 ! endif
2886 
2887 
2888  ENDIF ! END QCICE BLOCK
2889  ENDIF ! QTOT IF BLOCK
2890 
2891  ENDDO
2892  ENDDO
2893 !
2894  END SUBROUTINE rsipath
2895 
2896 
2925 !-----------------------------------
2926  subroutine rsipath2 &
2927 !...................................
2928 
2929 ! --- inputs:
2930  & ( plyr, plvl, tlyr, qlyr, qcwat, qcice, qrain, rrime, &
2931  & im, levs, iflip, flgmin, &
2932 ! --- outputs:
2933  & cwatp, cicep, rainp, snowp, recwat, rerain, resnow, snden &
2934  & )
2935 
2936 ! ================= subprogram documentation block ================ !
2937 ! !
2938 ! abstract: this program is a modified version of ferrier's original !
2939 ! "rsipath" subprogram. it computes layer's cloud liquid, ice, rain, !
2940 ! and snow water condensate path and the partical effective radius !
2941 ! for liquid droplet, rain drop, and snow flake. !
2942 ! !
2943 ! ==================== defination of variables ==================== !
2944 ! !
2945 ! input variables: !
2946 ! plyr (IM,LEVS) : model layer mean pressure in mb (100Pa) !
2947 ! plvl (IM,LEVS+1):model level pressure in mb (100Pa) !
2948 ! tlyr (IM,LEVS) : model layer mean temperature in k !
2949 ! qlyr (IM,LEVS) : layer specific humidity in gm/gm !
2950 ! qcwat (IM,LEVS) : layer cloud liquid water condensate amount !
2951 ! qcice (IM,LEVS) : layer cloud ice water condensate amount !
2952 ! qrain (IM,LEVS) : layer rain drop water amount !
2953 ! rrime (IM,LEVS) : mass ratio of total to unrimed ice ( >= 1 ) !
2954 ! IM : horizontal dimention !
2955 ! LEVS : vertical layer dimensions !
2956 ! iflip : control flag for in/out vertical indexing !
2957 ! =0: index from toa to surface !
2958 ! =1: index from surface to toa !
2959 ! flgmin : Minimum large ice fraction !
2960 ! lprnt : logical check print control flag !
2961 ! !
2962 ! output variables: !
2963 ! cwatp (IM,LEVS) : layer cloud liquid water path !
2964 ! cicep (IM,LEVS) : layer cloud ice water path !
2965 ! rainp (IM,LEVS) : layer rain water path !
2966 ! snowp (IM,LEVS) : layer snow water path !
2967 ! recwat(IM,LEVS) : layer cloud eff radius for liqid water (micron) !
2968 ! rerain(IM,LEVS) : layer rain water effective radius (micron) !
2969 ! resnow(IM,LEVS) : layer snow flake effective radius (micron) !
2970 ! snden (IM,LEVS) : 1/snow density !
2971 ! !
2972 ! !
2973 ! usage: call rsipath2 !
2974 ! !
2975 ! subroutines called: none !
2976 ! !
2977 ! program history log: !
2978 ! xx-xx-2001 b. ferrier - original program !
2979 ! xx-xx-2004 s. moorthi - modified for use in gfs model !
2980 ! 05-20-2004 y. hou - modified, added vertical index flag!
2981 ! to reduce data flipping, and rearrange code to !
2982 ! be comformable with radiation part programs. !
2983 ! !
2984 ! ==================== end of description ===================== !
2985 !
2986 
2987  implicit none
2988 
2989 ! --- constant parameter:
2990  real, parameter :: CEXP= 1.0/3.0
2991 
2992 ! --- inputs:
2993  real, dimension(:,:), intent(in) :: &
2994  & plyr, plvl, tlyr, qlyr, qcwat, qcice, qrain, rrime
2995 
2996  integer, intent(in) :: IM, LEVS, iflip
2997  real, dimension(:), intent(in) :: flgmin
2998 ! logical, intent(in) :: lprnt
2999 
3000 ! --- output:
3001  real, dimension(:,:), intent(out) :: &
3002  & cwatp, cicep, rainp, snowp, recwat, rerain, resnow, snden
3003 
3004 ! --- locals:
3005 ! real, dimension(IM,LEVS) :: delp, pp1, pp2
3006 
3007  real :: recw1, dsnow, qsnow, qqcice, flarge, xsimass, pfac, &
3008  & nlice, xli, nlimax, dum, tem, &
3009  & rho, cpath, rc, totcnd, tc
3010 
3011  integer :: i, k, indexs, ksfc, k1
3012 !
3013 !===> ... begin here
3014 !
3015  recw1 = 620.3505 / tnw**cexp ! cloud droplet effective radius
3016 
3017  do k = 1, levs
3018  do i = 1, im
3019  !--- hydrometeor's optical path
3020  cwatp(i,k) = 0.0
3021  cicep(i,k) = 0.0
3022  rainp(i,k) = 0.0
3023  snowp(i,k) = 0.0
3024  snden(i,k) = 0.0
3025  !--- hydrometeor's effective radius
3026  recwat(i,k) = recwmin
3027  rerain(i,k) = rerainmin
3028  resnow(i,k) = resnowmin
3029  enddo
3030  enddo
3031 
3034 ! --- set up pressure related arrays, convert unit from mb to cb (10Pa)
3035 ! cause the rest part uses cb in computation
3036 
3037  if (iflip == 0) then ! data from toa to sfc
3038  ksfc = levs + 1
3039  k1 = 0
3040  else ! data from sfc to top
3041  ksfc = 1
3042  k1 = 1
3043  endif ! end_if_iflip
3044 !
3045  do k = 1, levs
3046  do i = 1, im
3047  totcnd = qcwat(i,k) + qcice(i,k) + qrain(i,k)
3048  qsnow = 0.0
3049  if(totcnd > epsq) then
3051 ! --- air density (rho), model mass thickness (cpath), temperature in c (tc)
3052 
3053  rho = 0.1 * plyr(i,k) &
3054  & / (rd* tlyr(i,k) * (1.0 + eps1*qlyr(i,k)))
3055  cpath = abs(plvl(i,k+1) - plvl(i,k)) * (100000.0 / grav)
3056  tc = tlyr(i,k) - t0c
3057 
3060 ! cloud water
3061 !
3062 ! --- effective radius (recwat) & total water path (cwatp):
3063 ! assume monodisperse distribution of droplets (no factor of 1.5)
3064 
3065  if (qcwat(i,k) > 0.0) then
3066  recwat(i,k) = max(recwmin,recw1*(rho*qcwat(i,k))**cexp)
3067  cwatp(i,k) = cpath * qcwat(i,k) ! cloud water path
3068 ! tem = 5.0*(1.0 + max(0.0, min(1.0,-0.05*tc)))
3069 ! recwat(i,k) = max(recwat(i,k), tem)
3070  endif
3074 ! rain
3075 !
3076 ! --- effective radius (rerain) & total water path (rainp):
3077 ! factor of 1.5 accounts for r**3/r**2 moments for exponentially
3078 ! distributed drops in effective radius calculations
3079 ! (from m.d. chou's code provided to y.-t. hou)
3080 
3081  if (qrain(i,k) > 0.0) then
3082  tem = cn0r0 * sqrt(sqrt(rho*qrain(i,k)))
3083  rerain(i,k) = 1.5 * max(xmrmin, min(xmrmax, tem))
3084  rainp(i,k) = cpath * qrain(i,k) ! rain water path
3085  endif
3086 
3087 ! snow (large ice) & cloud ice
3088 !
3089 ! --- effective radius (resnow) & total ice path (snowp) for snow, and
3090 ! total ice path (cicep) for cloud ice:
3091 ! factor of 1.5 accounts for r**3/r**2 moments for exponentially
3092 ! distributed ice particles in effective radius calculations
3093 ! separation of cloud ice & "snow" uses algorithm from subroutine gsmcolumn
3094 
3095 ! pfac = max(0.5, sqrt(sqrt(min(1.0, pp1(i,k)*0.00004))))
3096 !go pfac = max(0.5, (sqrt(min(1.0, pp1(i,k)*0.000025))))
3097  pfac = 1.0
3098 
3099  if (qcice(i,k) > 0.0) then
3100 
3106 ! --- mean particle size following houze et al. (jas, 1979, p. 160),
3107 ! converted from fig. 5 plot of lamdas. an analogous set of
3108 ! relationships also shown by fig. 8 of ryan (bams, 1996, p. 66),
3109 ! but with a variety of different relationships that parallel
3110 ! the houze curves.
3111 
3112 ! dum = max(0.05, min(1.0, exp(0.0536*tc) ))
3113  dum = max(0.05, min(1.0, exp(0.0564*tc) ))
3114  indexs = min(mdimax, max(mdimin, int(xmimax*dum) ))
3115  dum=max(flgmin(i)*pfac, dum)
3116 
3117 !!\n - assumed number of fraction of large ice to total (large+small) ice
3118 !! particles, which is based on a general impression of the literature. small ice
3119 !! are assumed to have a mean diameter of 50 microns
3120 ! --- assumed number fraction of large ice to total (large & small) ice
3121 ! particles, which is based on a general impression of the literature.
3122 ! small ice are assumed to have a mean diameter of 50 microns.
3123 
3124  if (tc >= 0.0) then
3125  flarge = flg1p0
3126  else
3127  flarge = dum
3128 ! flarge = max(FLGmin*pfac, dum)
3129  endif
3130 !------------------------commented by moorthi -----------------------------
3131 ! elseif (tc >= -25.0) then
3132 !
3133 ! --- note that absence of cloud water (qcwat) is used as a quick
3134 ! substitute for calculating water subsaturation as in gsmcolumn
3135 !
3136 ! if (qcwat(i,k) <= 0.0 .or. tc < -8.0 &
3137 ! & .or. tc > -3.0) then
3138 ! flarge = FLG0P2
3139 ! else
3140 !
3141 ! --- parameterize effects of rime splintering by increasing
3142 ! number of small ice particles
3143 !
3144 ! flarge = FLG0P1
3145 ! endif
3146 ! elseif (tc <= -50.0) then
3147 ! flarge = 0.01
3148 ! else
3149 ! flarge = 0.2 * exp(0.1198*(tc+25.0))
3150 ! endif
3151 !____________________________________________________________________________
3152 
3153  xsimass = massi(mdimin) * (1.0 - flarge) / flarge
3154 ! nlimax = 20.0e3 !- ver3
3155 ! NLImax=50.E3 !- Ver3 => comment this line out
3156  nlimax=10.e3/sqrt(dum) !- Ver3
3157 ! NLImax=5.E3/sqrt(DUM) !- Ver3
3158 ! NLImax=6.E3/sqrt(DUM) !- Ver3
3159 ! NLImax=7.5E3/sqrt(DUM) !- Ver3
3160 
3161 ! indexs = min(MDImax, max(MDImin, int(XMImax*dum) ))
3162 !moorthi dsnow = XMImax * exp(0.0536*tc)
3163 !moorthi indexs = max(INDEXSmin, min(MDImax, int(dsnow)))
3164 
3165 ! if (lprnt) print *,' rrime=',rrime,' xsimass=',xsimass, &
3166 ! & ' indexs=',indexs,' massi=',massi(indexs),' flarge=',flarge
3167 
3168  tem = rho * qcice(i,k)
3169  nlice = tem / (xsimass +rrime(i,k)*massi(indexs))
3170 
3171 ! --- from subroutine gsmcolumn:
3172 ! minimum number concentration for large ice of NLImin=10/m**3
3173 ! at t>=0c. done in order to prevent unrealistically small
3174 ! melting rates and tiny amounts of snow from falling to
3175 ! unrealistically warm temperatures.
3176 
3177  if (tc >= 0.0) then
3178 
3179  nlice = max(nlimin, nlice)
3180 
3181  elseif (nlice > nlimax) then
3182 
3183 ! --- ferrier 6/13/01: prevent excess accumulation of ice
3184 
3185  xli = (tem/nlimax - xsimass) / rrime(i,k)
3186 
3187  if (xli <= massi(450) ) then
3188  dsnow = 9.5885e5 * xli**0.42066
3189  else
3190  dsnow = 3.9751e6 * xli** 0.49870
3191  endif
3192 
3193  indexs = min(mdimax, max(indexs, int(dsnow)))
3194  nlice = tem / (xsimass + rrime(i,k)*massi(indexs))
3195 
3196  endif ! end if_tc block
3197 
3198 ! if (abs(plvl(i,ksfc)-plvl(i,k+k1)) < 300.0 &
3199 ! if (abs(plvl(i,ksfc)-plvl(i,k+k1)) < 400.0 &
3200 ! if (plvl(i,k+k1) > 600.0 &
3201 ! & .and. indexs >= INDEXSmin) then
3202 ! if (tc > -20.0 .and. indexs >= indexsmin) then
3203  if (plvl(i,ksfc) > 850.0 .and. &
3204 ! & plvl(i,k+k1) > 600.0 .and. indexs >= indexsmin) then
3205  & plvl(i,k+k1) > 700.0 .and. indexs >= indexsmin) then ! 20060516
3206 !! if (plvl(i,ksfc) > 800.0 .and. &
3207 !! & plvl(i,k+k1) > 700.0 .and. indexs >= indexsmin) then
3208 ! if (plvl(i,ksfc) > 700.0 .and. &
3209 ! & plvl(i,k+k1) > 600.0 .and. indexs >= indexsmin) then
3210  qsnow = min( qcice(i,k), &
3211  & nlice*rrime(i,k)*massi(indexs)/rho )
3212  endif
3213 
3214  qqcice = max(0.0, qcice(i,k)-qsnow)
3215  cicep(i,k) = cpath * qqcice ! cloud ice path
3216  resnow(i,k) = 1.5 * float(indexs)
3217  snden(i,k) = sdens(indexs) / rrime(i,k) ! 1/snow density
3218  snowp(i,k) = cpath*qsnow ! snow path
3219 ! snowp (i,k) = cpath*qsnow*snden(i,k) ! snow path / snow density
3220 
3221 ! if (lprnt .and. i == ipr) then
3222 ! if (i == 2) then
3223 ! print *,' L=',k,' snowp=',snowp(i,k),' cpath=',cpath, &
3224 ! & ' qsnow=',qsnow,' sden=',snden(i,k),' rrime=',rrime(i,k),&
3225 ! & ' indexs=',indexs,' sdens=',sdens(indexs),' resnow=', &
3226 ! & resnow(i,k),' qcice=',qqcice,' cicep=',cicep(i,k)
3227 ! endif
3228 
3229  endif ! end if_qcice block
3230  endif ! end if_totcnd block
3231 
3232  enddo
3233  enddo
3234 !
3235 !...................................
3236  end subroutine rsipath2
3237 !-----------------------------------
3239  end MODULE module_microphysics
3240 
real, parameter, private rerainmin
real, parameter, private dmimax
integer, parameter itlo
real function condense(PP, QW, RHgrd, TK, WV)
integer, dimension(itlo:ithi, 4) nstats
subroutine rsipath2
This subroutine is a modified version of ferrier&#39;s original "rsipath" subprogram. It computes layer&#39;s...
integer, private mic_step
real(kind=kind_phys), parameter con_t0c
temp at 0C (K)
Definition: physcons.f:89
real(kind=kind_phys), parameter con_hfus
lat heat H2O fusion (J/kg)
Definition: physcons.f:85
real(kind=kind_phys), parameter con_hvap
lat heat H2O cond (J/kg)
Definition: physcons.f:83
integer, parameter, private nrime
real, parameter, private tnw
This module contains some the most frequently used math and physics constants for gcm models...
Definition: physcons.f:29
real function deposit(PP, RHgrd, Tdum, WVdum)
real, dimension(2:9, 0:nrime), private vel_rf
subroutine rsipath(im, ix, ix2, levs, prsl, prsi, t, q, clw
real, dimension(itlo:ithi, 5) qmax
real, parameter, private thom
subroutine my_growth_rates(DTPH)
integer, parameter, private mdimin
integer, parameter, private my_t2
real(kind=kind_phys), parameter con_rd
gas constant air (J/kg/K)
Definition: physcons.f:69
real, parameter, private dmrmax
subroutine gsmconst(DTPG, mype, first)
integer, parameter, private indexsmin
real, dimension(itlo:ithi, 22) qtot
integer, parameter, private my_t1
subroutine gsmcolumn(ARAING, ASNOWG, DTPG, I_index, J_index,
integer, parameter, private mdrmin
real, parameter, private dmrmin
integer, parameter ithi
real(kind=kind_phys), parameter con_rv
gas constant H2O (J/kg/K)
Definition: physcons.f:71
real(kind=kind_phys), parameter con_g
gravity (m/s2)
Definition: physcons.f:52
real(kind=kind_phys), parameter con_pi
Definition: physcons.f:41
integer, parameter, private mdrmax
real, parameter, private dmimin
real, dimension(my_t1:my_t2), private my_growth
integer, parameter, private mdimax
real(kind=kind_phys), parameter con_cp
spec heat air at p (J/kg/K)
Definition: physcons.f:73
real, parameter, private toler
real(kind=kind_phys), parameter con_fvirt
Definition: physcons.f:107
real(kind=kind_phys), parameter con_epsm1
Definition: physcons.f:109
real(kind=kind_phys), parameter con_eps
Definition: physcons.f:108